Skip to contents

Creates a SpatialExperiment from the Space Ranger output directories for 10x Genomics Visium spatial gene expression data.

Usage

read10xVisium(
  samples = "",
  sample_id = paste0("sample", sprintf("%02d", seq_along(samples))),
  type = c("HDF5", "sparse"),
  data = c("filtered", "raw"),
  images = "lowres",
  load = TRUE
)

Arguments

samples

a character vector specifying one or more directories, each corresponding to a 10x Genomics Visium sample (see Details); if provided, names will be used as sample identifiers

sample_id

character string specifying unique sample identifiers, one for each directory specified via samples; ignored if !is.null(names(samples))

type

character string specifying the type of format to read count data from (see read10xCounts)

data

character string specifying whether to read in filtered (spots mapped to tissue) or raw data (all spots).

images

character vector specifying which images to include. Valid values are "lowres", "hires", "fullres", "detected", "aligned"

load

logical; should the image(s) be loaded into memory as a grob? If FALSE, will store the path/URL instead.

Value

a SpatialExperiment object

Details

The constructor assumes data from each sample are located in a single output directory as returned by Space Ranger, thus having the following file organization (where "raw/filtered" refers to either "raw" or "filtered" to match the `data` argument.) The base directory "outs/" from Space Ranger can either be included manually in the paths provided in `samples`, or can be ignored; if ignored, it will be added automatically. The `.h5` files are used if `type = "HDF5"`. (Note that `tissue_positions.csv` was renamed in Space Ranger v2.0.0.)

sample
· | — outs
· · | — raw/filtered_feature_bc_matrix.h5
· · | — raw/filtered_feature_bc_matrix
· · · | — barcodes.tsv.gz
· · · | — features.tsv.gz
· · · | — matrix.mtx.gz
· · | — spatial
· · · | — scalefactors_json.json
· · · | — tissue_lowres_image.png
· · · | — tissue_positions.csv

Author

Helena L. Crowell

Examples

dir <- system.file(
  file.path("extdata", "10xVisium"), 
  package = "SpatialExperiment")
  
sample_ids <- c("section1", "section2")
samples <- file.path(dir, sample_ids, "outs")
  
list.files(samples[1])
#> [1] "raw_feature_bc_matrix" "spatial"              
list.files(file.path(samples[1], "spatial"))
#> [1] "scalefactors_json.json"    "tissue_lowres_image.png"  
#> [3] "tissue_positions_list.csv"
file.path(samples[1], "raw_feature_bc_matrix")
#> [1] "/tmp/RtmppdZi9A/temp_libpathcf263c794d1/SpatialExperiment/extdata/10xVisium/section1/outs/raw_feature_bc_matrix"

(spe <- read10xVisium(samples, sample_ids, 
  type = "sparse", data = "raw", 
  images = "lowres", load = FALSE))
#> class: SpatialExperiment 
#> dim: 50 99 
#> metadata(0):
#> assays(1): counts
#> rownames(50): ENSMUSG00000051951 ENSMUSG00000089699 ...
#>   ENSMUSG00000005886 ENSMUSG00000101476
#> 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

# base directory 'outs/' from Space Ranger can also be omitted
samples2 <- file.path(dir, sample_ids)
(spe2 <- read10xVisium(samples2, sample_ids, 
  type = "sparse", data = "raw", 
  images = "lowres", load = FALSE))
#> class: SpatialExperiment 
#> dim: 50 99 
#> metadata(0):
#> assays(1): counts
#> rownames(50): ENSMUSG00000051951 ENSMUSG00000089699 ...
#>   ENSMUSG00000005886 ENSMUSG00000101476
#> 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

# tabulate number of spots mapped to tissue
cd <- colData(spe)
table(
  in_tissue = cd$in_tissue, 
  sample_id = cd$sample_id)
#>          sample_id
#> in_tissue section1 section2
#>     FALSE       28       27
#>     TRUE        22       22

# view available images
imgData(spe)
#> DataFrame with 2 rows and 4 columns
#>     sample_id    image_id   data scaleFactor
#>   <character> <character> <list>   <numeric>
#> 1    section1      lowres   ####   0.0510334
#> 2    section2      lowres   ####   0.0510334