Skip to contents

The SpatialExperiment class provides methods for getting or setting named assays. For example, molecules(spe) will get or set an assay named molecules from object spe, equivalent to assay(spe, i = "molecules"). This provides a convenient interface for users and encourages standardization of assay names across packages.

Available methods

In the following code, spe is a SpatialExperiment object, value is a BumpyMatrix-like object with the same dimensions as spe, and ... are further arguments passed to assay (for the getter) or assay<- (for the setter).

molecules(x, ...), molecules(x, ...) <- value:

Get or set an assay named molecules, which is usually assumed to be a BumpyMatrix-formatted object containing spatial coordinates (and any other information) of the individual molecules per gene per cell.

See also

assay and assay<-

Author

Dario Righelli

Examples

example(SpatialExperiment)
#> 
#> SptlEx> #########################################################
#> SptlEx> # Example 1: Spot-based ST (10x Visium) using constructor
#> SptlEx> #########################################################
#> SptlEx> 
#> SptlEx> dir <- system.file(
#> SptlEx+     file.path("extdata", "10xVisium", "section1", "outs"),
#> SptlEx+     package = "SpatialExperiment")
#> 
#> SptlEx> # read in counts
#> SptlEx> fnm <- file.path(dir, "raw_feature_bc_matrix")
#> 
#> SptlEx> sce <- DropletUtils::read10xCounts(fnm)
#> 
#> SptlEx> # read in image data
#> SptlEx> img <- readImgData(
#> SptlEx+     path = file.path(dir, "spatial"),
#> SptlEx+     sample_id="foo")
#> 
#> SptlEx> # read in spatial coordinates
#> SptlEx> fnm <- file.path(dir, "spatial", "tissue_positions_list.csv")
#> 
#> SptlEx> xyz <- read.csv(fnm, header = FALSE,
#> SptlEx+     col.names = c(
#> SptlEx+         "barcode", "in_tissue", "array_row", "array_col",
#> SptlEx+         "pxl_row_in_fullres", "pxl_col_in_fullres"))
#> 
#> SptlEx> # construct observation & feature metadata
#> SptlEx> rd <- S4Vectors::DataFrame(
#> SptlEx+     symbol = rowData(sce)$Symbol)
#> 
#> SptlEx> # construct 'SpatialExperiment'
#> SptlEx> (spe <- SpatialExperiment(
#> SptlEx+     assays = list(counts = assay(sce)),
#> SptlEx+     rowData = rd,
#> SptlEx+     colData = DataFrame(xyz),
#> SptlEx+     spatialCoordsNames = c("pxl_col_in_fullres", "pxl_row_in_fullres"),
#> SptlEx+     imgData = img,
#> SptlEx+     sample_id = "foo"))
#> class: SpatialExperiment 
#> dim: 50 50 
#> metadata(0):
#> assays(1): counts
#> rownames(50): ENSMUSG00000051951 ENSMUSG00000089699 ...
#>   ENSMUSG00000005886 ENSMUSG00000101476
#> rowData names(1): symbol
#> colnames: NULL
#> colData names(5): barcode 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
#> 
#> SptlEx> #############################################################
#> SptlEx> # Example 2: Spot-based ST (10x Visium) using 'read10xVisium'
#> SptlEx> #############################################################
#> SptlEx> 
#> SptlEx> # see ?read10xVisium for details
#> SptlEx> example(read10xVisium)
#> 
#> rd10xV> dir <- system.file(
#> rd10xV+   file.path("extdata", "10xVisium"), 
#> rd10xV+   package = "SpatialExperiment")
#> 
#> rd10xV> sample_ids <- c("section1", "section2")
#> 
#> rd10xV> samples <- file.path(dir, sample_ids, "outs")
#> 
#> rd10xV> list.files(samples[1])
#> [1] "raw_feature_bc_matrix" "spatial"              
#> 
#> rd10xV> list.files(file.path(samples[1], "spatial"))
#> [1] "scalefactors_json.json"    "tissue_lowres_image.png"  
#> [3] "tissue_positions_list.csv"
#> 
#> rd10xV> file.path(samples[1], "raw_feature_bc_matrix")
#> [1] "/tmp/RtmppdZi9A/temp_libpathcf263c794d1/SpatialExperiment/extdata/10xVisium/section1/outs/raw_feature_bc_matrix"
#> 
#> rd10xV> (spe <- read10xVisium(samples, sample_ids, 
#> rd10xV+   type = "sparse", data = "raw", 
#> rd10xV+   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
#> 
#> rd10xV> # base directory 'outs/' from Space Ranger can also be omitted
#> rd10xV> samples2 <- file.path(dir, sample_ids)
#> 
#> rd10xV> (spe2 <- read10xVisium(samples2, sample_ids, 
#> rd10xV+   type = "sparse", data = "raw", 
#> rd10xV+   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
#> 
#> rd10xV> # tabulate number of spots mapped to tissue
#> rd10xV> cd <- colData(spe)
#> 
#> rd10xV> table(
#> rd10xV+   in_tissue = cd$in_tissue, 
#> rd10xV+   sample_id = cd$sample_id)
#>          sample_id
#> in_tissue section1 section2
#>     FALSE       28       27
#>     TRUE        22       22
#> 
#> rd10xV> # view available images
#> rd10xV> 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
#> 
#> SptlEx> ##############################
#> SptlEx> # Example 3: Molecule-based ST
#> SptlEx> ##############################
#> SptlEx> 
#> SptlEx> # create simulated data
#> SptlEx> n <- 1000; ng <- 50; nc <- 20
#> 
#> SptlEx> # sample xy-coordinates in [0,1]
#> SptlEx> x <- runif(n)
#> 
#> SptlEx> y <- runif(n)
#> 
#> SptlEx> # assign each molecule to some gene-cell pair
#> SptlEx> gs <- paste0("gene", seq(ng))
#> 
#> SptlEx> cs <- paste0("cell", seq(nc))
#> 
#> SptlEx> gene <- sample(gs, n, TRUE)
#> 
#> SptlEx> cell <- sample(cs, n, TRUE)
#> 
#> SptlEx> # construct data.frame of molecule coordinates
#> SptlEx> df <- data.frame(gene, cell, x, y)
#> 
#> SptlEx> # (assure gene & cell are factor so that
#> SptlEx> # missing observations aren't dropped)
#> SptlEx> df$gene <- factor(df$gene, gs)
#> 
#> SptlEx> df$cell <- factor(df$cell, cs)
#> 
#> SptlEx> # construct BumpyMatrix
#> SptlEx> mol <- BumpyMatrix::splitAsBumpyMatrix(
#> SptlEx+     df[, c("x", "y")],
#> SptlEx+     row = df$gene, column = df$cell)
#> 
#> SptlEx> # get count matrix
#> SptlEx> y <- with(df, table(gene, cell))
#> 
#> SptlEx> y <- as.matrix(unclass(y))
#> 
#> SptlEx> # construct SpatialExperiment
#> SptlEx> (spe_mol <- SpatialExperiment(
#> SptlEx+     assays = list(
#> SptlEx+         counts = y, 
#> SptlEx+         molecules = mol)))
#> class: SpatialExperiment 
#> dim: 50 20 
#> metadata(0):
#> assays(2): counts molecules
#> rownames(50): gene1 gene2 ... gene49 gene50
#> rowData names(0):
#> colnames(20): cell1 cell2 ... cell19 cell20
#> colData names(1): sample_id
#> reducedDimNames(0):
#> mainExpName: NULL
#> altExpNames(0):
#> spatialCoords names(0) :
#> imgData names(0):
molecules(spe_mol)
#> 50 x 20 BumpyDataFrameMatrix
#> rownames: gene1 gene2 ... gene49 gene50 
#> colnames: cell1 cell2 ... cell19 cell20 
#> preview [1,1]:
#>   DataFrame with 1 row and 2 columns
#>             x         y
#>     <numeric> <numeric>
#>   1  0.990712  0.517485