This function loads a SpatialExperiment object (or a
list with several SpatialExperiment objects) into a
SpatialDDLS object.
Usage
loadSTProfiles(
object,
st.data,
st.spot.ID.column,
st.gene.ID.column,
st.min.counts = 0,
st.min.spots = 0,
st.n.slides = 3,
verbose = TRUE
)Arguments
- object
A
SpatialDDLSobject.- st.data
A
SpatialExperimentobject (or a list with severalSpatialExperimentobjects) to be deconvoluted.- st.spot.ID.column
Name or number of the column in spots metadata corresponding to spot names in the expression matrix.
- st.gene.ID.column
Name or number of the column in genes metadata corresponding to names used for features/genes.
- st.min.counts
Minimum gene counts to filter (0 by default).
- st.min.spots
Minimum of spots with more than
min.counts(0 by default).- st.n.slides
Minimum number of slides (
SpatialExperimentobjects) in which a gene has to be expressed in order to keep it. This parameter is applicable only when multipleSpatialExperimentobjects are provided. Genes not present in at leastst.n.slideswill be discarded. If no filtering is desired, setst.n.slides = 1.- verbose
Show informative messages during execution (
TRUEby default).
Value
A SpatialDDLS object with the provided spatial
trainscriptomics data loaded into the spatial.experiments slot.
Details
It is recommended to perform this step when creating the
SpatialDDLS object using the
createSpatialDDLSobject function in order to only keep genes
shared between the spatial transcriptomics and the single-cell
transcriptomics data used as reference. In addition, please, make sure the
gene identifiers used the spatial and single-cell transcriptomics data are
consistent.
Examples
# \donttest{
set.seed(123)
sce <- SingleCellExperiment::SingleCellExperiment(
assays = list(
counts = matrix(
rpois(100, lambda = 5), nrow = 40, ncol = 30,
dimnames = list(paste0("Gene", seq(40)), paste0("RHC", seq(30)))
)
),
colData = data.frame(
Cell_ID = paste0("RHC", seq(30)),
Cell_Type = sample(x = paste0("CellType", seq(4)), size = 30,
replace = TRUE)
),
rowData = data.frame(
Gene_ID = paste0("Gene", seq(40))
)
)
SDDLS <- createSpatialDDLSobject(
sc.data = sce,
sc.cell.ID.column = "Cell_ID",
sc.gene.ID.column = "Gene_ID",
sc.filt.genes.cluster = FALSE
)
#> === Spatial transcriptomics data not provided
#> === Processing single-cell data
#> - Filtering features:
#> - Selected features: 40
#> - Discarded features: 0
#>
#> === No mitochondrial genes were found by using ^mt- as regrex
#>
#> === Final number of dimensions for further analyses: 40
## simulating a SpatialExperiment object
counts <- matrix(rpois(30, lambda = 5), ncol = 6)
rownames(counts) <- paste0("Gene", 1:5)
colnames(counts) <- paste0("Spot", 1:6)
coordinates <- matrix(
c(1, 2, 3, 1, 2, 3, 1, 2, 3, 1, 2, 3), ncol = 2
)
ste <- SpatialExperiment::SpatialExperiment(
assays = list(counts = as.matrix(counts)),
rowData = data.frame(Gene_ID = paste0("Gene", 1:5)),
colData = data.frame(Cell_ID = paste0("Spot", 1:6)),
spatialCoords = coordinates
)
## previous SpatialDDLS object
SDDLS <- loadSTProfiles(
object = SDDLS,
st.data = ste,
st.spot.ID.column = "Cell_ID",
st.gene.ID.column = "Gene_ID"
)
#> === 1 SpatialExperiment objects provided
#> === Processing spatial transcriptomics data
#> - Filtering features:
#> - Selected features: 5
#> - Discarded features: 0
#>
# }
