The SpatialExperiment
class is designed to represent spatially
resolved transcriptomics (ST) data. It inherits from the
SingleCellExperiment
class and is used in the same manner.
In addition, the class supports storage of spatial information
via spatialData
and spatialCoords
, and storage of
images via imgData
.
For demonstration of the general class structure, we load an example
SpatialExperiment
(abbreviated as SPE) object (variable spe
):
library(SpatialExperiment)
## Loading required package: SingleCellExperiment
## Loading required package: SummarizedExperiment
## Loading required package: MatrixGenerics
## Loading required package: matrixStats
##
## Attaching package: 'MatrixGenerics'
## The following objects are masked from 'package:matrixStats':
##
## colAlls, colAnyNAs, colAnys, colAvgsPerRowSet, colCollapse,
## colCounts, colCummaxs, colCummins, colCumprods, colCumsums,
## colDiffs, colIQRDiffs, colIQRs, colLogSumExps, colMadDiffs,
## colMads, colMaxs, colMeans2, colMedians, colMins, colOrderStats,
## colProds, colQuantiles, colRanges, colRanks, colSdDiffs, colSds,
## colSums2, colTabulates, colVarDiffs, colVars, colWeightedMads,
## colWeightedMeans, colWeightedMedians, colWeightedSds,
## colWeightedVars, rowAlls, rowAnyNAs, rowAnys, rowAvgsPerColSet,
## rowCollapse, rowCounts, rowCummaxs, rowCummins, rowCumprods,
## rowCumsums, rowDiffs, rowIQRDiffs, rowIQRs, rowLogSumExps,
## rowMadDiffs, rowMads, rowMaxs, rowMeans2, rowMedians, rowMins,
## rowOrderStats, rowProds, rowQuantiles, rowRanges, rowRanks,
## rowSdDiffs, rowSds, rowSums2, rowTabulates, rowVarDiffs, rowVars,
## rowWeightedMads, rowWeightedMeans, rowWeightedMedians,
## rowWeightedSds, rowWeightedVars
## Loading required package: GenomicRanges
## Loading required package: stats4
## Loading required package: BiocGenerics
##
## Attaching package: 'BiocGenerics'
## The following objects are masked from 'package:stats':
##
## IQR, mad, sd, var, xtabs
## The following objects are masked from 'package:base':
##
## Filter, Find, Map, Position, Reduce, anyDuplicated, append,
## as.data.frame, basename, cbind, colnames, dirname, do.call,
## duplicated, eval, evalq, get, grep, grepl, intersect, is.unsorted,
## lapply, mapply, match, mget, order, paste, pmax, pmax.int, pmin,
## pmin.int, rank, rbind, rownames, sapply, setdiff, sort, table,
## tapply, union, unique, unsplit, which.max, which.min
## Loading required package: S4Vectors
##
## Attaching package: 'S4Vectors'
## The following objects are masked from 'package:base':
##
## I, expand.grid, unname
## Loading required package: IRanges
## Loading required package: GenomeInfoDb
## Loading required package: Biobase
## Welcome to Bioconductor
##
## Vignettes contain introductory material; view with
## 'browseVignettes()'. To cite Bioconductor, see
## 'citation("Biobase")', and for packages 'citation("pkgname")'.
##
## Attaching package: 'Biobase'
## The following object is masked from 'package:MatrixGenerics':
##
## rowMedians
## The following objects are masked from 'package:matrixStats':
##
## anyMissing, rowMedians
example(read10xVisium, echo = FALSE)
spe
## 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(1): sample_id
## reducedDimNames(0):
## mainExpName: NULL
## altExpNames(0):
## spatialData names(3) : in_tissue array_row array_col
## spatialCoords names(2) : pxl_col_in_fullres pxl_row_in_fullres
## imgData names(4): sample_id image_id data scaleFactor
spatialData
& -Coords
In addition to observation metadata stored inside the colData
slot of
a SingleCellExperiment
, the SpatialExperiment
class can accommodate:
spatialData
, a DataFrame
containing spatial metadataspatialCoords
, a numeric matrix of spatial coordinates (e.g. x
and y
)Both spatialData
and spatialCoords
are stored
separately inside the int_colData
slot.
Note that the colData
, spatialData
, and spatialCoords
slots follow
a hierarchical structure where colData
> spatialData
> spatialCoords
.
Here, each accessor function allows joint accession of the target slot, and (optionally)
any slot(s) that precedes it.
Specifically, the following commands are supported that may be used to access specific subsets of (spatial) metadata associated with each column (observation, e.g. spots or cells) in a SPE:
spatialCoords(spe)
spatialData(spe)
spatialData(spe, spatialCoords = TRUE)
colData(spe, spatialData = TRUE)
colData(spe, spatialCoords = TRUE)
colData(spe, spatialData = TRUE, spatialCoords = TRUE)
imgData
All image related data are stored inside the int_metadata
’s
imgData
field as a DataFrame
of the following structure:
sample_id
the image belongs toimage_id
in order to accommodate multiple imagesdata
(a SpatialImage
object)scaleFactor
that adjusts pixel positions of the original,The imgData()
accessor can be used to retrieve
the image data stored within the object:
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
SpatialImage
classImages are stored inside the data
field of the imgData
as a list of
SpatialImage
s. Each image may be of one of the following sub-classes:
LoadedSpatialImage
raster
object@image
contains a raster
object: a matrix
of RGB colors for each pixel in the imageStoredSpatialImage
@path
specifies a local file from which to retrieve the imageRemoteSpatialImage
@url
specifies where to retrieve the image fromA SpatialImage
can be accessed using getImg()
,
or retrieved directly from the imgData()
:
(spi <- getImg(spe))
## A 600 x 576 StoredSpatialImage
## imgSource():
## /tmp/RtmpGjXjlV/Rinst21541f790c1ff1/Spat
## ialExperiment/extdata/10xVisium/section1
## /spatial/tissue_lowres_image.png
identical(spi, imgData(spe)$data[[1]])
## [1] TRUE
Data available in an object of class SpatialImage
may be
accessed via the imgRaster()
and imgSource()
accessors:
plot(imgRaster(spe))
Images entries may be added or removed from a SpatialExperiment
’s
imgData
DataFrame
using addImg()
and rmvImg()
, respectively.
Besides a path or URL to source the image from and a numeric scale factor,
addImg()
requires specification of the sample_id
the new image belongs to,
and an image_id
that is not yet in use for that sample:
url <- "https://i.redd.it/3pw5uah7xo041.jpg"
spe <- addImg(spe,
sample_id = "section1",
image_id = "pomeranian",
imageSource = url,
scaleFactor = NA_real_,
load = TRUE)
img <- imgRaster(spe,
sample_id = "section1",
image_id = "pomeranian")
plot(img)
The rmvImg()
function is more flexible in the specification
of the sample_id
and image_id
arguments. Specifically:
TRUE
is equivalent to all, e.g.sample_id = "<sample>"
, image_id = TRUE
NULL
defaults to the first entry available, e.g.sample_id = "<sample>"
, image_id = NULL
For example, sample_id = TRUE
, image_id = TRUE
will specify all images;
sample_id = NULL
, image_id = NULL
corresponds to the first image entry in the imgData
;
sample_id = TRUE
, image_id = NULL
equals the first image for all samples; and
sample_id = NULL
, image_id = TRUE
matches all images for the first sample.
Here, we remove section1
’s pomeranian
image added in the previous
code chunk; the image is now completely gone from the imgData
:
imgData(spe <- rmvImg(spe, "section1", "pomeranian"))
## 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
The SpatialExperiment
constructor provides several arguments
to give maximum flexibility to the user.
In particular, these include:
spatialData
, a DataFrame
with all the data
associated with the spatial information (optionally
including spatial coordinates to use as spatialCoords
)spatialCoords
, a numeric matrix
containing spatial coordinatesspatialDataNames
, a character vector specifyingcolData
fields correspond to spatial metadataspatialCoordsNames
, a character vector specifying whichcolData
or spatialData
fields correspond to spatial coordinatesBoth spatialData
and SpatialCoords
can be supplied directly via colData
by specifying the column names that correspond to metadata and spatial coordinates,
respectively, via spatialDataNames
and spatialCoordsNames
:
n <- length(z <- letters)
y <- matrix(nrow = n, ncol = n)
cd <- DataFrame(x = seq(n), y = seq(n), z)
spe1 <- SpatialExperiment(
assay = y,
colData = cd,
spatialDataNames = "z",
spatialCoordsNames = c("x", "y"))
Alternatively, spatialData
and spatialCoords
may be supplied separately
directly as a DataFrame
and matrix
, e.g.:
xy <- as.matrix(cd[, c("x", "y")])
spe2 <- SpatialExperiment(
assay = y,
spatialData = cd["z"],
spatialCoords = xy)
Or, one of spatialData
or spatialCoords
can be supplied, while the other may be
extracted from the input colData
according to spatialData
or spatialCoordsNames
:
spe3 <- SpatialExperiment(
assay = y,
colData = cd[-3],
spatialData = cd["z"],
spatialCoordsNames = c("x", "y"))
spe4 <- SpatialExperiment(
assay = y,
colData = cd[-c(1, 2)],
spatialCoords = xy,
spatialDataNames = "z")
Importantly, all of the above SpatialExperiment()
function calls
lead to construction of the exact same object:
all(identical(spe1, spe2),
identical(spe1, spe3),
identical(spe1, spe4),
identical(spe2, spe3),
identical(spe2, spe4),
identical(spe3, spe4))
## [1] TRUE
Finally, all of spatialData/Coords(Names)
are optional. I.e.,
we can construct a SPE using only a subset of the above arguments:
spe <- SpatialExperiment(
assays = y,
spatialCoords = xy)
isEmpty(spatialData(spe))
## [1] TRUE
In general, spatialData/CoordsNames
take precedence over spatialData/Coords
,
i.e., if both are supplied, the latter will be ignored. In other words,
spatialData/Coords
are preferentially extracted from the DataFrame
s
provided via spatial/colData
. E.g., in the following function call,
spatialCoords
will be ignored, and xy-coordinates are instead extracted
from the input colData
according to the specified spatialCoordsNames
:
n <- 10; m <- 20
y <- matrix(nrow = n, ncol = m)
cd <- DataFrame(x = seq(m), y = seq(m))
xy <- matrix(nrow = m, ncol = 2)
colnames(xy) <- c("x", "y")
SpatialExperiment(
assay = y,
colData = cd,
spatialCoordsNames = c("x", "y"),
spatialCoords = xy)
## both 'spatialCoords' and 'spatialCoordsNames' have been supplied; using 'spatialCoords'. Set either to NULL to suppress this message
When working with spot-based ST data, such as 10x Genomics Visium or other
platforms providing images, it is possible to store the image information
in the dedicated imgData
structure.
Also, the SpatialExperiment
class stores a sample_id
value in the
spatialData
structure, which is possible to set with the sample_id
argument (default is “sample_01”).
Here we show how to load the default Space Ranger data files from a
10x Genomics Visium experiment, and build a SpatialExperiment
object.
In particular, the readImgData()
function is used to build an imgData
DataFrame
to be passed to the SpatialExperiment
constructor.
The sample_id
used to build the imgData
object must be the same one
used to build the SpatialExperiment
object, otherwise an error is returned.
dir <- system.file(
file.path("extdata", "10xVisium", "section1"),
package = "SpatialExperiment")
# read in counts
fnm <- file.path(dir, "raw_feature_bc_matrix")
sce <- DropletUtils::read10xCounts(fnm)
# read in image data
img <- readImgData(
path = file.path(dir, "spatial"),
sample_id = "foo")
# read in spatial coordinates
fnm <- file.path(dir, "spatial", "tissue_positions_list.csv")
xyz <- read.csv(fnm, header = FALSE,
col.names = c(
"barcode", "in_tissue", "array_row", "array_col",
"pxl_row_in_fullres", "pxl_col_in_fullres"))
# construct observation & feature metadata
rd <- S4Vectors::DataFrame(
symbol = rowData(sce)$Symbol)
# construct 'SpatialExperiment'
(spe <- SpatialExperiment(
assays = list(counts = assay(sce)),
rowData = rd,
colData = colData(sce),
imgData = img,
spatialData = DataFrame(xyz),
spatialCoordsNames = c("pxl_col_in_fullres", "pxl_row_in_fullres"),
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(3): Sample Barcode sample_id
## reducedDimNames(0):
## mainExpName: NULL
## altExpNames(0):
## spatialData names(4) : barcode in_tissue array_row array_col
## spatialCoords names(2) : pxl_col_in_fullres pxl_row_in_fullres
## imgData names(4): sample_id image_id data scaleFactor
Alternatively, the read10xVisium()
function facilitates the import of
10x Genomics Visium data to handle one or more samples organized in
folders reflecting the default Space Ranger folder tree organization:
sample
. |—outs
· · |—raw/filtered_feature_bc_matrix.h5
· · |—raw/filtered_feature_bc_matrix
· · · · |—barcodes.tsv
· · · · |—features.tsv
· · · · |—matrix.mtx
· · |—spatial
· · · · |—scalefactors_json.json
· · · · |—tissue_lowres_image.png
· · · · |—tissue_positions_list.csv
Using read10xVisium()
, the above code to construct the same SPE is reduced to:
dir <- system.file(
file.path("extdata", "10xVisium"),
package = "SpatialExperiment")
sample_ids <- c("section1", "section2")
samples <- file.path(dir, sample_ids)
(spe10x <- 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(1): sample_id
## reducedDimNames(0):
## mainExpName: NULL
## altExpNames(0):
## spatialData names(3) : in_tissue array_row array_col
## spatialCoords names(2) : pxl_col_in_fullres pxl_row_in_fullres
## imgData names(4): sample_id image_id data scaleFactor
To demonstrate how to accommodate molecule-based ST data
(e.g. seqFISH platform) inside a SpatialExperiment
object,
we generate some mock data of 1000 molecule coordinates across
50 genes and 20 cells. These should be formatted into a data.frame
where each row corresponds to a molecule, and columns specify the
xy-positions as well as which gene/cell the molecule has been assigned to:
n <- 1e3 # number of molecules
ng <- 50 # number of genes
nc <- 20 # number of cells
# sample xy-coordinates in [0, 1]
x <- runif(n)
y <- runif(n)
# assign each molecule to some gene-cell pair
gs <- paste0("gene", seq(ng))
cs <- paste0("cell", seq(nc))
gene <- sample(gs, n, TRUE)
cell <- sample(cs, n, TRUE)
# assure gene & cell are factors so that
# missing observations aren't dropped
gene <- factor(gene, gs)
cell <- factor(cell, cs)
# construct data.frame of molecule coordinates
df <- data.frame(gene, cell, x, y)
head(df)
## gene cell x y
## 1 gene48 cell11 0.03539863 0.64465423
## 2 gene15 cell12 0.78304050 0.27896636
## 3 gene40 cell9 0.38427813 0.27187383
## 4 gene31 cell2 0.33349229 0.67031502
## 5 gene45 cell6 0.55584953 0.98536730
## 6 gene23 cell20 0.20483172 0.09803537
Next, it is possible to re-shape the above table into a
BumpyMatrix using splitAsBumpyMatrix()
, which takes
as input the xy-coordinates, as well as arguments specifying the row and column
index of each observation:
# construct 'BumpyMatrix'
library(BumpyMatrix)
mol <- splitAsBumpyMatrix(
df[, c("x", "y")],
row = gene, col = cell)
Finally, it is possible to construct a SpatialExperiment
object with two data
slots:
counts
assay stores the number of molecules per gene and cellmolecules
assay holds the spatial molecule positions (xy-coordinates)speEach entry in the molecules
assay is a DFrame
that contains the positions
of all molecules from a given gene that have been assigned to a given cell.
# get count matrix
y <- with(df, table(gene, cell))
y <- as.matrix(unclass(y))
y[1:5, 1:5]
## cell
## gene cell1 cell2 cell3 cell4 cell5
## gene1 0 2 2 1 0
## gene2 2 1 2 1 0
## gene3 0 1 2 0 1
## gene4 2 0 0 1 2
## gene5 0 1 1 0 1
# construct SpatialExperiment
spe <- SpatialExperiment(
assays = list(
counts = y,
molecules = mol))
spe
## 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):
## spatialData names(0) :
## spatialCoords names(0) :
## imgData names(0):
The BumpyMatrix
of molecule locations can be accessed
using the dedicated molecules()
accessor:
molecules(spe)
## 50 x 20 BumpyDataFrameMatrix
## rownames: gene1 gene2 ... gene49 gene50
## colnames: cell1 cell2 ... cell19 cell20
## preview [1,1]:
## DataFrame with 0 rows and 2 columns
Subsetting objects is automatically defined to synchronize across all attributes, as for any other Bioconductor Experiment class.
For example, it is possible to subset
by sample_id
as follows:
sub <- spe10x[, spe10x$sample_id == "section1"]
Or to retain only observations that map to tissue via:
sub <- spe10x[, spatialData(spe10x)$in_tissue]
sum(spatialData(spe10x)$in_tissue) == ncol(sub)
## [1] TRUE
To work with multiple samples, the SpatialExperiment
class provides the cbind
method, which assumes unique sample_id
(s) are provided for each sample.
In case the sample_id
(s) are duplicated across multiple samples, the cbind
method takes care of this by appending indices to create unique sample identifiers.
spe1 <- spe2 <- spe
spe3 <- cbind(spe1, spe2)
## 'sample_id's are duplicated across 'SpatialExperiment' objects to cbind; appending sample indices.
unique(spe3$sample_id)
## [1] "sample01.1" "sample01.2"
Alternatively (and preferentially), we can create unique
sample_id
(s) prior to cbind
ing as follows:
# make sample identifiers unique
spe1 <- spe2 <- spe
spe1$sample_id <- paste(spe1$sample_id, "A", sep = ".")
spe2$sample_id <- paste(spe2$sample_id, "B", sep = ".")
# combine into single object
spe3 <- cbind(spe1, spe2)
In particular, when trying to replace the sample_id
(s) of a SpatialExperiment
object, these must map uniquely with the already existing ones, otherwise an
error is returned.
new <- spe3$sample_id
new[1] <- "section2.A"
spe3$sample_id <- new
## Error in .local(x, ..., value): Number of unique 'sample_id's is 2, but 3 were provided.
new[1] <- "third.one.of.two"
spe3$sample_id <- new
## Error in .local(x, ..., value): Number of unique 'sample_id's is 2, but 3 were provided.
Importantly, the sample_id
colData
field is protected, i.e.,
it will be retained upon attempted removal (= replacement by NULL
):
# backup original sample IDs
tmp <- spe$sample_id
# try to remove sample IDs
spe$sample_id <- NULL
# sample IDs remain unchanged
identical(tmp, spe$sample_id)
## [1] TRUE
sessionInfo()
## R version 4.1.1 (2021-08-10)
## Platform: x86_64-pc-linux-gnu (64-bit)
## Running under: Ubuntu 20.04.3 LTS
##
## Matrix products: default
## BLAS: /home/biocbuild/bbs-3.14-bioc/R/lib/libRblas.so
## LAPACK: /home/biocbuild/bbs-3.14-bioc/R/lib/libRlapack.so
##
## 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
##
## attached base packages:
## [1] stats4 stats graphics grDevices utils datasets methods
## [8] base
##
## other attached packages:
## [1] BumpyMatrix_1.2.0 SpatialExperiment_1.4.0
## [3] SingleCellExperiment_1.16.0 SummarizedExperiment_1.24.0
## [5] Biobase_2.54.0 GenomicRanges_1.46.0
## [7] GenomeInfoDb_1.30.0 IRanges_2.28.0
## [9] S4Vectors_0.32.0 BiocGenerics_0.40.0
## [11] MatrixGenerics_1.6.0 matrixStats_0.61.0
## [13] BiocStyle_2.22.0
##
## loaded via a namespace (and not attached):
## [1] Rcpp_1.0.7 locfit_1.5-9.4
## [3] lattice_0.20-45 digest_0.6.28
## [5] R6_2.5.1 evaluate_0.14
## [7] highr_0.9 sparseMatrixStats_1.6.0
## [9] zlibbioc_1.40.0 rlang_0.4.12
## [11] curl_4.3.2 jquerylib_0.1.4
## [13] magick_2.7.3 R.utils_2.11.0
## [15] R.oo_1.24.0 Matrix_1.3-4
## [17] rmarkdown_2.11 BiocParallel_1.28.0
## [19] stringr_1.4.0 RCurl_1.98-1.5
## [21] beachmat_2.10.0 DelayedArray_0.20.0
## [23] HDF5Array_1.22.0 compiler_4.1.1
## [25] xfun_0.27 DropletUtils_1.14.0
## [27] htmltools_0.5.2 GenomeInfoDbData_1.2.7
## [29] bookdown_0.24 edgeR_3.36.0
## [31] codetools_0.2-18 rhdf5filters_1.6.0
## [33] bitops_1.0-7 R.methodsS3_1.8.1
## [35] grid_4.1.1 jsonlite_1.7.2
## [37] formatR_1.11 magrittr_2.0.1
## [39] dqrng_0.3.0 stringi_1.7.5
## [41] scuttle_1.4.0 XVector_0.34.0
## [43] limma_3.50.0 bslib_0.3.1
## [45] DelayedMatrixStats_1.16.0 Rhdf5lib_1.16.0
## [47] rjson_0.2.20 tools_4.1.1
## [49] parallel_4.1.1 fastmap_1.1.0
## [51] yaml_2.2.1 rhdf5_2.38.0
## [53] BiocManager_1.30.16 knitr_1.36
## [55] sass_0.4.0