SpatialExperiment-coercion | R Documentation |
The SpatialExperiment
class inherits from the
SingleCellExperiment
class making it necessary to coerce between these
classes.
To do so, we designed two different methods: the traditional as
method
and the toSpatialExperiment
function (recommended).
The as
method checks if the SingleCellExperiment
object has
already populated int_colData
with three elements:
spatialData
, spatialCoords
, and imgData
.
It also checks if colData
already contains a sample_id
.
In case these checks pass the new SpatialExperiment
will have the same
values as the SingleCellExperiment
passed object.
Otherwise a SpatialExperiment
with default values for these slots
will be created.
The toSpatialExperiment
method expects a SingleCellExperiment
object
and additional arguments as explained in the related section of this
documentation. In case the SingleCellExperiment
object has already
populated int_colData
with spatialData
and/or
spatialCoords
and/or imgData
, these will be respectively
overwritten in case the arguments spatialData
/spatialDataNames
and/or spatialCoords
/spatialCoordsNames
and/or imgData
are not NULL
.
sce |
A |
sample_id |
A |
spatialCoordsNames |
A |
spatialCoords |
A numeric matrix containing columns of spatial
coordinates, which will be accessible with |
scaleFactors |
Optional scale factors associated with the image(s). This
can be provided as a numeric value, numeric vector, list, or file path to a
JSON file for the 10x Genomics Visium platform. For 10x Genomics Visium,
the correct scale factor will automatically be selected depending on the
resolution of the image from |
imgData |
Optional |
imageSources |
Optional file path(s) or URL(s) for one or more image sources. |
image_id |
Optional character vector (same length as
|
loadImage |
Logical indicating whether to load image into memory.
Default = |
spatialDataNames |
(Deprecated) A |
spatialData |
(Deprecated) A |
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 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"))
# read in image data
img <- readImgData(
path = file.path(dir, "spatial"),
sample_id = "sample01")
## as method
(spe <- as(sce, "SpatialExperiment"))
colData(sce) <- DataFrame(xyz[,c(1:4)])
int_colData(sce)$spatialCoords <- as.matrix(xyz[,c(5,6)])
## Coercing an sce without imgData
(spe <- as(sce, "SpatialExperiment"))
## Coercing an sce with imgData
int_colData(sce)$imgData <- img
(spe <- as(sce, "SpatialExperiment"))
## toSpatialExperiment method
colData(sce) <- DataFrame(xyz)
(spe <- toSpatialExperiment(sce,
imgData = img,
spatialCoordsNames = c("pxl_col_in_fullres", "pxl_row_in_fullres"),
sample_id = "sample01"))
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.