1 Overview

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

2 CosMx (Nanostring)

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/RtmpBasoI1/Rinst33659443778b2d/SpaceTrooper/extdata/CosMx_DBKero_Tiny"

2.1 Route A — Direct loading with SpaceTrooper

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/RtmpBasoI1/Rinst33659443778b2d/SpaceTrooper/extdata/CosMx_DBKero_Tiny/DBKero-polygons.csv"

2.1.1 CosMx protein

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

2.2 Route B — Via SpatialExperimentIO, then standardize

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

3 Xenium (10x Genomics)

A small Xenium example is also included for demonstration.

xepath <- system.file("extdata", "Xenium_small", package = "SpaceTrooper")
xepath
#> [1] "/tmp/RtmpBasoI1/Rinst33659443778b2d/SpaceTrooper/extdata/Xenium_small"

3.1 Route A — Direct loading with SpaceTrooper

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/RtmpBasoI1/Rinst33659443778b2d/SpaceTrooper/extdata/Xenium_small/cell_boundaries.parquet"
metadata(spe_xen_a)$technology
#> [1] "10X_Xenium"

3.2 Route B — Via SpatialExperimentIO, then standardize

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

4 File layout notes

  • CosMx: per-cell metadata (e.g. 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.
  • Xenium: unzipped bundle contains an 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.

5 Session info

sessionInfo()
#> R version 4.5.1 Patched (2025-08-23 r88802)
#> Platform: x86_64-pc-linux-gnu
#> Running under: Ubuntu 24.04.3 LTS
#> 
#> Matrix products: default
#> BLAS:   /home/biocbuild/bbs-3.22-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] SpaceTrooper_1.0.0          SpatialExperiment_1.20.0   
#>  [3] SingleCellExperiment_1.32.0 SummarizedExperiment_1.40.0
#>  [5] Biobase_2.70.0              GenomicRanges_1.62.0       
#>  [7] Seqinfo_1.0.0               IRanges_2.44.0             
#>  [9] S4Vectors_0.48.0            BiocGenerics_0.56.0        
#> [11] generics_0.1.4              MatrixGenerics_1.22.0      
#> [13] matrixStats_1.5.0           BiocStyle_2.38.0           
#> 
#> loaded via a namespace (and not attached):
#>   [1] DBI_1.2.3                 gridExtra_2.3            
#>   [3] rlang_1.1.6               magrittr_2.0.4           
#>   [5] scater_1.38.0             e1071_1.7-16             
#>   [7] compiler_4.5.1            DelayedMatrixStats_1.32.0
#>   [9] SpatialExperimentIO_1.2.0 sfheaders_0.4.4          
#>  [11] vctrs_0.6.5               shape_1.4.6.1            
#>  [13] pkgconfig_2.0.3           fastmap_1.2.0            
#>  [15] backports_1.5.0           magick_2.9.0             
#>  [17] XVector_0.50.0            labeling_0.4.3           
#>  [19] scuttle_1.20.0            rmarkdown_2.30           
#>  [21] ggbeeswarm_0.7.2          tinytex_0.57             
#>  [23] purrr_1.1.0               bit_4.6.0                
#>  [25] glmnet_4.1-10             xfun_0.53                
#>  [27] cachem_1.1.0              beachmat_2.26.0          
#>  [29] jsonlite_2.0.0            rhdf5filters_1.22.0      
#>  [31] DelayedArray_0.36.0       Rhdf5lib_1.32.0          
#>  [33] BiocParallel_1.44.0       irlba_2.3.5.1            
#>  [35] broom_1.0.10              parallel_4.5.1           
#>  [37] R6_2.6.1                  bslib_0.9.0              
#>  [39] RColorBrewer_1.1-3        limma_3.66.0             
#>  [41] car_3.1-3                 jquerylib_0.1.4          
#>  [43] iterators_1.0.14          Rcpp_1.1.0               
#>  [45] bookdown_0.45             assertthat_0.2.1         
#>  [47] knitr_1.50                R.utils_2.13.0           
#>  [49] splines_4.5.1             Matrix_1.7-4             
#>  [51] tidyselect_1.2.1          viridis_0.6.5            
#>  [53] dichromat_2.0-0.1         abind_1.4-8              
#>  [55] yaml_2.3.10               codetools_0.2-20         
#>  [57] lattice_0.22-7            tibble_3.3.0             
#>  [59] withr_3.0.2               S7_0.2.0                 
#>  [61] evaluate_1.0.5            sf_1.0-21                
#>  [63] survival_3.8-3            units_1.0-0              
#>  [65] proxy_0.4-27              pillar_1.11.1            
#>  [67] BiocManager_1.30.26       ggpubr_0.6.2             
#>  [69] carData_3.0-5             KernSmooth_2.23-26       
#>  [71] foreach_1.5.2             ggplot2_4.0.0            
#>  [73] sparseMatrixStats_1.22.0  scales_1.4.0             
#>  [75] class_7.3-23              glue_1.8.0               
#>  [77] tools_4.5.1               BiocNeighbors_2.4.0      
#>  [79] robustbase_0.99-6         data.table_1.17.8        
#>  [81] ScaledMatrix_1.18.0       locfit_1.5-9.12          
#>  [83] ggsignif_0.6.4            cowplot_1.2.0            
#>  [85] rhdf5_2.54.0              grid_4.5.1               
#>  [87] tidyr_1.3.1               DropletUtils_1.30.0      
#>  [89] edgeR_4.8.0               beeswarm_0.4.0           
#>  [91] BiocSingular_1.26.0       HDF5Array_1.38.0         
#>  [93] vipor_0.4.7               rsvd_1.0.5               
#>  [95] Formula_1.2-5             cli_3.6.5                
#>  [97] viridisLite_0.4.2         S4Arrays_1.10.0          
#>  [99] arrow_22.0.0              dplyr_1.1.4              
#> [101] DEoptimR_1.1-4            gtable_0.3.6             
#> [103] R.methodsS3_1.8.2         rstatix_0.7.3            
#> [105] sass_0.4.10               digest_0.6.37            
#> [107] classInt_0.4-11           ggrepel_0.9.6            
#> [109] SparseArray_1.10.0        dqrng_0.4.1              
#> [111] rjson_0.2.23              farver_2.1.2             
#> [113] htmltools_0.5.8.1         R.oo_1.27.1              
#> [115] lifecycle_1.0.4           h5mread_1.2.0            
#> [117] statmod_1.5.1             bit64_4.6.0-1