R/patchDetection.R
patchDetection.Rd
Function to detect spatial clusters of defined types of cells. By defining a certain distance threshold, all cells within the vicinity of these clusters are detected as well.
patchDetection(
object,
patch_cells,
colPairName,
min_patch_size = 1,
name = "patch_id",
expand_by = 0,
coords = c("Pos_X", "Pos_Y"),
convex = FALSE,
img_id = NULL,
BPPARAM = SerialParam()
)
a SingleCellExperiment
or SpatialExperiment
object
logical vector of length equal to the number of cells
contained in object
. TRUE
entries define the cells to consider
for patch detection (see Details).
single character indicating the colPair(object)
entry containing the neighbor information.
single integer indicating the minimum number of connected cells that make up a patch before expansion.
single character specifying the colData
entry storing
the patch IDs in the returned object.
single numeric indicating in which vicinity range cells should be considered as belonging to the patch (see Details).
character vector of length 2 specifying the names of the
colData
(for a SingleCellExperiment
object) or the
spatialCoords
entries of the cells' x and y locations.
should the convex hull be computed before expansion? Default: the concave hull is computed.
single character indicating the colData(object)
entry
containing the unique image identifiers.
a BiocParallelParam-class
object
defining how to parallelize computations.
An object of class(object)
containing a patch ID for each
cell in colData(object)[[name]]
. If expand_by > 0
, cells in the
output object are grouped by entries in img_id
.
This function works as follows:
1. Only cells defined by patch_cells
are considered for patch
detection.
2. Patches of connected cells are detected. Here, cell-to-cell connections
are defined by the interaction graph stored in
colPair(object, colPairName)
. At this point, patches that contain
fewer than min_patch_size
cells are removed.
3. If expand_by > 0
, a concave (default) or convex hull is constructed
around each patch. This is is then expanded by expand_by
and cells
within the expanded hull are detected and assigned to the patch. This
expansion only works if a patch contains at least 3 cells.
The returned object contains an additional entry
colData(object)[[name]]
, which stores the patch ID per cell. NA
indicate cells that are not part of a patch.
If expand_by > 0
, the patchDetection
function operates on individual images.
Therefore the returned object is grouped by entries in img_id
.
This means all cells of a given image are grouped together in the object.
The ordering of cells within each individual image is the same as the ordering
of these cells in the input object.
If expand_by = 0
, the ordering of cells in the output object is the same as
in the input object.
library(cytomapper)
data(pancreasSCE)
# Visualize cell types
plotSpatial(pancreasSCE,
img_id = "ImageNb",
node_color_by = "CellType",
scales = "free")
# Build interaction graph
pancreasSCE <- buildSpatialGraph(pancreasSCE, img_id = "ImageNb",
type = "expansion", threshold = 20)
#> The returned object is ordered by the 'ImageNb' entry.
# Detect patches of "celltype_B" cells
pancreasSCE <- patchDetection(pancreasSCE,
patch_cells = pancreasSCE$CellType == "celltype_B",
colPairName = "expansion_interaction_graph")
plotSpatial(pancreasSCE,
img_id = "ImageNb",
node_color_by = "patch_id",
scales = "free")
# Include cells in vicinity
pancreasSCE <- patchDetection(pancreasSCE,
patch_cells = pancreasSCE$CellType == "celltype_B",
colPairName = "expansion_interaction_graph",
expand_by = 20,
img_id = "ImageNb")
#> The returned object is ordered by the 'ImageNb' entry.
plotSpatial(pancreasSCE,
img_id = "ImageNb",
node_color_by = "patch_id",
scales = "free")