View source: R/spatialEnhance.R
spatialEnhance | R Documentation |
Enhanced clustering of a spatial expression dataset to subspot resolution.
spatialEnhance(
sce,
q,
platform = c("Visium", "VisiumHD", "ST"),
use.dimred = "PCA",
d = 15,
nsubspots.per.edge = 3,
init = NULL,
init.method = c("spatialCluster", "mclust", "kmeans"),
model = c("t", "normal"),
nrep = 1e+05,
gamma = NULL,
mu0 = NULL,
lambda0 = NULL,
alpha = 1,
beta = 0.01,
save.chain = FALSE,
chain.fname = NULL,
burn.in = 10000,
thin = 100,
jitter.scale = 5,
jitter.prior = 0.3,
adapt.before = burn.in,
cores = 1,
verbose = FALSE
)
coreTune(sce, test.cores = detectCores(), test.times = 1, ...)
adjustClusterLabels(sce, burn.in)
sce |
A SingleCellExperiment object containing the spatial data. |
q |
The number of clusters. |
platform |
Spatial transcriptomic platform. Specify 'Visium' for hex
lattice geometry or 'ST' and 'VisiumHD' for square lattice geometry.
Specifying this parameter is optional when analyzing SingleCellExperiments
processed using |
use.dimred |
Name of a reduced dimensionality result in
|
d |
Number of top principal components to use when clustering. |
nsubspots.per.edge |
Number of subspots per edge of the square. Only
valid when |
init |
Initial cluster assignments for spots. |
init.method |
If |
model |
Error model. ('normal' or 't') |
nrep |
The number of MCMC iterations. |
gamma |
Smoothing parameter. (Values in range of 1-3 seem to work well.) |
mu0 |
Prior mean hyperparameter for mu. If not provided, mu0 is set to the mean of PCs over all spots. |
lambda0 |
Prior precision hyperparam for mu. If not provided, lambda0
is set to a diagonal matrix |
alpha |
Hyperparameter for Wishart distributed precision lambda. |
beta |
Hyperparameter for Wishart distributed precision lambda. |
save.chain |
If true, save the MCMC chain to an HDF5 file. |
chain.fname |
File path for saved chain. Tempfile used if not provided. |
burn.in |
Number of iterations to exclude as burn-in period. The MCMC
iterations are currently thinned to every 100; accordingly |
thin |
Thinning rate. |
jitter.scale |
Controls the amount of jittering. Small amounts of
jittering are more likely to be accepted but result in exploring the space
more slowly. We suggest tuning |
jitter.prior |
Scale factor for the prior variance, parameterized as the
proportion (default = 0.3) of the mean variance of the PCs.
We suggest making |
adapt.before |
Adapting the MCMC chain before the specified number
or proportion of iterations (by default equal to |
cores |
The number of threads to use. The results are invariate to the
value of |
verbose |
Log progress to stderr. |
test.cores |
Either a list of, or a maximum number of cores to test. In the latter case, a list of values (power of 2) will be created |
test.times |
Times to repeat the benchmarking with microbenchmark. |
... |
Arguments for |
The enhanced SingleCellExperiment
has most of the properties of the
input SCE - rowData
, colData
, reducedDims
- but does
not include expression data in counts
or logcounts
. To impute
enhanced expression vectors, please use [enhanceFeatures()] after
running spatialEnhance
.
The colData
of the enhanced SingleCellExperiment
includes the
following columns to permit referencing the subspots in spatial context and
linking back to the original spots:
spot.idx
: Index of the spot this subspot belongs to (with
respect to the input SCE).
subspot.idx
: Index of the subspot within its parent spot.
spot.row
: Array row of the subspot's parent spot.
spot.col
: Array col of the subspot's parent spot.
array_row
: Array row of the subspot. This is the parent spot's row
plus an offset based on the subspot's position within the spot.
array_col
: Array col of the subspot. This is the parent spot's col
plus an offset based on the subspot's position within the spot.
pxl_row_in_fullres
: Pixel row of the subspot. This is the parent spot's
row plus an offset based on the subspot's position within the spot.
pxl_col_in_fullres
: Pixel col of the subspot. This is the parent spot's
col plus an offset based on the subspot's position within the spot.
spatialEnhance
returns a new SingleCellExperiment object.
By default, the assays
of this object are empty, and the enhanced
resolution PCs are stored as a reduced dimensionality result accessible
with reducedDim(sce, 'PCA')
.
coresTune
returns the output of microbenchmark
.
adjustClusterLabels
adjusts the cluster labels from the MCMC samples
via burn.in
, the percentage of samples to drop. The MCMC chain
must be retained.
spatialCluster
for clustering at the spot level
before enhancing, clusterPlot
for visualizing the cluster
assignments, enhanceFeatures
for imputing enhanced
expression, and mcmcChain
for examining the full MCMC chain
associated with the enhanced clustering.
.
set.seed(149)
sce <- exampleSCE()
sce <- spatialCluster(sce, 7, nrep = 100, burn.in = 10)
enhanced <- spatialEnhance(sce, 7, nrep = 100, burn.in = 10)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.