Contents

Introduction

The VisiumIO package provides a set of functions to import 10X Genomics Visium experiment data into a SpatialExperiment object. The package makes use of the SpatialExperiment data structure, which provides a set of classes and methods to handle spatially resolved transcriptomics data.

TENxIO Supported Formats

Extension Class Imported as
.h5 TENxH5 SingleCellExperiment w/ TENxMatrix
.mtx / .mtx.gz TENxMTX SummarizedExperiment w/ dgCMatrix
.tar.gz TENxFileList SingleCellExperiment w/ dgCMatrix
peak_annotation.tsv TENxPeaks GRanges
fragments.tsv.gz TENxFragments RaggedExperiment
.tsv / .tsv.gz TENxTSV tibble

VisiumIO Supported Formats

Extension Class Imported as
spatial.tar.gz TENxSpatialList inter. DataFrame list

Installation

if (!require("BiocManager", quietly = TRUE))
    install.packages("BiocManager")
BiocManager::install("VisiumIO")

Loading package

library(VisiumIO)

TENxVisium

The TENxVisium class is used to import a single sample of 10X Visium data. The TENxVisium constructor function takes the following arguments:

TENxVisium(
    resources = "path/to/10x/visium/file.tar.gz",
    spatialResource = "path/to/10x/visium/spatial/file.spatial.tar.gz",
    spacerangerOut = "path/to/10x/visium/sample/folder",
    sample_id = "sample01",
    images = c("lowres", "hires", "detected", "aligned"),
    jsonFile = "scalefactors_json.json",
    tissuePattern = "tissue_positions.*\\.csv",
    spatialCoordsNames = c("pxl_col_in_fullres", "pxl_row_in_fullres")
)

The resource argument is the path to the 10X Visium file. The spatialResource argument is the path to the 10X Visium spatial file. It usually ends in spatial.tar.gz.

Example from SpatialExperiment

Note that we use the image = "lowres" and processing = "raw" arguments based on the name of the tissue_*_image.png file and *_feature_bc_matrix folder in the spaceranger output. The directory structure for a single sample is shown below:

    section1
    └── outs
        ├── spatial
        │   ├── tissue_lowres_image.png
        │   └── tissue_positions_list.csv
        └── raw_feature_bc_matrix
            ├── barcodes.tsv
            ├── features.tsv
            └── matrix.mtx

Creating a TENxVisium instance

Using the example data in SpatialExperiment, we can load the section1 sample using TENxVisium.

sample_dir <- system.file(
    file.path("extdata", "10xVisium", "section1"),
    package = "SpatialExperiment"
)

vis <- TENxVisium(
    spacerangerOut = sample_dir, processing = "raw", images = "lowres"
)
vis
## An object of class "TENxVisium"
## Slot "resources":
## TENxFileList of length 3
## names(3): barcodes.tsv features.tsv matrix.mtx
## 
## Slot "spatialList":
## TENxSpatialList of length 3
## names(3): scalefactors_json.json tissue_lowres_image.png tissue_positions_list.csv
## 
## Slot "coordNames":
## [1] "pxl_col_in_fullres" "pxl_row_in_fullres"
## 
## Slot "sampleId":
## [1] "sample01"

The show method of the TENxVisium class displays the object’s metadata.

Importing into SpatialExperiment

The TEnxVisium object can be imported into a SpatialExperiment object using the import function.

import(vis)
## class: SpatialExperiment 
## dim: 50 50 
## metadata(0):
## assays(1): counts
## rownames: NULL
## rowData names(1): Symbol
## colnames(50): AAACAACGAATAGTTC-1 AAACAAGTATCTCCCA-1 ...
##   AAAGTCGACCCTCAGT-1 AAAGTGCCATCAATTA-1
## colData names(4): in_tissue array_row array_col sample_id
## reducedDimNames(0):
## mainExpName: NULL
## altExpNames(0):
## spatialCoords names(2) : pxl_col_in_fullres pxl_row_in_fullres
## imgData names(4): sample_id image_id data scaleFactor

TENxVisiumList

The TENxVisiumList class is used to import multiple samples of 10X Visium. The interface is a bit more simple in that you only need to provide the space ranger output folder as input to the function.

TENxVisiumList(
    sampleFolders = "path/to/10x/visium/sample/folder",
    sample_ids = c("sample01", "sample02"),
    ...
)

The sampleFolders argument is a character vector of paths to the spaceranger output folder. Note that each folder must contain an outs directory. The sample_ids argument is a character vector of sample ids.

Example from SpatialExperiment

The directory structure for multiple samples (section1 and section2) is shown below:

    section1
    └── outs
    |   ├── spatial
    |   └── raw_feature_bc_matrix
    section2
    └── outs
        ├── spatial
        └── raw_feature_bc_matrix

Creating a TENxVisiumList

The main inputs to TENxVisiumList are the sampleFolders and sample_ids. These correspond to the spaceranger output sample folders and a vector of sample identifiers, respectively.

sample_dirs <- list.dirs(
    system.file(
        file.path("extdata", "10xVisium"), package = "SpatialExperiment"
    ),
    recursive = FALSE, full.names = TRUE
)
    
vlist <- TENxVisiumList(
    sampleFolders = sample_dirs,
    sample_ids = basename(sample_dirs),
    processing = "raw",
    image = "lowres"
)
vlist
## An object of class "TENxVisiumList"
## Slot "VisiumList":
## List of length 2

Importing into SpatialExperiment

The import method combines both SingleCellExperiment objects along with the spatial information into a single SpatialExperiment object. The number of columns in the SpatialExperiment object is equal to the number of cells across both samples (section1 and section2).

import(vlist)
## class: SpatialExperiment 
## dim: 50 99 
## metadata(0):
## assays(1): counts
## rownames: NULL
## rowData names(1): Symbol
## colnames(99): AAACAACGAATAGTTC-1 AAACAAGTATCTCCCA-1 ...
##   AAAGTCGACCCTCAGT-1 AAAGTGCCATCAATTA-1
## colData names(4): in_tissue array_row array_col sample_id
## reducedDimNames(0):
## mainExpName: NULL
## altExpNames(0):
## spatialCoords names(2) : pxl_col_in_fullres pxl_row_in_fullres
## imgData names(4): sample_id image_id data scaleFactor

Session Info

sessionInfo()
## R version 4.4.0 RC (2024-04-16 r86468)
## Platform: x86_64-pc-linux-gnu
## Running under: Ubuntu 22.04.4 LTS
## 
## Matrix products: default
## BLAS:   /home/biocbuild/bbs-3.20-bioc/R/lib/libRblas.so 
## LAPACK: /usr/lib/x86_64-linux-gnu/lapack/liblapack.so.3.10.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] VisiumIO_1.1.0              TENxIO_1.7.0               
##  [3] SingleCellExperiment_1.27.0 SummarizedExperiment_1.35.0
##  [5] Biobase_2.65.0              GenomicRanges_1.57.0       
##  [7] GenomeInfoDb_1.41.0         IRanges_2.39.0             
##  [9] S4Vectors_0.43.0            BiocGenerics_0.51.0        
## [11] MatrixGenerics_1.17.0       matrixStats_1.3.0          
## [13] BiocStyle_2.33.0           
## 
## loaded via a namespace (and not attached):
##  [1] rjson_0.2.21             xfun_0.43                bslib_0.7.0             
##  [4] lattice_0.22-6           tzdb_0.4.0               vctrs_0.6.5             
##  [7] tools_4.4.0              parallel_4.4.0           tibble_3.2.1            
## [10] fansi_1.0.6              pkgconfig_2.0.3          BiocBaseUtils_1.7.0     
## [13] Matrix_1.7-0             lifecycle_1.0.4          GenomeInfoDbData_1.2.12 
## [16] compiler_4.4.0           codetools_0.2-20         htmltools_0.5.8.1       
## [19] sass_0.4.9               yaml_2.3.8               pillar_1.9.0            
## [22] crayon_1.5.2             jquerylib_0.1.4          DelayedArray_0.31.0     
## [25] cachem_1.0.8             magick_2.8.3             abind_1.4-5             
## [28] tidyselect_1.2.1         digest_0.6.35            bookdown_0.39           
## [31] fastmap_1.1.1            grid_4.4.0               archive_1.1.8           
## [34] cli_3.6.2                SparseArray_1.5.0        magrittr_2.0.3          
## [37] S4Arrays_1.5.0           utf8_1.2.4               readr_2.1.5             
## [40] UCSC.utils_1.1.0         bit64_4.0.5              rmarkdown_2.26          
## [43] XVector_0.45.0           httr_1.4.7               bit_4.0.5               
## [46] hms_1.1.3                SpatialExperiment_1.15.0 evaluate_0.23           
## [49] knitr_1.46               BiocIO_1.15.0            rlang_1.1.3             
## [52] Rcpp_1.0.12              glue_1.7.0               BiocManager_1.30.22     
## [55] vroom_1.6.5              jsonlite_1.8.8           R6_2.5.1                
## [58] zlibbioc_1.51.0