Combining SpatialExperiment objects
Source:R/SpatialExperiment-combine.R
SpatialExperiment-combine.Rd
The SpatialExperiment
class provides modified methods to
combine multiple SpatialExperiment
objects by column, for example from
multiple samples. These methods ensure that all data fields remain
synchronized when samples are added or removed.
Usage
# S4 method for class 'SpatialExperiment'
cbind(..., deparse.level = 1)
Arguments
- ...
a list of
SpatialExperiment
objects- deparse.level
refer to
?rbind
Value
a combined SpatialExperiment
object
Combining
The ...
argument is assumed to contain one or more
SpatialExperiment
objects.
cbind(..., deparse.level=1)
:Returns a
SpatialExperiment
where all objects in...
are combined column-wise, i.e., columns in successive objects are appended to the first object.Each
SpatialExperiment
object in...
must have the samecolData
(with the samespatialCoords
). If multiple objects use the samesample_id
, the method will proceed by assigning uniquesample_id
s.Additionally, the method combines
imgData
by row usingrbind
.Refer to
?"cbind,SingleCellExperiment-method"
for details on how metadata and other inherited attributes are combined in the output object.Refer to
?cbind
for the interpretation ofdeparse.level
.
Examples
example(read10xVisium, echo = FALSE)
# merging with duplicated 'sample_id's
# will automatically assign unique 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] "section1.1" "section2.1" "section1.2" "section2.2"
# assign unique sample identifiers
spe1 <- spe2 <- spe
spe1$sample_id <- paste(spe1$sample_id, "sample1", sep = ".")
spe2$sample_id <- paste(spe2$sample_id, "sample2", sep = ".")
# combine into single object
spe <- cbind(spe1, spe2)
# view joint 'imgData'
imgData(spe)
#> DataFrame with 4 rows and 4 columns
#> sample_id image_id data scaleFactor
#> <character> <character> <list> <numeric>
#> 1 section1.sample1 lowres #### 0.0510334
#> 2 section2.sample1 lowres #### 0.0510334
#> 3 section1.sample2 lowres #### 0.0510334
#> 4 section2.sample2 lowres #### 0.0510334
# 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.sample1 section1.sample2 section2.sample1 section2.sample2
#> FALSE 28 28 27 27
#> TRUE 22 22 22 22