The SpatialExperiment
class is an R/Bioconductor S4 class for storing
data from spatial -omics experiments. The class
extends the SingleCellExperiment
class for single-cell data to support
storage and retrieval of additional information from spot-based and
molecule-based platforms, including spatial coordinates, images, and
image metadata. A specialized constructor function is included for data
from the 10x Genomics Visium platform.
The following schematic illustrates the SpatialExperiment
class
structure.
As shown, an object consists of: (i) assays
containing expression counts,
(ii) rowData
containing information on features, i.e. genes, (iii)
colData
containing information on spots or cells, including nonspatial
and spatial metadata, (iv) spatialCoords
containing spatial coordinates,
and (v) imgData
containing image data. For spot-based ST data (e.g. 10x
Genomics Visium), a single assay
named counts
is used. For molecule-based
ST data (e.g. seqFISH), two assays
named counts
and molecules
can be used.
Additional details on the class structure are provided in our preprint.
For demonstration of the general class structure, we load an example
SpatialExperiment
(abbreviated as SPE) object (variable spe
):
library(SpatialExperiment)
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(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
spatialCoords
In addition to observation metadata stored inside the colData
slot,
the SpatialExperiment
class stores spatial coordinates as:
spatialCoords
, a numeric matrix of spatial coordinates (e.g. x
and y
)spatialCoords
are stored inside the int_colData
, and are directly
accessible via the corresponding accessor:
head(spatialCoords(spe))
## pxl_col_in_fullres pxl_row_in_fullres
## AAACAACGAATAGTTC-1 2312 1252
## AAACAAGTATCTCCCA-1 8230 7237
## AAACAATCTACTAGCA-1 4170 1611
## AAACACCAATAACTGC-1 2519 8315
## AAACAGAGCGACTCCT-1 7679 2927
## AAACAGCTTTCAGAAG-1 1831 6400
The corresponding column names can be also accessed with spatialCoordsNames()
:
spatialCoordsNames(spe)
## [1] "pxl_col_in_fullres" "pxl_row_in_fullres"
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))
## 576 x 600 (width x height) StoredSpatialImage
## imgSource():
## /home/biocbuild/bbs-3.20-bioc/tmpdir/RtmpH7Cpoy/Rinste5c022a45aa32/SpatialExperi
## ment/extdata/10xVisium/section1/outs/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))
Image 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:
spatialCoords
, a numeric matrix
containing spatial coordinatesspatialCoordsNames
, a character vector specifying whichcolData
fields correspond to spatial coordinatesspatialCoords
can be supplied via colData
by specifying the column names that correspond to spatial coordinates
with 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,
spatialCoordsNames = c("x", "y"))
Alternatively, spatialCoords
may be supplied separately
as a matrix
, e.g.:
xy <- as.matrix(cd[, c("x", "y")])
spe2 <- SpatialExperiment(
assay = y,
colData = cd["z"],
spatialCoords = xy)
Importantly, both of the above SpatialExperiment()
function calls
lead to construction of the exact same object:
identical(spe1, spe2)
## [1] TRUE
Finally, spatialCoords(Names)
are optional, i.e.,
we can construct a SPE using only a subset of the above arguments:
spe <- SpatialExperiment(
assays = y)
isEmpty(spatialCoords(spe))
## [1] TRUE
In general, spatialCoordsNames
takes precedence over spatialCoords
,
i.e., if both are supplied, the latter will be ignored. In other words,
spatialCoords
are preferentially extracted from the DataFrame
provided via 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
.
In this case, a message is also provided:
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
colData
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", "outs"),
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 = DataFrame(xyz),
spatialCoordsNames = c("pxl_col_in_fullres", "pxl_row_in_fullres"),
imgData = img,
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
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,
as illustrated below (where “raw/filtered” refers to either “raw” or
“filtered” to match the data
argument). Note that the base directory
“outs/” from Space Ranger can either be included manually in the paths
provided in the samples
argument, 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
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, "outs")
(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(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
Or alternatively, omitting the base directory outs/
from Space Ranger:
samples2 <- file.path(dir, sample_ids)
(spe10x2 <- 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
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 gene14 cell4 0.3828213 0.3917898
## 2 gene6 cell5 0.1932988 0.3310752
## 3 gene17 cell3 0.5967864 0.6500672
## 4 gene9 cell19 0.7395685 0.4661870
## 5 gene7 cell2 0.1744469 0.2706218
## 6 gene11 cell14 0.7718313 0.4556108
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)Each 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 2 1 1 1 0
## gene2 0 2 2 2 1
## gene3 1 0 0 4 1
## gene4 0 1 3 1 1
## gene5 0 1 1 1 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):
## 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 2 rows and 2 columns
## x y
## <numeric> <numeric>
## 1 0.0487952 0.707750
## 2 0.3682675 0.130663
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[, colData(spe10x)$in_tissue]
sum(colData(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
Both the SpatialImage
(SpI) and SpatialExperiment
(SpE) class currently support
two basic image transformations, namely, rotation (via rotateImg()
) and
mirroring (via mirrorImg()
). Specifically, for a SpI/E x
:
rotateImg(x, degrees)
expects as degrees
a single numeric in +/-[0,90,…,360].mirrorImg(x, axis)
expects as axis
a character string specifying"h"
) or vertically ("v"
).Here, we apply various transformations using both a SpI (spi
) and SpE (spe
)
as input, and compare the resulting images to the original:
# extract first image
spi <- getImg(spe10x)
# apply counter-/clockwise rotation
spi1 <- rotateImg(spi, -90)
spi2 <- rotateImg(spi, +90)
# visual comparison
par(mfrow = c(1, 3))
plot(as.raster(spi))
plot(as.raster(spi1))
plot(as.raster(spi2))
# specify sample & image identifier
sid <- "section1"
iid <- "lowres"
# counter-clockwise rotation
tmp <- rotateImg(spe10x,
sample_id = sid,
image_id = iid,
degrees = -90)
# visual comparison
par(mfrow = c(1, 2))
plot(imgRaster(spe10x, sid, iid))
plot(imgRaster(tmp, sid, iid))
# extract first image
spi <- getImg(spe10x)
# mirror horizontally/vertically
spi1 <- mirrorImg(spi, "h")
spi2 <- mirrorImg(spi, "v")
# visual comparison
par(mfrow = c(1, 3))
plot(as.raster(spi))
plot(as.raster(spi1))
plot(as.raster(spi2))
# specify sample & image identifier
sid <- "section2"
iid <- "lowres"
# mirror horizontally
tmp <- mirrorImg(spe10x,
sample_id = sid,
image_id = iid,
axis = "h")
# visual comparison
par(mfrow = c(1, 2))
plot(imgRaster(spe10x, sid, iid))
plot(imgRaster(tmp, sid, iid))
sessionInfo()
## R version 4.4.1 (2024-06-14)
## Platform: x86_64-pc-linux-gnu
## Running under: Ubuntu 24.04.1 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.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] BumpyMatrix_1.14.0 SpatialExperiment_1.16.0
## [3] SingleCellExperiment_1.28.0 SummarizedExperiment_1.36.0
## [5] Biobase_2.66.0 GenomicRanges_1.58.0
## [7] GenomeInfoDb_1.42.0 IRanges_2.40.0
## [9] S4Vectors_0.44.0 BiocGenerics_0.52.0
## [11] MatrixGenerics_1.18.0 matrixStats_1.4.1
## [13] BiocStyle_2.34.0
##
## loaded via a namespace (and not attached):
## [1] rjson_0.2.23 xfun_0.48
## [3] bslib_0.8.0 rhdf5_2.50.0
## [5] lattice_0.22-6 rhdf5filters_1.18.0
## [7] tools_4.4.1 curl_5.2.3
## [9] parallel_4.4.1 highr_0.11
## [11] R.oo_1.26.0 Matrix_1.7-1
## [13] sparseMatrixStats_1.18.0 dqrng_0.4.1
## [15] lifecycle_1.0.4 GenomeInfoDbData_1.2.13
## [17] compiler_4.4.1 tinytex_0.53
## [19] statmod_1.5.0 codetools_0.2-20
## [21] htmltools_0.5.8.1 sass_0.4.9
## [23] yaml_2.3.10 crayon_1.5.3
## [25] jquerylib_0.1.4 R.utils_2.12.3
## [27] BiocParallel_1.40.0 DelayedArray_0.32.0
## [29] cachem_1.1.0 limma_3.62.0
## [31] magick_2.8.5 abind_1.4-8
## [33] locfit_1.5-9.10 digest_0.6.37
## [35] bookdown_0.41 fastmap_1.2.0
## [37] grid_4.4.1 cli_3.6.3
## [39] SparseArray_1.6.0 magrittr_2.0.3
## [41] S4Arrays_1.6.0 edgeR_4.4.0
## [43] DelayedMatrixStats_1.28.0 UCSC.utils_1.2.0
## [45] rmarkdown_2.28 XVector_0.46.0
## [47] httr_1.4.7 DropletUtils_1.26.0
## [49] R.methodsS3_1.8.2 beachmat_2.22.0
## [51] HDF5Array_1.34.0 evaluate_1.0.1
## [53] knitr_1.48 rlang_1.1.4
## [55] Rcpp_1.0.13 scuttle_1.16.0
## [57] formatR_1.14 BiocManager_1.30.25
## [59] jsonlite_1.8.9 R6_2.5.1
## [61] Rhdf5lib_1.28.0 zlibbioc_1.52.0