This vignette shows how to build a SpatialExperiment (SPE) from: - Nanostring CosMx (RNA/protein) outputs - 10x Genomics Xenium outputs
For each technology, we demonstrate: - Route A: SpaceTrooper’s high-level reader (readCosmxSPE, readCosmxProteinSPE, readXeniumSPE) - Route B: read with SpatialExperimentIO, then standardize with updateCosmxSPE / updateXeniumSPE
We begin with CosMx. The package ships with a small CosMx example for demonstration.
cospath <- system.file(file.path("extdata", "CosMx_DBKero_Tiny"), package="SpaceTrooper")
cospath
#> [1] "/tmp/RtmpaWTLmh/Rinst2adb0d52563f17/SpaceTrooper/extdata/CosMx_DBKero_Tiny"
Use readCosmxSPE() to construct an SPE from CosMx outputs; it also normalizes names/metadata and records polygons/FOV info if present.
spe_cos <- readCosmxSPE(
dirName=cospath,
sampleName="DBKero_Tiny",
coordNames=c("CenterX_global_px", "CenterY_global_px"),
countMatFPattern="exprMat_file.csv",
metadataFPattern="metadata_file.csv",
polygonsFPattern="polygons.csv",
fovPosFPattern="fov_positions_file.csv",
fovdims=c(xdim=4256, ydim=4256)
)
spe_cos
#> class: SpatialExperiment
#> dim: 1010 905
#> metadata(4): fov_positions fov_dim polygons technology
#> assays(1): counts
#> rownames(1010): RAMP2 CD83 ... NegPrb09 NegPrb10
#> rowData names(0):
#> colnames(905): f16_c1 f16_c10 ... f16_c98 f16_c99
#> colData names(20): fov cellID ... sample_id cell_id
#> reducedDimNames(0):
#> mainExpName: NULL
#> altExpNames(0):
#> spatialCoords names(2) : CenterX_global_px CenterY_global_px
#> imgData names(1): sample_id
Inspect essentials:
assayNames(spe_cos)
#> [1] "counts"
dim(spe_cos)
#> [1] 1010 905
head(colnames(spatialCoords(spe_cos)))
#> [1] "CenterX_global_px" "CenterY_global_px"
metadata(spe_cos)$technology
#> [1] "Nanostring_CosMx"
metadata(spe_cos)$polygons
#> [1] "/tmp/RtmpaWTLmh/Rinst2adb0d52563f17/SpaceTrooper/extdata/CosMx_DBKero_Tiny/DBKero-polygons.csv"
If working with CosMx Protein, use the convenience wrapper:
protfolder <- system.file("extdata", "S01_prot", package = "SpaceTrooper")
spe_cos_prot <- readCosmxProteinSPE(
dirName=protfolder,
sampleName="cosmx_prots",
coordNames=c("CenterX_global_px", "CenterY_global_px"),
countMatFPattern="exprMat_file.csv",
metadataFPattern="metadata_file.csv",
polygonsFPattern="polygons.csv",
fovPosFPattern="fov_positions_file.csv",
fovdims=c(xdim=4256, ydim=4256)
)
metadata(spe_cos_prot)$technology
If you prefer to read with SpatialExperimentIO first, upgrade the object with updateCosmxSPE() to harmonize names/metadata and attach polygons.
spe_cos_raw <- SpatialExperimentIO::readCosmxSXE(
dirName=cospath,
returnType="SPE",
countMatPattern="exprMat_file.csv",
metaDataPattern="metadata_file.csv",
coordNames=c("CenterX_global_px", "CenterY_global_px"),
addFovPos=TRUE,
fovPosPattern="fov_positions_file.csv",
altExps=NULL,
addParquetPaths=FALSE
)
spe_cos_std <- updateCosmxSPE(
spe=spe_cos_raw,
dirName=cospath,
sampleName="DBKero_Tiny",
polygonsFPattern="polygons.csv",
fovdims=c(xdim=4256, ydim=4256)
)
identical(spe_cos_std, spe_cos)
#> [1] TRUE
A small Xenium example is also included for demonstration.
xepath <- system.file("extdata", "Xenium_small", package = "SpaceTrooper")
xepath
#> [1] "/tmp/RtmpaWTLmh/Rinst2adb0d52563f17/SpaceTrooper/extdata/Xenium_small"
readXeniumSPE() builds the SPE from a Xenium Output Bundle (root or outs/) and standardizes metadata.
Key options: - type: "HDF5" or "sparse" (feature matrix) - boundariesType: "parquet" or "csv" (cell boundaries) - computeMissingMetrics: compute QC metrics (area/aspect ratio) if needed - keepPolygons: append polygons to colData - addFOVs: derive FOV IDs from transcript parquet
spe_xen_a <- readXeniumSPE(
dirName=xepath,
type="HDF5",
coordNames=c("x_centroid", "y_centroid"),
boundariesType="parquet",
computeMissingMetrics=TRUE,
keepPolygons=TRUE,
countsFilePattern="cell_feature_matrix",
metadataFPattern="cells.csv.gz",
polygonsFPattern="cell_boundaries",
polygonsCol="polygons",
txPattern="transcripts",
addFOVs=FALSE
)
#> Computing missing metrics, this could take some time...
spe_xen_a
#> class: SpatialExperiment
#> dim: 4 6
#> metadata(2): polygons technology
#> assays(1): counts
#> rownames(4): ABCC11 ACTA2 ACTG2 ADAM9
#> rowData names(3): ID Symbol Type
#> colnames(6): 1 2 ... 5 6
#> colData names(11): X cell_id ... Area_um polygons
#> reducedDimNames(0):
#> mainExpName: NULL
#> altExpNames(0):
#> spatialCoords names(2) : x_centroid y_centroid
#> imgData names(1): sample_id
Quick checks:
assayNames(spe_xen_a)
#> [1] "counts"
dim(spe_xen_a)
#> [1] 4 6
colnames(spatialCoords(spe_xen_a))
#> [1] "x_centroid" "y_centroid"
metadata(spe_xen_a)$polygons
#> [1] "/tmp/RtmpaWTLmh/Rinst2adb0d52563f17/SpaceTrooper/extdata/Xenium_small/cell_boundaries.parquet"
metadata(spe_xen_a)$technology
#> [1] "10X_Xenium"
Read with SpatialExperimentIO and then pass through updateXeniumSPE() for SpaceTrooper-standardized metadata and optional metrics/FOV extraction.
spe_xen_b <- SpatialExperimentIO::readXeniumSXE(
dirName=xepath,
countMatPattern="cell_feature_matrix.h5",
metaDataPattern="cells.csv.gz",
coordNames=c("x_centroid", "y_centroid"),
returnType="SPE",
addExperimentXenium=FALSE,
altExps=NULL,
addParquetPaths=FALSE
)
spe_xen_b <- updateXeniumSPE(
spe=spe_xen_b,
dirName=xepath,
boundariesType="parquet",
computeMissingMetrics=TRUE,
keepPolygons=TRUE,
polygonsFPattern="cell_boundaries",
polygonsCol="polygons",
txPattern="transcripts",
addFOVs=FALSE
)
#> Computing missing metrics, this could take some time...
spe_xen_b
#> class: SpatialExperiment
#> dim: 4 6
#> metadata(2): polygons technology
#> assays(1): counts
#> rownames(4): ABCC11 ACTA2 ACTG2 ADAM9
#> rowData names(3): ID Symbol Type
#> colnames(6): 1 2 ... 5 6
#> colData names(11): X cell_id ... Area_um polygons
#> reducedDimNames(0):
#> mainExpName: NULL
#> altExpNames(0):
#> spatialCoords names(2) : x_centroid y_centroid
#> imgData names(1): sample_id
Validate:
identical(metadata(spe_xen_b)$technology, "10X_Xenium")
#> [1] TRUE
identical(spe_xen_a, spe_xen_b)
#> [1] TRUE
metadata_file.csv), expression matrix (exprMat_file.csv), optional polygon CSV (polygons.csv), and FOV positions. updateCosmxSPE() also fixes common field names and records FOV dims in metadata.outs/ folder. Feature matrix may be cell_feature_matrix.h5 (HDF5) or sparse folder; cells metadata cells.csv.gz; boundaries as .parquet or .csv.gz; transcript parquet for FOV attribution. readXeniumSPE() auto-detects outs/ if you pass the bundle root.sessionInfo()
#> R Under development (unstable) (2026-01-15 r89304)
#> Platform: x86_64-pc-linux-gnu
#> Running under: Ubuntu 24.04.3 LTS
#>
#> Matrix products: default
#> BLAS: /home/biocbuild/bbs-3.23-bioc/R/lib/libRblas.so
#> LAPACK: /usr/lib/x86_64-linux-gnu/lapack/liblapack.so.3.12.0 LAPACK version 3.12.0
#>
#> locale:
#> [1] LC_CTYPE=en_US.UTF-8 LC_NUMERIC=C
#> [3] LC_TIME=en_GB LC_COLLATE=C
#> [5] LC_MONETARY=en_US.UTF-8 LC_MESSAGES=en_US.UTF-8
#> [7] LC_PAPER=en_US.UTF-8 LC_NAME=C
#> [9] LC_ADDRESS=C LC_TELEPHONE=C
#> [11] LC_MEASUREMENT=en_US.UTF-8 LC_IDENTIFICATION=C
#>
#> time zone: America/New_York
#> tzcode source: system (glibc)
#>
#> attached base packages:
#> [1] stats4 stats graphics grDevices utils datasets methods
#> [8] base
#>
#> other attached packages:
#> [1] ggplot2_4.0.1 SpaceTrooper_1.1.4
#> [3] SpatialExperiment_1.21.0 SingleCellExperiment_1.33.0
#> [5] SummarizedExperiment_1.41.0 Biobase_2.71.0
#> [7] GenomicRanges_1.63.1 Seqinfo_1.1.0
#> [9] IRanges_2.45.0 S4Vectors_0.49.0
#> [11] BiocGenerics_0.57.0 generics_0.1.4
#> [13] MatrixGenerics_1.23.0 matrixStats_1.5.0
#> [15] BiocStyle_2.39.0
#>
#> loaded via a namespace (and not attached):
#> [1] RColorBrewer_1.1-3 jsonlite_2.0.0
#> [3] shape_1.4.6.1 magrittr_2.0.4
#> [5] ggbeeswarm_0.7.3 magick_2.9.0
#> [7] farver_2.1.2 rmarkdown_2.30
#> [9] vctrs_0.7.1 DelayedMatrixStats_1.33.0
#> [11] tinytex_0.58 rstatix_0.7.3
#> [13] htmltools_0.5.9 S4Arrays_1.11.1
#> [15] BiocNeighbors_2.5.2 broom_1.0.11
#> [17] Rhdf5lib_1.33.0 SparseArray_1.11.10
#> [19] Formula_1.2-5 rhdf5_2.55.12
#> [21] sass_0.4.10 KernSmooth_2.23-26
#> [23] bslib_0.10.0 cachem_1.1.0
#> [25] lifecycle_1.0.5 iterators_1.0.14
#> [27] pkgconfig_2.0.3 rsvd_1.0.5
#> [29] Matrix_1.7-4 R6_2.6.1
#> [31] fastmap_1.2.0 digest_0.6.39
#> [33] scater_1.39.2 dqrng_0.4.1
#> [35] irlba_2.3.5.1 ggpubr_0.6.2
#> [37] beachmat_2.27.2 labeling_0.4.3
#> [39] SpatialExperimentIO_1.3.0 abind_1.4-8
#> [41] compiler_4.6.0 proxy_0.4-29
#> [43] bit64_4.6.0-1 withr_3.0.2
#> [45] S7_0.2.1 backports_1.5.0
#> [47] BiocParallel_1.45.0 carData_3.0-5
#> [49] viridis_0.6.5 DBI_1.2.3
#> [51] HDF5Array_1.39.0 R.utils_2.13.0
#> [53] ggsignif_0.6.4 DelayedArray_0.37.0
#> [55] rjson_0.2.23 classInt_0.4-11
#> [57] tools_4.6.0 units_1.0-0
#> [59] vipor_0.4.7 otel_0.2.0
#> [61] beeswarm_0.4.0 R.oo_1.27.1
#> [63] glue_1.8.0 h5mread_1.3.1
#> [65] rhdf5filters_1.23.3 grid_4.6.0
#> [67] sf_1.0-24 gtable_0.3.6
#> [69] R.methodsS3_1.8.2 class_7.3-23
#> [71] tidyr_1.3.2 data.table_1.18.0
#> [73] BiocSingular_1.27.1 ScaledMatrix_1.19.0
#> [75] car_3.1-3 XVector_0.51.0
#> [77] ggrepel_0.9.6 foreach_1.5.2
#> [79] pillar_1.11.1 limma_3.67.0
#> [81] robustbase_0.99-6 splines_4.6.0
#> [83] dplyr_1.1.4 lattice_0.22-7
#> [85] survival_3.8-6 bit_4.6.0
#> [87] tidyselect_1.2.1 locfit_1.5-9.12
#> [89] scuttle_1.21.0 sfheaders_0.4.5
#> [91] knitr_1.51 gridExtra_2.3
#> [93] bookdown_0.46 edgeR_4.9.2
#> [95] xfun_0.56 statmod_1.5.1
#> [97] DropletUtils_1.31.0 DEoptimR_1.1-4
#> [99] yaml_2.3.12 evaluate_1.0.5
#> [101] codetools_0.2-20 tibble_3.3.1
#> [103] BiocManager_1.30.27 cli_3.6.5
#> [105] arrow_23.0.0 jquerylib_0.1.4
#> [107] dichromat_2.0-0.1 Rcpp_1.1.1
#> [109] parallel_4.6.0 assertthat_0.2.1
#> [111] sparseMatrixStats_1.23.0 glmnet_4.1-10
#> [113] viridisLite_0.4.2 scales_1.4.0
#> [115] e1071_1.7-17 purrr_1.2.1
#> [117] rlang_1.1.7 cowplot_1.2.0