Calculate SpatialFeatures from a MoleculeExperiment object
Source:R/spatialFeatures.R
spatialFeatures.Rd
This function takes a MoleculeExperiment (ME) object as input and calculates new molecule-based spatial features according to the following segmentations:
Sub-sector polygons
Sub-concentric polygons
Super-sector polygons
Super-concentric polygons
The function has 3 key steps: 1. adds new boundaries to the ME object according to the subcellular and supercellular segmentations, 2. calculates entropy for each of the cell and featureType combinations, and 3. combines the entropy values into a SingleCellExperiment object.
Usage
spatialFeatures(
me,
featureTypes = c("subsector", "subconcentric", "supersector", "superconcentric"),
k = 5,
nCores = 1,
includeCounts = FALSE,
concatenateFeatures = FALSE,
...
)
Arguments
- me
A MoleculeExperiment (ME) object
- featureTypes
A character string specifying the feature type. Supported values include "subsector", "subconcentric", "supersector", and "superconcentric".
- k
A numeric value indicating the number of polygons to calculate entropy (default 5)
- nCores
number of cores for parallel processing (default 1)
- includeCounts
logical (default FALSE) whether to include gene counts as features
- concatenateFeatures
logical whether to concatenate all the features into a single assay (default FALSE). If FALSE the output SE object has multiple assays
- ...
arguments passed to
loadBoundaries
andEntropyMatrix
Examples
data(example_me)
se <- spatialFeatures(me)
se
#> class: SingleCellExperiment
#> dim: 178 14
#> metadata(0):
#> assays(4): subsector subconcentric supersector superconcentric
#> rownames(178): 2010300C02Rik Acsbg1 ... Zfp536 Zfpm2
#> rowData names(0):
#> colnames(14): sample1.67500 sample1.67512 ... sample2.65070
#> sample2.65071
#> colData names(5): Cell Sample_id x_central y_central boundaries
#> reducedDimNames(0):
#> mainExpName: NULL
#> altExpNames(0):