Title: | Methods for imaging mass cytometry data analysis |
---|---|
Description: | This R package supports the handling and analysis of imaging mass cytometry and other highly multiplexed imaging data. The main functionality includes reading in single-cell data after image segmentation and measurement, data formatting to perform channel spillover correction and a number of spatial analysis approaches. First, cell-cell interactions are detected via spatial graph construction; these graphs can be visualized with cells representing nodes and interactions representing edges. Furthermore, per cell, its direct neighbours are summarized to allow spatial clustering. Per image/grouping level, interactions between types of cells are counted, averaged and compared against random permutations. In that way, types of cells that interact more (attraction) or less (avoidance) frequently than expected by chance are detected. |
Authors: | Nils Eling [aut], Tobias Hoch [ctb], Vito Zanotelli [ctb], Jana Fischer [ctb], Daniel Schulz [ctb, cre] , Lasse Meyer [ctb] |
Maintainer: | Daniel Schulz <[email protected]> |
License: | GPL-3 |
Version: | 1.13.0 |
Built: | 2024-11-29 08:30:59 UTC |
Source: | https://github.com/bioc/imcRtools |
Function to summarize categorical or expression values of all neighbors of each cell.
aggregateNeighbors( object, colPairName, aggregate_by = c("metadata", "expression"), count_by = NULL, proportions = TRUE, assay_type = NULL, subset_row = NULL, statistic = c("mean", "median", "sd", "var"), name = NULL )
aggregateNeighbors( object, colPairName, aggregate_by = c("metadata", "expression"), count_by = NULL, proportions = TRUE, assay_type = NULL, subset_row = NULL, statistic = c("mean", "median", "sd", "var"), name = NULL )
object |
a |
colPairName |
single character indicating the |
aggregate_by |
character specifying whether the neighborhood should be
summarized by cellular features stored in |
count_by |
for |
proportions |
single logical indicating whether aggregated metadata should be returned in form of proportions instead of absolute counts. |
assay_type |
for |
subset_row |
for |
statistic |
for |
name |
single character specifying the name of the data frame to be
saved in the |
returns an object of class(object)
containing the aggregated
values in form of a DataFrame
object in
colData(object)[[name]]
.
Daniel Schulz ([email protected])
library(cytomapper) data(pancreasSCE) sce <- buildSpatialGraph(pancreasSCE, img_id = "ImageNb", type = "knn", k = 3) # Aggregating neighboring cell-types sce <- aggregateNeighbors(sce, colPairName = "knn_interaction_graph", aggregate_by = "metadata", count_by = "CellType") sce$aggregatedNeighbors # Aggregating neighboring expression values sce <- aggregateNeighbors(sce, colPairName = "knn_interaction_graph", aggregate_by = "expression", assay_type = "exprs", statistic = "mean") sce$mean_aggregatedExpression
library(cytomapper) data(pancreasSCE) sce <- buildSpatialGraph(pancreasSCE, img_id = "ImageNb", type = "knn", k = 3) # Aggregating neighboring cell-types sce <- aggregateNeighbors(sce, colPairName = "knn_interaction_graph", aggregate_by = "metadata", count_by = "CellType") sce$aggregatedNeighbors # Aggregating neighboring expression values sce <- aggregateNeighbors(sce, colPairName = "knn_interaction_graph", aggregate_by = "expression", assay_type = "exprs", statistic = "mean") sce$mean_aggregatedExpression
Helper function for estimating the spillover matrix. Per metal spot, consecutive pixels a aggregated (default: summed).
binAcrossPixels( object, bin_size, spot_id = "sample_id", assay_type = "counts", statistic = "sum", ... )
binAcrossPixels( object, bin_size, spot_id = "sample_id", assay_type = "counts", statistic = "sum", ... )
object |
a |
bin_size |
single numeric indicating how many consecutive pixels per spot should be aggregated. |
spot_id |
character string indicating which |
assay_type |
character string indicating which assay to use. |
statistic |
character string indicating the statistic to use for aggregating consecutive pixels. |
... |
additional arguments passed to |
returns the binned pixel intensities in form of a
SingleCellExperiment
object
Nils Eling ([email protected])
aggregateAcrossCells
for the aggregation
function
path <- system.file("extdata/spillover", package = "imcRtools") # Read in .txt files sce <- readSCEfromTXT(path) dim(sce) # Visualizes heatmap before aggregation plotSpotHeatmap(sce) # Sum consecutive pixels sce <- binAcrossPixels(sce, bin_size = 10) dim(sce) # Visualizes heatmap after aggregation plotSpotHeatmap(sce)
path <- system.file("extdata/spillover", package = "imcRtools") # Read in .txt files sce <- readSCEfromTXT(path) dim(sce) # Visualizes heatmap before aggregation plotSpotHeatmap(sce) # Sum consecutive pixels sce <- binAcrossPixels(sce, bin_size = 10) dim(sce) # Visualizes heatmap after aggregation plotSpotHeatmap(sce)
Function to define cell-cell interactions via distance-based expansion, delaunay triangulation or k nearest neighbor detection.
buildSpatialGraph( object, img_id, type = c("expansion", "knn", "delaunay"), k = NULL, directed = TRUE, max_dist = NULL, threshold = NULL, coords = c("Pos_X", "Pos_Y"), name = NULL, BNPARAM = KmknnParam(), BPPARAM = SerialParam(), ... )
buildSpatialGraph( object, img_id, type = c("expansion", "knn", "delaunay"), k = NULL, directed = TRUE, max_dist = NULL, threshold = NULL, coords = c("Pos_X", "Pos_Y"), name = NULL, BNPARAM = KmknnParam(), BPPARAM = SerialParam(), ... )
object |
a |
img_id |
single character indicating the |
type |
single character specifying the type of graph to be build.
Supported entries are |
k |
(when |
directed |
(when |
max_dist |
(when |
threshold |
(when |
coords |
character vector of length 2 specifying the names of the
|
name |
single character specifying the name of the graph. |
BNPARAM |
a |
BPPARAM |
a |
... |
additional parameters passed to the
|
returns a SpatialExperiment
or SingleCellExperiment
containing the graph in form of a SelfHits
object in
colPair(object, name)
. The object is grouped by entries in the
img_id
slot.
This function defines interacting cells in different ways. They are based on the cells' centroids and do not incorporate cell shape or area.
1. When type = "expansion"
, all cells within the radius
threshold
are considered interacting cells.
2. When type = "delaunay"
, interacting cells are found via a delaunay
triangulation of the cells' centroids.
3. When type = "knn"
, interacting cells are defined as the k
nearest neighbors in the 2D spatial plane.
The directed
parameter only affects graph construction via k nearest
neighbor search. For directed = FALSE
, each interaction will be
stored as mutual edge (e.g. node 2 is connected to node 10 and vise versa).
For type = "expansion"
and type = "delaunay"
, each edge is
stored as mutual edge by default.
The graph is stored in form of a SelfHits
object in
colPair(object, name)
. This object can be regarded as an edgelist
and coerced to an igraph
object via
graph_from_edgelist(as.matrix(colPair(object, name)))
.
When finding interactions via expansion
or knn
, the
findNeighbors
or
findKNN
functions are used, respectively. Both
functions accept the BNPARAM
parameter via which the graph
construction method can be defined (default
KmknnParam
). For an overview on the different
algorithms, see
BiocNeighborParam
. Within the
BiocNeighborParam
object, distance
can be set to
"Euclidean"
(default), "Manhattan"
or "Cosine"
.
The buildSpatialGraph
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.
Nils Eling ([email protected])
findNeighbors
for the function finding
interactions via expansion
findKNN
for the function finding interactions
via nearest neighbor search
triangulate
for the function finding interactions
via delaunay triangulation
plotSpatial
for visualizing spatial graphs
path <- system.file("extdata/mockData/steinbock", package = "imcRtools") spe <- read_steinbock(path) # Constructing a graph via expansion spe <- buildSpatialGraph(spe, img_id = "sample_id", type = "expansion", threshold = 10) colPair(spe, "expansion_interaction_graph") # Constructing a graph via delaunay triangulation spe <- buildSpatialGraph(spe, img_id = "sample_id", type = "delaunay") colPair(spe, "delaunay_interaction_graph") # Constructing a graph via k nearest neighbor search spe <- buildSpatialGraph(spe, img_id = "sample_id", type = "knn", k = 5) colPair(spe, "knn_interaction_graph")
path <- system.file("extdata/mockData/steinbock", package = "imcRtools") spe <- read_steinbock(path) # Constructing a graph via expansion spe <- buildSpatialGraph(spe, img_id = "sample_id", type = "expansion", threshold = 10) colPair(spe, "expansion_interaction_graph") # Constructing a graph via delaunay triangulation spe <- buildSpatialGraph(spe, img_id = "sample_id", type = "delaunay") colPair(spe, "delaunay_interaction_graph") # Constructing a graph via k nearest neighbor search spe <- buildSpatialGraph(spe, img_id = "sample_id", type = "knn", k = 5) colPair(spe, "knn_interaction_graph")
Function to calculate the average number of neighbors B that a cell of type A has using different approaches.
countInteractions( object, group_by, label, colPairName, method = c("classic", "histocat", "patch"), patch_size = NULL )
countInteractions( object, group_by, label, colPairName, method = c("classic", "histocat", "patch"), patch_size = NULL )
object |
a |
group_by |
a single character indicating the |
label |
single character specifying the |
colPairName |
single character indicating the |
method |
which cell-cell interaction counting method to use (see details) |
patch_size |
if |
a DataFrame containing one row per group_by
entry and unique
label
entry combination (from_label
, to_label
). The
ct
entry stores the interaction count as described in the details.
NA
is returned if a certain label is not present in this grouping
level.
In principle, the countInteractions
function counts the number
of edges (interactions) between each set of unique entries in
colData(object)[[label]]
. Simplified, it counts for each cell of
type A the number of neighbors of type B.
This count is averaged within each unique entry
colData(object)[[group_by]]
in three different ways:
1. method = "classic"
: The count is divided by the total number of
cells of type A. The final count can be interpreted as "How many neighbors
of type B does a cell of type A have on average?"
2. method = "histocat"
: The count is divided by the number of cells
of type A that have at least one neighbor of type B. The final count can be
interpreted as "How many many neighbors of type B has a cell of type A on
average, given it has at least one neighbor of type B?"
3. method = "patch"
: For each cell, the count is binarized to 0
(less than patch_size
neighbors of type B) or 1 (more or equal to
patch_size
neighbors of type B). The binarized counts are averaged
across all cells of type A. The final count can be interpreted as "What
fraction of cells of type A have at least a given number of neighbors of
type B?"
Vito Zanotelli
Jana Fischer
adapted by Nils Eling ([email protected])
testInteractions
for testing cell-cell
interactions per grouping level.
library(cytomapper) data(pancreasSCE) pancreasSCE <- buildSpatialGraph(pancreasSCE, img_id = "ImageNb", type = "knn", k = 3) # Classic style calculation (out <- countInteractions(pancreasSCE, group_by = "ImageNb", label = "CellType", method = "classic", colPairName = "knn_interaction_graph")) # Histocat style calculation (out <- countInteractions(pancreasSCE, group_by = "ImageNb", label = "CellType", method = "histocat", colPairName = "knn_interaction_graph")) # Patch style calculation (out <- countInteractions(pancreasSCE, group_by = "ImageNb", label = "CellType", method = "patch", patch_size = 3, colPairName = "knn_interaction_graph"))
library(cytomapper) data(pancreasSCE) pancreasSCE <- buildSpatialGraph(pancreasSCE, img_id = "ImageNb", type = "knn", k = 3) # Classic style calculation (out <- countInteractions(pancreasSCE, group_by = "ImageNb", label = "CellType", method = "classic", colPairName = "knn_interaction_graph")) # Histocat style calculation (out <- countInteractions(pancreasSCE, group_by = "ImageNb", label = "CellType", method = "histocat", colPairName = "knn_interaction_graph")) # Patch style calculation (out <- countInteractions(pancreasSCE, group_by = "ImageNb", label = "CellType", method = "patch", patch_size = 3, colPairName = "knn_interaction_graph"))
Function to detect the spatial community of each cell as proposed by Jackson et al., The single-cell pathology landscape of breast cancer, Nature, 2020. Each cell is clustered based on its interactions as defined by a spatial object graph.
detectCommunity( object, colPairName, size_threshold = 0, group_by = NULL, name = "spatial_community", cluster_fun = "louvain", BPPARAM = SerialParam() )
detectCommunity( object, colPairName, size_threshold = 0, group_by = NULL, name = "spatial_community", cluster_fun = "louvain", BPPARAM = SerialParam() )
object |
a |
colPairName |
single character indicating the |
size_threshold |
single positive numeric that specifies the minimum number of cells per community. Defaults to 0. |
group_by |
single character indicating that spatial community
detection will be performed separately for all unique entries to
|
name |
single character specifying the name of the output saved in
|
cluster_fun |
single character specifying the function to use for
community detection. Options are all strings that contain the suffix of an
|
BPPARAM |
a |
returns an object of class(object)
containing a new column
entry to colData(object)[[name]]
.
1. Create an igraph object from the edge list stored in
colPair(object, colPairName)
.
2. Perform community detection using the specified cluster_fun
algorithm.
3. Store the community IDs in a vector and replace all communities with
a size smaller than size_threshold
by NA.
Optional steps: Specify group_by
to perform spatial community
detection separately for all unique entries to colData(object)[,group_by]
e.g. for tumor and non-tumor cells.
Lasse Meyer ([email protected])
Jackson et al., The single-cell pathology landscape of breast cancer, Nature, 2020
library(cytomapper) library(BiocParallel) data(pancreasSCE) sce <- buildSpatialGraph(pancreasSCE, img_id = "ImageNb", type = "expansion", name = "neighborhood", threshold = 20) ## Detect spatial community set.seed(22) sce <- detectCommunity(sce, colPairName = "neighborhood", size_threshold = 10) plotSpatial(sce, img_id = "ImageNb", node_color_by = "spatial_community", scales = "free") ## Detect spatial community - specify group_by sce <- detectCommunity(sce, colPairName = "neighborhood", group_by = "CellType", size_threshold = 10, BPPARAM = SerialParam(RNGseed = 22)) plotSpatial(sce, img_id = "ImageNb", node_color_by = "spatial_community", scales = "free")
library(cytomapper) library(BiocParallel) data(pancreasSCE) sce <- buildSpatialGraph(pancreasSCE, img_id = "ImageNb", type = "expansion", name = "neighborhood", threshold = 20) ## Detect spatial community set.seed(22) sce <- detectCommunity(sce, colPairName = "neighborhood", size_threshold = 10) plotSpatial(sce, img_id = "ImageNb", node_color_by = "spatial_community", scales = "free") ## Detect spatial community - specify group_by sce <- detectCommunity(sce, colPairName = "neighborhood", group_by = "CellType", size_threshold = 10, BPPARAM = SerialParam(RNGseed = 22)) plotSpatial(sce, img_id = "ImageNb", node_color_by = "spatial_community", scales = "free")
Function to detect the spatial context (SC) of each cell. Based on its sorted (high-to-low) cellular neighborhood (CN) fractions in a spatial interaction graph, the SC of each cell is assigned as the set of CNs that cumulatively exceed a user-defined fraction threshold.
The term was coined by Bhate S. et al., Tissue schematics map the specialization of immune tissue motifs and their appropriation by tumors, Cell Systems, 2022 and describes tissue regions in which distinct CNs may be interacting.
detectSpatialContext( object, entry = "aggregatedNeighbors", threshold = 0.9, name = "spatial_context" )
detectSpatialContext( object, entry = "aggregatedNeighbors", threshold = 0.9, name = "spatial_context" )
object |
a |
entry |
single character specifying the |
threshold |
single numeric between 0 and 1 that specifies the fraction threshold for SC assignment. Defaults to 0.9. |
name |
single character specifying the name of the output saved in
|
returns an object of class(object)
containing a new column
entry to colData(object)[[name]]
The function relies on CN fractions for each cell in a spatial interaction graph (originally a k-nearest neighbor (KNN) graph).
We can retrieve the CN fractions using the buildSpatialGraph and aggregateNeighbors functions.
The window size (k for KNN) for buildSpatialGraph should reflect a length scale on which biological signals can be exchanged and depends, among others, on cell density and tissue area. In view of their divergent functionality, we recommend to use a larger window size for SC (interaction between local processes) than for CN (local processes) detection.
Subsequently, the CN fractions are sorted from high-to-low and the SC of each cell is assigned the minimal combination of SCs that additively surpass a user-defined threshold. The default threshold of 0.9 aims to represent the dominant CNs, hence the most prevalent signals, in a given window.
For more details, please refer to: Bhate S. et al., Tissue schematics map the specialization of immune tissue motifs and their appropriation by tumors, Cell Systems, 2022.
Lasse Meyer ([email protected])
filterSpatialContext
for the function to filter
spatial contexts
plotSpatialContext
for the function to plot
spatial context graphs
set.seed(22) library(cytomapper) data(pancreasSCE) ## 1. Cellular neighborhood (CN) sce <- buildSpatialGraph(pancreasSCE, img_id = "ImageNb", type = "knn", name = "knn_cn_graph", k = 5) sce <- aggregateNeighbors(sce, colPairName = "knn_cn_graph", aggregate_by = "metadata", count_by = "CellType", name = "aggregatedCellTypes") cur_cluster <- kmeans(sce$aggregatedCellTypes, centers = 3) sce$cellular_neighborhood <- factor(cur_cluster$cluster) plotSpatial(sce, img_id = "ImageNb", colPairName = "knn_cn_graph", node_color_by = "cellular_neighborhood", scales = "free") ## 2. Spatial context (SC) sce <- buildSpatialGraph(sce, img_id = "ImageNb", type = "knn", name = "knn_sc_graph", k = 15) sce <- aggregateNeighbors(sce, colPairName = "knn_sc_graph", aggregate_by = "metadata", count_by = "cellular_neighborhood", name = "aggregatedNeighborhood") # Detect spatial context sce <- detectSpatialContext(sce, entry = "aggregatedNeighborhood", threshold = 0.9) plotSpatial(sce, img_id = "ImageNb", colPairName = "knn_sc_graph", node_color_by = "spatial_context", scales = "free")
set.seed(22) library(cytomapper) data(pancreasSCE) ## 1. Cellular neighborhood (CN) sce <- buildSpatialGraph(pancreasSCE, img_id = "ImageNb", type = "knn", name = "knn_cn_graph", k = 5) sce <- aggregateNeighbors(sce, colPairName = "knn_cn_graph", aggregate_by = "metadata", count_by = "CellType", name = "aggregatedCellTypes") cur_cluster <- kmeans(sce$aggregatedCellTypes, centers = 3) sce$cellular_neighborhood <- factor(cur_cluster$cluster) plotSpatial(sce, img_id = "ImageNb", colPairName = "knn_cn_graph", node_color_by = "cellular_neighborhood", scales = "free") ## 2. Spatial context (SC) sce <- buildSpatialGraph(sce, img_id = "ImageNb", type = "knn", name = "knn_sc_graph", k = 15) sce <- aggregateNeighbors(sce, colPairName = "knn_sc_graph", aggregate_by = "metadata", count_by = "cellular_neighborhood", name = "aggregatedNeighborhood") # Detect spatial context sce <- detectSpatialContext(sce, entry = "aggregatedNeighborhood", threshold = 0.9) plotSpatial(sce, img_id = "ImageNb", colPairName = "knn_sc_graph", node_color_by = "spatial_context", scales = "free")
Helper function for estimating the spillover matrix. After assigning each pixel to a spotted mass, this function will filter incorrectly assigned pixels and remove small pixel sets.
filterPixels( object, bc_id = "bc_id", spot_mass = "sample_mass", minevents = 40, correct_pixels = TRUE )
filterPixels( object, bc_id = "bc_id", spot_mass = "sample_mass", minevents = 40, correct_pixels = TRUE )
object |
a |
bc_id |
character string indicating which |
spot_mass |
character string indicating which |
minevents |
single numeric indicating the threshold under which pixel sets are excluded from spillover estimation. |
correct_pixels |
logical indicating if incorrectly assigned pixels should be excluded from spillover estimation. |
returns a SingleCellExperiment object in which
colData(object)$bc_id
has been adjusted based on the filter
criteria.
Vito Zanotelli, adapted by Nils Eling ([email protected])
path <- system.file("extdata/spillover", package = "imcRtools") sce <- readSCEfromTXT(path) assay(sce, "exprs") <- asinh(counts(sce)/5) # Pre-process via CATALYST library(CATALYST) bc_key <- as.numeric(unique(sce$sample_mass)) sce <- assignPrelim(sce, bc_key = bc_key) sce <- estCutoffs(sce) sce <- applyCutoffs(sce) sce <- filterPixels(sce) table(sce$sample_mass, sce$bc_id)
path <- system.file("extdata/spillover", package = "imcRtools") sce <- readSCEfromTXT(path) assay(sce, "exprs") <- asinh(counts(sce)/5) # Pre-process via CATALYST library(CATALYST) bc_key <- as.numeric(unique(sce$sample_mass)) sce <- assignPrelim(sce, bc_key = bc_key) sce <- estCutoffs(sce) sce <- applyCutoffs(sce) sce <- filterPixels(sce) table(sce$sample_mass, sce$bc_id)
Function to filter detected spatial contexts (SCs) based on a user-defined threshold for number of group entries and/or cells.
filterSpatialContext( object, entry = "spatial_context", group_by = "sample_id", group_threshold = NULL, cells_threshold = NULL, name = "spatial_context_filtered" )
filterSpatialContext( object, entry = "spatial_context", group_by = "sample_id", group_threshold = NULL, cells_threshold = NULL, name = "spatial_context_filtered" )
object |
a |
entry |
a single character specifying the |
group_by |
a single character indicating the |
group_threshold |
a single numeric specifying the minimum number of group entries in which a SC is detected. |
cells_threshold |
a single numeric specifying the minimum total number of cells in a SC. |
name |
a single character specifying the name of the output saved in
|
returns an object of class(object)
containing a new column
entry to colData(object)[[name]]
and a new data.frame
entry to
metadata(object)[["filterSpatialContext"]]
containing the group and
cell counts per SC.
Lasse Meyer ([email protected])
detectSpatialContext
for the function to detect
spatial contexts
plotSpatialContext
for the function to plot
spatial context graphs
set.seed(22) library(cytomapper) data(pancreasSCE) ## 1. Cellular neighborhood (CN) sce <- buildSpatialGraph(pancreasSCE, img_id = "ImageNb", type = "knn", name = "knn_cn_graph", k = 5) sce <- aggregateNeighbors(sce, colPairName = "knn_cn_graph", aggregate_by = "metadata", count_by = "CellType", name = "aggregatedCellTypes") cur_cluster <- kmeans(sce$aggregatedCellTypes, centers = 3) sce$cellular_neighborhood <- factor(cur_cluster$cluster) plotSpatial(sce, img_id = "ImageNb", colPairName = "knn_cn_graph", node_color_by = "cellular_neighborhood", scales = "free") ## 2. Spatial context (SC) sce <- buildSpatialGraph(sce, img_id = "ImageNb", type = "knn", name = "knn_sc_graph", k = 15) sce <- aggregateNeighbors(sce, colPairName = "knn_sc_graph", aggregate_by = "metadata", count_by = "cellular_neighborhood", name = "aggregatedNeighborhood") # Detect spatial context sce <- detectSpatialContext(sce, entry = "aggregatedNeighborhood", threshold = 0.9) plotSpatial(sce, img_id = "ImageNb", colPairName = "knn_sc_graph", node_color_by = "spatial_context", scales = "free") # Filter spatial context # By group sce <- filterSpatialContext(sce, group_by = "ImageNb", group_threshold = 2) plotSpatial(sce, img_id = "ImageNb", colPairName = "knn_sc_graph", node_color_by = "spatial_context_filtered", scales = "free") # By cells sce <- filterSpatialContext(sce, group_by = "ImageNb", cells_threshold = 15) plotSpatial(sce, img_id = "ImageNb", colPairName = "knn_sc_graph", node_color_by = "spatial_context_filtered", scales = "free")
set.seed(22) library(cytomapper) data(pancreasSCE) ## 1. Cellular neighborhood (CN) sce <- buildSpatialGraph(pancreasSCE, img_id = "ImageNb", type = "knn", name = "knn_cn_graph", k = 5) sce <- aggregateNeighbors(sce, colPairName = "knn_cn_graph", aggregate_by = "metadata", count_by = "CellType", name = "aggregatedCellTypes") cur_cluster <- kmeans(sce$aggregatedCellTypes, centers = 3) sce$cellular_neighborhood <- factor(cur_cluster$cluster) plotSpatial(sce, img_id = "ImageNb", colPairName = "knn_cn_graph", node_color_by = "cellular_neighborhood", scales = "free") ## 2. Spatial context (SC) sce <- buildSpatialGraph(sce, img_id = "ImageNb", type = "knn", name = "knn_sc_graph", k = 15) sce <- aggregateNeighbors(sce, colPairName = "knn_sc_graph", aggregate_by = "metadata", count_by = "cellular_neighborhood", name = "aggregatedNeighborhood") # Detect spatial context sce <- detectSpatialContext(sce, entry = "aggregatedNeighborhood", threshold = 0.9) plotSpatial(sce, img_id = "ImageNb", colPairName = "knn_sc_graph", node_color_by = "spatial_context", scales = "free") # Filter spatial context # By group sce <- filterSpatialContext(sce, group_by = "ImageNb", group_threshold = 2) plotSpatial(sce, img_id = "ImageNb", colPairName = "knn_sc_graph", node_color_by = "spatial_context_filtered", scales = "free") # By cells sce <- filterSpatialContext(sce, group_by = "ImageNb", cells_threshold = 15) plotSpatial(sce, img_id = "ImageNb", colPairName = "knn_sc_graph", node_color_by = "spatial_context_filtered", scales = "free")
Detection of cells close to the image border for subsequent exclusion from downstream analyses.
findBorderCells(object, img_id, border_dist, coords = c("Pos_X", "Pos_Y"))
findBorderCells(object, img_id, border_dist, coords = c("Pos_X", "Pos_Y"))
object |
a |
img_id |
single character indicating the |
border_dist |
single numeric defining the distance to the image border. The image border here is defined as the minimum and maximum among the cells' x and y location. |
coords |
character vector of length 2 specifying the names of the
|
an object of class(object)
containing the logical
border_cells
entry in the colData
slot.
library(cytomapper) data("pancreasSCE") sce <- findBorderCells(pancreasSCE, img_id = "ImageNb", border_dist = 10) plotSpatial(sce, img_id = "ImageNb", node_color_by = "border_cells", scales = "free")
library(cytomapper) data("pancreasSCE") sce <- findBorderCells(pancreasSCE, img_id = "ImageNb", border_dist = 10) plotSpatial(sce, img_id = "ImageNb", node_color_by = "border_cells", scales = "free")
Function to return the distance of the closest cell of interest for each cell in the data. In the case of patched/clustered cells negative distances are returned by default which indicate the distance of the cells of interest to the closest cell that is not of the type of cells of interest.
minDistToCells( object, x_cells, img_id, name = "distToCells", coords = c("Pos_X", "Pos_Y"), return_neg = TRUE, BPPARAM = SerialParam() )
minDistToCells( object, x_cells, img_id, name = "distToCells", coords = c("Pos_X", "Pos_Y"), return_neg = TRUE, BPPARAM = SerialParam() )
object |
a |
x_cells |
logical vector of length equal to the number of cells
contained in |
img_id |
single character indicating the |
name |
character specifying the name of the |
coords |
character vector of length 2 specifying the names of the
|
return_neg |
logical indicating whether negative distances are to be returned for the distances of patched/spatially clustered cells. |
BPPARAM |
a |
returns an object of class(object)
containing a new column
entry to colData(object)[[name]]
. Cells in the object are grouped
by entries in img_id
.
The minDistToCells
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.
Daniel Schulz ([email protected])
library(cytomapper) data(pancreasSCE) # Build interaction graph pancreasSCE <- buildSpatialGraph(pancreasSCE, img_id = "ImageNb", type = "expansion",threshold = 20) # Detect patches of "celltype_B" cells pancreasSCE <- patchDetection(pancreasSCE, img_id = "ImageNb", patch_cells = pancreasSCE$CellType == "celltype_B", colPairName = "expansion_interaction_graph", min_patch_size = 20, expand_by = 1) plotSpatial(pancreasSCE, img_id = "ImageNb", node_color_by = "patch_id", scales = "free") # Distance to celltype_B patches pancreasSCE <- minDistToCells(pancreasSCE, x_cells = !is.na(pancreasSCE$patch_id), coords = c("Pos_X","Pos_Y"), img_id = "ImageNb") plotSpatial(pancreasSCE, img_id = "ImageNb", node_color_by = "distToCells", scales = "free")
library(cytomapper) data(pancreasSCE) # Build interaction graph pancreasSCE <- buildSpatialGraph(pancreasSCE, img_id = "ImageNb", type = "expansion",threshold = 20) # Detect patches of "celltype_B" cells pancreasSCE <- patchDetection(pancreasSCE, img_id = "ImageNb", patch_cells = pancreasSCE$CellType == "celltype_B", colPairName = "expansion_interaction_graph", min_patch_size = 20, expand_by = 1) plotSpatial(pancreasSCE, img_id = "ImageNb", node_color_by = "patch_id", scales = "free") # Distance to celltype_B patches pancreasSCE <- minDistToCells(pancreasSCE, x_cells = !is.na(pancreasSCE$patch_id), coords = c("Pos_X","Pos_Y"), img_id = "ImageNb") plotSpatial(pancreasSCE, img_id = "ImageNb", node_color_by = "distToCells", scales = "free")
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() )
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() )
object |
a |
patch_cells |
logical vector of length equal to the number of cells
contained in |
colPairName |
single character indicating the |
min_patch_size |
single integer indicating the minimum number of connected cells that make up a patch before expansion. |
name |
single character specifying the |
expand_by |
single numeric indicating in which vicinity range cells should be considered as belonging to the patch (see Details). |
coords |
character vector of length 2 specifying the names of the
|
convex |
should the convex hull be computed before expansion? Default: the concave hull is computed. |
img_id |
single character indicating the |
BPPARAM |
a |
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.
Tobias Hoch
adapted by Nils Eling ([email protected])
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) # 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") plotSpatial(pancreasSCE, img_id = "ImageNb", node_color_by = "patch_id", scales = "free")
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) # 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") plotSpatial(pancreasSCE, img_id = "ImageNb", node_color_by = "patch_id", scales = "free")
This function constructs polygons around patch cells and computes their area.
patchSize( object, patch_name = "patch_id", coords = c("Pos_X", "Pos_Y"), convex = FALSE )
patchSize( object, patch_name = "patch_id", coords = c("Pos_X", "Pos_Y"), convex = FALSE )
object |
a |
patch_name |
single character indicating the |
coords |
character vector of length 2 specifying the names of the
|
convex |
should the convex hull be computed to construct the polygon? Default: the concave hull is computed. |
A DataFrame object containing the patch identifier, the constructed polygon and the polygon size.
Nils Eling ([email protected])
library(cytomapper) data(pancreasSCE) # Build interaction graph pancreasSCE <- buildSpatialGraph(pancreasSCE, img_id = "ImageNb", type = "expansion", threshold = 20) # Detect patches of "celltype_B" cells pancreasSCE <- patchDetection(pancreasSCE, patch_cells = pancreasSCE$CellType == "celltype_B", expand_by = 5, img_id = "ImageNb", colPairName = "expansion_interaction_graph") # Compute the patch area patchSize(pancreasSCE)
library(cytomapper) data(pancreasSCE) # Build interaction graph pancreasSCE <- buildSpatialGraph(pancreasSCE, img_id = "ImageNb", type = "expansion", threshold = 20) # Detect patches of "celltype_B" cells pancreasSCE <- patchDetection(pancreasSCE, patch_cells = pancreasSCE$CellType == "celltype_B", expand_by = 5, img_id = "ImageNb", colPairName = "expansion_interaction_graph") # Compute the patch area patchSize(pancreasSCE)
A general function to plot spatial locations of cells while specifying color, shape, size. Cell-cell interactions can be visualized in form of edges between points.
plotSpatial( object, img_id, coords = c("Pos_X", "Pos_Y"), node_color_by = NULL, node_shape_by = NULL, node_size_by = NULL, node_color_fix = NULL, node_shape_fix = NULL, node_size_fix = NULL, assay_type = NULL, draw_edges = FALSE, directed = TRUE, edge_color_by = NULL, edge_width_by = NULL, edge_color_fix = NULL, edge_width_fix = NULL, arrow = NULL, end_cap = NULL, colPairName = NULL, nodes_first = TRUE, ncols = NULL, nrows = NULL, scales = "fixed", flip_x = FALSE, flip_y = TRUE, aspect_ratio = "auto" )
plotSpatial( object, img_id, coords = c("Pos_X", "Pos_Y"), node_color_by = NULL, node_shape_by = NULL, node_size_by = NULL, node_color_fix = NULL, node_shape_fix = NULL, node_size_fix = NULL, assay_type = NULL, draw_edges = FALSE, directed = TRUE, edge_color_by = NULL, edge_width_by = NULL, edge_color_fix = NULL, edge_width_fix = NULL, arrow = NULL, end_cap = NULL, colPairName = NULL, nodes_first = TRUE, ncols = NULL, nrows = NULL, scales = "fixed", flip_x = FALSE, flip_y = TRUE, aspect_ratio = "auto" )
object |
a |
img_id |
single character indicating the |
coords |
character vector of length 2 specifying the names of the
|
node_color_by |
single character indicating the |
node_shape_by |
single character indicating the |
node_size_by |
single character indicating the |
node_color_fix |
single character or numeric specifying the color of all nodes. |
node_shape_fix |
single numeric or character specifying the shape of all nodes. |
node_size_fix |
single numeric specifying the size of all nodes |
assay_type |
single character indicating the assay slot from which
to extract the expression data when |
draw_edges |
should cell-cell interactions be drawn as edges between nodes? |
directed |
should cell-cell interactions be handled as a directed graph? |
edge_color_by |
single character indicating by which to color the edges. See details for more information. |
edge_width_by |
single character determining the size of the edges. See details for more information. |
edge_color_fix |
single character or numeric specifying the color of all edges. |
edge_width_fix |
single numeric specifying the size of all edges. |
arrow |
an |
end_cap |
a |
colPairName |
single character specifying the
|
nodes_first |
should the nodes be plotted first and then the edges? |
ncols |
number of columns of the grid to arrange individual images. |
nrows |
number of rows of the grid to arrange individual images. |
scales |
one of |
flip_x |
flip the x-axis? |
flip_y |
flip the y-axis? |
aspect_ratio |
single numeric, "auto" or NULL to define the relative
ratio between the physical units of the x and y axis. If |
returns a ggplot
object.
By default, the cells' locations are visualized in form of points (here also
referred to as "nodes") on a 2-dimensional plane. The cells' coordinates are
extracted either from colData(object)
slot (for a
SingleCellExperiment
input object) or from the
spatialCoords(object)
slot (for a SpatialExperiment
input
object). Node aesthetics are controlled by setting node_color_by
,
node_shape_by
and node_size_by
for associating the aesthetics
with variables. If node aesthetics should be the same for all nodes,
node_color_fix
, node_shape_fix
and node_size_fix
can
be set.
When draw_edges = TRUE
, cell-cell interactions are visualized in form
of edges between nodes. For this, object
needs to contain
column pairings in colPair(object, colPairName)
. Edge color and size
can be set by specifying either an entry in
mcols(colPair(object, colPairName))
(edge attributes) or in
colData(object)
. In the latter case, edges are colored by attributes
associated to the "from" node. Variable aesthetics can be set using
edge_color_by
and edge_width_by
. If all edges should have
the same width or color, edge_color_fix
and edge_width_fix
can be set.
Arrows for displaying directed graphs can be drawn by supplying a
arrow
object. Arrow attributes can be set within this
class. To cap the edge before it reaches the next node, the end_cap
parameter can be used.
Nils Eling ([email protected])
buildSpatialGraph
for constructing interaction graphs
ggraph
for handling graph aesthetics
library(cytomapper) data(pancreasSCE) sce <- buildSpatialGraph(pancreasSCE, img_id = "ImageNb", type = "knn", k = 3, directed = FALSE) # Only nodes plotSpatial(sce, img_id = "ImageNb", node_color_by = "CellType", node_shape_by = "ImageNb", node_size_by = "Area", scales = "free") # With edges and nodes colored by expression plotSpatial(sce, img_id = "ImageNb", node_color_by = "PIN", assay_type = "exprs", node_shape_by = "ImageNb", node_size_by = "Area", draw_edges = TRUE, colPairName = "knn_interaction_graph", edge_color_by = "Pattern", scales = "free") # With arrows plotSpatial(sce, img_id = "ImageNb", node_color_by = "CellType", node_shape_by = "ImageNb", node_size_by = "Area", draw_edges = TRUE, colPairName = "knn_interaction_graph", edge_color_fix = "green", arrow = grid::arrow(length = grid::unit(0.1, "inch")), end_cap = ggraph::circle(0.2, "cm"), scales = "free")
library(cytomapper) data(pancreasSCE) sce <- buildSpatialGraph(pancreasSCE, img_id = "ImageNb", type = "knn", k = 3, directed = FALSE) # Only nodes plotSpatial(sce, img_id = "ImageNb", node_color_by = "CellType", node_shape_by = "ImageNb", node_size_by = "Area", scales = "free") # With edges and nodes colored by expression plotSpatial(sce, img_id = "ImageNb", node_color_by = "PIN", assay_type = "exprs", node_shape_by = "ImageNb", node_size_by = "Area", draw_edges = TRUE, colPairName = "knn_interaction_graph", edge_color_by = "Pattern", scales = "free") # With arrows plotSpatial(sce, img_id = "ImageNb", node_color_by = "CellType", node_shape_by = "ImageNb", node_size_by = "Area", draw_edges = TRUE, colPairName = "knn_interaction_graph", edge_color_fix = "green", arrow = grid::arrow(length = grid::unit(0.1, "inch")), end_cap = ggraph::circle(0.2, "cm"), scales = "free")
Function to plot directed spatial context graphs based on symbolic edge-lists and vertex metadata, which operates on the cohort-level. The user can specify node, node_label and edge aesthetics.
plotSpatialContext( object, entry = "spatial_context", group_by = "sample_id", node_color_by = NULL, node_size_by = NULL, node_color_fix = NULL, node_size_fix = NULL, node_label_repel = TRUE, node_label_color_by = NULL, node_label_color_fix = NULL, draw_edges = TRUE, edge_color_fix = NULL, return_data = FALSE )
plotSpatialContext( object, entry = "spatial_context", group_by = "sample_id", node_color_by = NULL, node_size_by = NULL, node_color_fix = NULL, node_size_fix = NULL, node_label_repel = TRUE, node_label_color_by = NULL, node_label_color_fix = NULL, draw_edges = TRUE, edge_color_fix = NULL, return_data = FALSE )
object |
a |
entry |
single character specifying the |
group_by |
a single character indicating the |
node_color_by |
single character either
|
node_size_by |
single character either |
node_color_fix |
single character specifying the color of all nodes. |
node_size_fix |
single numeric specifying the size of all nodes. |
node_label_repel |
should nodes be labelled? Defaults to TRUE. |
node_label_color_by |
single character either
|
node_label_color_fix |
single character specifying the color of all node labels. |
draw_edges |
should edges be drawn between nodes? Defaults to TRUE. |
edge_color_fix |
single character specifying the color of all edges. |
return_data |
should the edge list and vertex metadata for graph
construction be returned as a |
returns a ggplot
object or a list
of two
data.frames
.
Lasse Meyer ([email protected])
detectSpatialContext
for the function to detect
spatial contexts
filterSpatialContext
for the function to filter
spatial contexts
set.seed(22) library(cytomapper) data(pancreasSCE) ## 1. Cellular neighborhood (CN) sce <- buildSpatialGraph(pancreasSCE, img_id = "ImageNb", type = "knn", name = "knn_cn_graph", k = 5) sce <- aggregateNeighbors(sce, colPairName = "knn_cn_graph", aggregate_by = "metadata", count_by = "CellType", name = "aggregatedCellTypes") cur_cluster <- kmeans(sce$aggregatedCellTypes, centers = 3) sce$cellular_neighborhood <- factor(cur_cluster$cluster) plotSpatial(sce, img_id = "ImageNb", colPairName = "knn_cn_graph", node_color_by = "cellular_neighborhood", scales = "free") ## 2. Spatial context (SC) sce <- buildSpatialGraph(sce, img_id = "ImageNb", type = "knn", name = "knn_sc_graph", k = 15) sce <- aggregateNeighbors(sce, colPairName = "knn_sc_graph", aggregate_by = "metadata", count_by = "cellular_neighborhood", name = "aggregatedNeighborhood") # Detect spatial context sce <- detectSpatialContext(sce, entry = "aggregatedNeighborhood", threshold = 0.9) plotSpatial(sce, img_id = "ImageNb", colPairName = "knn_sc_graph", node_color_by = "spatial_context", scales = "free") # Plot spatial context - default plotSpatialContext(sce, group_by = "ImageNb") # Plot spatial context - adjust aesthetics plotSpatialContext(sce, group_by = "ImageNb", node_color_by = "name", node_size_by = "n_cells", node_label_color_by = "name") plotSpatialContext(sce, group_by = "ImageNb", node_color_by = "n_cells", node_size_by = "n_group") # Plot spatial context - return data plotSpatialContext(sce, group_by = "ImageNb", return_data = TRUE)
set.seed(22) library(cytomapper) data(pancreasSCE) ## 1. Cellular neighborhood (CN) sce <- buildSpatialGraph(pancreasSCE, img_id = "ImageNb", type = "knn", name = "knn_cn_graph", k = 5) sce <- aggregateNeighbors(sce, colPairName = "knn_cn_graph", aggregate_by = "metadata", count_by = "CellType", name = "aggregatedCellTypes") cur_cluster <- kmeans(sce$aggregatedCellTypes, centers = 3) sce$cellular_neighborhood <- factor(cur_cluster$cluster) plotSpatial(sce, img_id = "ImageNb", colPairName = "knn_cn_graph", node_color_by = "cellular_neighborhood", scales = "free") ## 2. Spatial context (SC) sce <- buildSpatialGraph(sce, img_id = "ImageNb", type = "knn", name = "knn_sc_graph", k = 15) sce <- aggregateNeighbors(sce, colPairName = "knn_sc_graph", aggregate_by = "metadata", count_by = "cellular_neighborhood", name = "aggregatedNeighborhood") # Detect spatial context sce <- detectSpatialContext(sce, entry = "aggregatedNeighborhood", threshold = 0.9) plotSpatial(sce, img_id = "ImageNb", colPairName = "knn_sc_graph", node_color_by = "spatial_context", scales = "free") # Plot spatial context - default plotSpatialContext(sce, group_by = "ImageNb") # Plot spatial context - adjust aesthetics plotSpatialContext(sce, group_by = "ImageNb", node_color_by = "name", node_size_by = "n_cells", node_label_color_by = "name") plotSpatialContext(sce, group_by = "ImageNb", node_color_by = "n_cells", node_size_by = "n_group") # Plot spatial context - return data plotSpatialContext(sce, group_by = "ImageNb", return_data = TRUE)
Helper function for estimating the spillover matrix. This function visualizes the median pixel intensities per spot (rows) and per channel (columns) in form of a heatmap.
plotSpotHeatmap( object, spot_id = "sample_id", channel_id = "channel_name", assay_type = "counts", statistic = "median", log = TRUE, threshold = NULL, order_metals = TRUE, color = viridis(100), breaks = NA, legend_breaks = NA, cluster_cols = FALSE, cluster_rows = FALSE, ... )
plotSpotHeatmap( object, spot_id = "sample_id", channel_id = "channel_name", assay_type = "counts", statistic = "median", log = TRUE, threshold = NULL, order_metals = TRUE, color = viridis(100), breaks = NA, legend_breaks = NA, cluster_cols = FALSE, cluster_rows = FALSE, ... )
object |
a |
spot_id |
character string indicating which |
channel_id |
character string indicating which |
assay_type |
character string indicating which assay to use (default
|
statistic |
the statistic to use when aggregating channels per spot
(default |
log |
should the aggregated pixel intensities be |
threshold |
single numeric indicating a threshold after pixel
aggregation. All aggregated values larger than |
order_metals |
should the metals be ordered based on spotted mass? |
color |
see parameter in |
breaks |
see parameter in |
legend_breaks |
see parameter in |
cluster_cols |
see parameter in |
cluster_rows |
see parameter in |
... |
other arguments passed to |
a pheatmap
object
Visualizing the aggregated pixel intensities serves two purposes:
Small median pixel intensities (< 200 counts) might hinder the robust
estimation of the channel spillover. In that case, consecutive pixels can be
summed (see binAcrossPixels
).
Each spotted metal (row) should show the highest median pixel intensity in its corresponding channel (column). If this is not the case, either the naming of the .txt files was incorrect or the incorrect metal was spotted.
By setting the threshold
parameter, the user can easily identify spots
where pixel intensities are too low for robust spillover estimation.
Nils Eling ([email protected])
pheatmap
for visual modifications
aggregateAcrossCells
for the aggregation
function
path <- system.file("extdata/spillover", package = "imcRtools") # Read in .txt files sce <- readSCEfromTXT(path) # Visualizes heatmap plotSpotHeatmap(sce) # Visualizes thresholding results plotSpotHeatmap(sce, log = FALSE, threshold = 200)
path <- system.file("extdata/spillover", package = "imcRtools") # Read in .txt files sce <- readSCEfromTXT(path) # Visualizes heatmap plotSpotHeatmap(sce) # Visualizes thresholding results plotSpotHeatmap(sce, log = FALSE, threshold = 200)
Reader function to generate a
SpatialExperiment
or
SingleCellExperiment
object from single-cell data
obtained by the
ImcSegmentationPipeline
pipeline.
read_cpout( path, object_file = "cell.csv", image_file = "Image.csv", panel_file = "panel.csv", graph_file = "Object relationships.csv", object_feature_file = "var_cell.csv", intensities = "Intensity_MeanIntensity_FullStack", extract_imgid_from = "ImageNumber", extract_cellid_from = "ObjectNumber", extract_coords_from = c("Location_Center_X", "Location_Center_Y"), extract_cellmetadata_from = c("AreaShape_Area", "Neighbors_NumberOfNeighbors_8", "AreaShape_Eccentricity", "AreaShape_MajorAxisLength", "AreaShape_MinorAxisLength", "AreaShape_MeanRadius"), extract_imagemetadata_from = c("Metadata_acname", "Metadata_acid", "Metadata_description"), extract_graphimageid_from = "First Image Number", extract_graphcellids_from = c("First Object Number", "Second Object Number"), extract_metal_from = "Metal Tag", scale_intensities = TRUE, extract_scalingfactor_from = "Scaling_FullStack", return_as = c("spe", "sce") )
read_cpout( path, object_file = "cell.csv", image_file = "Image.csv", panel_file = "panel.csv", graph_file = "Object relationships.csv", object_feature_file = "var_cell.csv", intensities = "Intensity_MeanIntensity_FullStack", extract_imgid_from = "ImageNumber", extract_cellid_from = "ObjectNumber", extract_coords_from = c("Location_Center_X", "Location_Center_Y"), extract_cellmetadata_from = c("AreaShape_Area", "Neighbors_NumberOfNeighbors_8", "AreaShape_Eccentricity", "AreaShape_MajorAxisLength", "AreaShape_MinorAxisLength", "AreaShape_MeanRadius"), extract_imagemetadata_from = c("Metadata_acname", "Metadata_acid", "Metadata_description"), extract_graphimageid_from = "First Image Number", extract_graphcellids_from = c("First Object Number", "Second Object Number"), extract_metal_from = "Metal Tag", scale_intensities = TRUE, extract_scalingfactor_from = "Scaling_FullStack", return_as = c("spe", "sce") )
path |
full path to the CellProfiler output folder. |
object_file |
single character indicating the file name storing the object/cell-specific intensities and metadata. |
image_file |
single character indicating the file name storing meta data
per image (can be |
panel_file |
single character indicating the file name storing the panel
information (can be |
graph_file |
single character indicating the file name storing the
object/cell interaction information (can be |
object_feature_file |
single character indicating the file name storing object/cell feature information. |
intensities |
single character indicating which column entries of the
|
extract_imgid_from |
single character indicating which column entries of
the |
extract_cellid_from |
single character indicating which column entry
of the |
extract_coords_from |
character vector indicating which column entries
of the |
extract_cellmetadata_from |
character vector indicating which additional
object/cell specific metadata to extract from the |
extract_imagemetadata_from |
character vector indicating which
additional image specific metadata to extract from the |
extract_graphimageid_from |
single character indicating which column
entries of the |
extract_graphcellids_from |
character vector indicating which column
entries of the |
extract_metal_from |
single character indicating which column entry of
the |
scale_intensities |
single logical. Should the measured intensity
features be scaled by |
extract_scalingfactor_from |
single character indicating which column
entries of the |
return_as |
should the object be returned as
|
returns a SpatialExperiment
or SingleCellExperiment
object with markers in rows and cells in columns.
In the case of both containers x
, intensity features (as selected
by the intensities
parameter) are stored in the counts(x)
slot.
Cell metadata (e.g morphological features) are stored in the
colData(x)
slot. The interaction graphs are stored as
SelfHits
object in the
colPair(x, "neighborhood")
slot.
Intensity features are extracted via partial string matching. Internally, the
read_cpout
function checks if per channel a single intensity feature
is read in (by checking the _cXY
ending where XY
is the
channel number).
In the case of a returned SpatialExperiment
object, the cell
coordinates are stored in the spatialCoords(x)
slot.
In the case of a returned SingleCellExperiment
object, the cell
coordinates are stored in the colData(x)
slot named as Pos_X
and Pos_Y
.
Tobias Hoch
Nils Eling ([email protected])
https://github.com/BodenmillerGroup/ImcSegmentationPipeline for the pipeline
read_steinbock
for reading in single-cell data as produced by
the steinbock pipeline
colPair
for information on how to work
with the cell-cell interaction graphs
path <- system.file("extdata/mockData/cpout", package = "imcRtools") # Read in as SpatialExperiment object x <- read_cpout(path, graph_file = "Object_relationships.csv") x # Read in as SingleCellExperiment object x <- read_cpout(path, graph_file = "Object_relationships.csv", return_as = "sce") x
path <- system.file("extdata/mockData/cpout", package = "imcRtools") # Read in as SpatialExperiment object x <- read_cpout(path, graph_file = "Object_relationships.csv") x # Read in as SingleCellExperiment object x <- read_cpout(path, graph_file = "Object_relationships.csv", return_as = "sce") x
Reader function to generate a
SpatialExperiment
or
SingleCellExperiment
object from single-cell data
obtained by the
steinbock
pipeline.
read_steinbock( path, intensities_folder = "intensities", regionprops_folder = "regionprops", graphs_folder = "neighbors", pattern = NULL, extract_cellid_from = "Object", extract_coords_from = c("centroid-1", "centroid-0"), image_file = "images.csv", extract_imagemetadata_from = c("width_px", "height_px"), panel_file = "panel.csv", extract_names_from = "name", return_as = c("spe", "sce"), BPPARAM = SerialParam() )
read_steinbock( path, intensities_folder = "intensities", regionprops_folder = "regionprops", graphs_folder = "neighbors", pattern = NULL, extract_cellid_from = "Object", extract_coords_from = c("centroid-1", "centroid-0"), image_file = "images.csv", extract_imagemetadata_from = c("width_px", "height_px"), panel_file = "panel.csv", extract_names_from = "name", return_as = c("spe", "sce"), BPPARAM = SerialParam() )
path |
full path to the steinbock output folder |
intensities_folder |
name of the folder containing the intensity measurements per image |
regionprops_folder |
name of the folder containing the cell-specific
morphology and spatial measurements per image. Can be set to |
graphs_folder |
name of the folder containing the spatial connectivity
graphs per image. Can be set to |
pattern |
regular expression specifying a subset of files that should be read in. |
extract_cellid_from |
single character indicating which column entry in the intensity files contains the integer cell id. |
extract_coords_from |
character vector indicating which column entries in the regionprops files contain the x (first entry) and y (second entry) coordinates. |
image_file |
single character indicating the file name storing meta data
per image (can be |
extract_imagemetadata_from |
character vector indicating which
additional image specific metadata to extract from the |
panel_file |
single character containing the name of the panel file. This can either be inside the steinbock path (recommended) or located somewhere else. |
extract_names_from |
single character indicating the column of the panel file containing the channel names. |
return_as |
should the object be returned as
|
BPPARAM |
parameters for parallelised processing. |
returns a SpatialExperiment
or SingleCellExperiment
object markers in rows and cells in columns.
In the case of both containers x
, intensity features are stored in
the counts(x)
slot. Morphological features are stored in the
colData(x)
slot. The graphs are stored as
SelfHits
object in the
colPair(x, "neighborhood")
slot.
In the case of a returned SpatialExperiment
object, the cell
coordinates are stored in the spatialCoords(x)
slot.
In the case of a returned SingleCellExperiment
object, the cell
coordinates are stored in the colData(x)
slot named as Pos_X
and Pos_Y
.
Nils Eling ([email protected])
https://github.com/BodenmillerGroup/steinbock for the pipeline
read_cpout
for reading in single-cell data as produced by the
ImcSegmentationPipeline
SingleCellExperiment
and
SpatialExperiment
for the constructor
functions.
colPair
for information on how to work
with the cell-cell interaction graphs
bpparam
for the parallelised backend
path <- system.file("extdata/mockData/steinbock", package = "imcRtools") # Read in as SpatialExperiment object x <- read_steinbock(path) x # Read in as SingleCellExperiment object x <- read_steinbock(path, return_as = "sce") x # Read in a subset of files x <- read_steinbock(path, pattern = "mockData1") x # Only read in intensities x <- read_steinbock(path, graphs_folder = NULL, regionprops_folder = NULL) x # Parallelisation x <- read_steinbock(path, BPPARAM = BiocParallel::bpparam())
path <- system.file("extdata/mockData/steinbock", package = "imcRtools") # Read in as SpatialExperiment object x <- read_steinbock(path) x # Read in as SingleCellExperiment object x <- read_steinbock(path, return_as = "sce") x # Read in a subset of files x <- read_steinbock(path, pattern = "mockData1") x # Only read in intensities x <- read_steinbock(path, graphs_folder = NULL, regionprops_folder = NULL) x # Parallelisation x <- read_steinbock(path, BPPARAM = BiocParallel::bpparam())
Reader function to generate Image
objects
in form of a CytoImageList
container from .txt files.
readImagefromTXT( path, pattern = ".txt$", channel_pattern = "[A-Za-z]{1,2}[0-9]{2,3}Di", index_names = c("X", "Y"), BPPARAM = SerialParam() )
readImagefromTXT( path, pattern = ".txt$", channel_pattern = "[A-Za-z]{1,2}[0-9]{2,3}Di", index_names = c("X", "Y"), BPPARAM = SerialParam() )
path |
Full path to where the individual .txt files are located. This is usualy the path where the .mcd file is located. |
pattern |
pattern to select which files should be read in (default
|
channel_pattern |
regular expression to select the channel names from the files. |
index_names |
exact names of the columns storing the x and y coordinates of the image |
BPPARAM |
parameters for parallelized reading in of images. This is only recommended for very large images. |
returns a CytoImageList
object containing one
Image
object per .txt file.
As part of the raw data folder, the Hyperion imaging system writes out one .txt file per acquisition. These files store the ion counts per pixel and channel.
This function reads these .txt files into a single
CytoImageList
object for downstream analysis. The
pattern
argument allows selection of all .txt files or a specific
subset of files. The channelNames
of the
CytoImageList
object are determined by the channel_pattern
argument.
Nils Eling ([email protected])
CytoImageList
for the container
MulticoreParam
for parallelized processing
Image
for the multi-channel image object
vignette("cytomapper")
for visualization of multi-channel images
path <- system.file("extdata/mockData/raw", package = "imcRtools") # Read in all images x <- readImagefromTXT(path) x # Read in specific files y <- readImagefromTXT(path, pattern = "ROI_002") y # Read in other channelNames z <- readImagefromTXT(path, channel_pattern = "[A-Za-z]{2}[0-9]{3}") z
path <- system.file("extdata/mockData/raw", package = "imcRtools") # Read in all images x <- readImagefromTXT(path) x # Read in specific files y <- readImagefromTXT(path, pattern = "ROI_002") y # Read in other channelNames z <- readImagefromTXT(path, channel_pattern = "[A-Za-z]{2}[0-9]{3}") z
Helper function to process raw .txt files acquired by the
Hyperion imaging system into a SingleCellExperiment
object. This function is mainly used to read-in data generated from a
"spillover slide". Here, each .txt file contains the measurements of
multiple pixels for a single stain across all open channels.
readSCEfromTXT( x, pattern = ".txt$", metadata_cols = c("Start_push", "End_push", "Pushes_duration", "X", "Y", "Z"), verbose = TRUE, read_metal_from_filename = TRUE )
readSCEfromTXT( x, pattern = ".txt$", metadata_cols = c("Start_push", "End_push", "Pushes_duration", "X", "Y", "Z"), verbose = TRUE, read_metal_from_filename = TRUE )
x |
input can be of different types:
|
pattern |
pattern to select which files should be read in (default
|
metadata_cols |
character vector indicating which column entries of the
.txt files should be saved in the |
verbose |
logical indicating if additional information regarding the spotted and acquired masses should be shown. |
read_metal_from_filename |
should the sample metal and mass be extracted from the file/object names? |
returns a SCE object where pixels are stored as columns and acquired channels are stored as rows.
As described in the original publication, single metal spots are acquired using the Hyperion imaging system. Each acquisition corresponds to one spot. All acquisitions are stored in a single .mcd file and individual acquisitions are stored in single .txt files.
This function aggregates these measurements into a single
SingleCellExperiment
object. For this, two inputs are possible:
x
is a path:
By default all .txt files are read in from the specified path. Here, the path
should indicate the location of the spillover slide measurement. The file
names of the .txt file must contain the spotted metal isotope name in the
format (mt)(mass)
(e.g. Sm152
for Samarium isotope with the
atomic mass 152). Internally, the last occurrence of such a pattern is read
in as the metal isotope name and stored in the colData(sce)$sample_id
slot.
x
is a named list:
If there are issues with reading in the metal isotope names from the .txt
file names, the user can provide a list for which each entry contains the
contents of a single .txt file. The names of the list must indicate the
spotted metal in the format (mt)(mass)
. These names will be stored in
the colData(sce)$sample_id
slot.
When read_metal_from_filename = FALSE
, the function will not attempt
to read in the spotted metal isotopes from the file or list names. Therefore,
only the sample_id
will be set based on the file/list names.
Nils Eling ([email protected])
# Read files from path path <- system.file("extdata/spillover", package = "imcRtools") sce <- readSCEfromTXT(path) sce # Read files as list cur_file_names <- list.files(path, pattern = ".txt", full.names = TRUE) cur_files <- lapply(cur_file_names, read.delim) names(cur_files) <- sub(".txt", "", basename(cur_file_names)) sce <- readSCEfromTXT(cur_files) sce
# Read files from path path <- system.file("extdata/spillover", package = "imcRtools") sce <- readSCEfromTXT(path) sce # Read files as list cur_file_names <- list.files(path, pattern = ".txt", full.names = TRUE) cur_files <- lapply(cur_file_names, read.delim) names(cur_files) <- sub(".txt", "", basename(cur_file_names)) sce <- readSCEfromTXT(cur_files) sce
Searchable datatable object of cell and image features as extracted by CellProfiler.
show_cpout_features( path, display = c("cell_features", "image_features"), cell_features = "var_cell.csv", image_features = "var_Image.csv" )
show_cpout_features( path, display = c("cell_features", "image_features"), cell_features = "var_cell.csv", image_features = "var_Image.csv" )
path |
full path to the CellProfiler output folder |
display |
single character indicating which features to display.
Accepted entries are |
cell_features |
single character indicating the name of the file storing the extracted cell features. |
image_features |
single character indicating the name of the file storing the extracted image features. |
a datatable
object
Nils Eling ([email protected])
read_cpout
for the CellProfiler reader function
path <- system.file("extdata/mockData/cpout", package = "imcRtools") # Display cell features show_cpout_features(path) # Display image features show_cpout_features(path, display = "image_features")
path <- system.file("extdata/mockData/cpout", package = "imcRtools") # Display cell features show_cpout_features(path) # Display image features show_cpout_features(path, display = "image_features")
Cell-cell interactions are summarized in different ways and the resulting count is compared to a distribution of counts arising from random permutations.
testInteractions( object, group_by, label, colPairName, method = c("classic", "histocat", "patch"), patch_size = NULL, iter = 1000, p_threshold = 0.01, return_samples = FALSE, tolerance = sqrt(.Machine$double.eps), BPPARAM = SerialParam() )
testInteractions( object, group_by, label, colPairName, method = c("classic", "histocat", "patch"), patch_size = NULL, iter = 1000, p_threshold = 0.01, return_samples = FALSE, tolerance = sqrt(.Machine$double.eps), BPPARAM = SerialParam() )
object |
a |
group_by |
a single character indicating the |
label |
single character specifying the |
colPairName |
single character indicating the |
method |
which cell-cell interaction counting method to use (see details) |
patch_size |
if |
iter |
single numeric specifying the number of permutations to perform |
p_threshold |
single numeric indicating the empirical p-value threshold at which interactions are considered to be significantly enriched or depleted per group. |
return_samples |
single logical indicating if the permuted interaction counts of all iterations should be returned. |
tolerance |
single numeric larger than 0. This parameter defines the
difference between the permuted count and the actual counts at which both
are regarded as equal. Default taken from |
BPPARAM |
parameters for parallelized processing. |
a DataFrame containing one row per group_by
entry and unique
label
entry combination (from_label
, to_label
). The
object contains following entries:
ct
: stores the interaction count as described in the details
p_gt
: stores the fraction of perturbations equal or greater
than ct
p_lt
: stores the fraction of perturbations equal or less than
ct
interaction
: is there the tendency for a positive interaction
(attraction) between from_label
and to_label
? Is p_lt
greater than p_gt
?
p
: the smaller value of p_gt
and p_lt
.
sig
: is p
smaller than p_threshold
?
sigval
: Combination of interaction
and sig
.
-1: interaction == FALSE
and sig == TRUE
0: sig == FALSE
1: interaction == TRUE
and sig == TRUE
NA
is returned if a certain label is not present in this grouping
level.
In principle, the countInteractions
function counts the number
of edges (interactions) between each set of unique entries in
colData(object)[[label]]
. Simplified, it counts for each cell of
type A the number of neighbors of type B.
This count is averaged within each unique entry
colData(object)[[group_by]]
in three different ways:
1. method = "classic"
: The count is divided by the total number of
cells of type A. The final count can be interpreted as "How many neighbors
of type B does a cell of type A have on average?"
2. method = "histocat"
: The count is divided by the number of cells
of type A that have at least one neighbor of type B. The final count can be
interpreted as "How many many neighbors of type B has a cell of type A on
average, given it has at least one neighbor of type B?"
3. method = "patch"
: For each cell, the count is binarized to 0
(less than patch_size
neighbors of type B) or 1 (more or equal to
patch_size
neighbors of type B). The binarized counts are averaged
across all cells of type A. The final count can be interpreted as "What
fraction of cells of type A have at least a given number of neighbors of
type B?"
Within each unique entry to
colData(object)[[group_by]]
, the entries of
colData(object)[[label]]
are randomized iter
times. For each
iteration, the interactions are counted as described above. The result is a
distribution of the interaction count under spatial randomness. The
observed interaction count is compared against this Null distribution to
derive empirical p-values:
p_gt
: fraction of perturbations equal or greater than the observed
count
p_lt
: fraction of perturbations equal or less than the observed count
Based on these empirical p-values, the interaction
score (attraction
or avoidance), overall p
value and significance by comparison to
p_treshold
(sig
and sigval
) are derived.
Vito Zanotelli
Jana Fischer
adapted by Nils Eling ([email protected])
countInteractions
for counting (but not testing) cell-cell
interactions per grouping level.
bpparam
for the parallelised backend
library(cytomapper) library(BiocParallel) data(pancreasSCE) pancreasSCE <- buildSpatialGraph(pancreasSCE, img_id = "ImageNb", type = "knn", k = 3) # Classic style calculation - setting the seed inside SerialParam for reproducibility (out <- testInteractions(pancreasSCE, group_by = "ImageNb", label = "CellType", method = "classic", colPairName = "knn_interaction_graph", iter = 1000, BPPARAM = SerialParam(RNGseed = 123))) # Histocat style calculation (out <- testInteractions(pancreasSCE, group_by = "ImageNb", label = "CellType", method = "histocat", colPairName = "knn_interaction_graph", iter = 1000, BPPARAM = SerialParam(RNGseed = 123))) # Patch style calculation (out <- testInteractions(pancreasSCE, group_by = "ImageNb", label = "CellType", method = "patch", patch_size = 3, colPairName = "knn_interaction_graph", iter = 1000, BPPARAM = SerialParam(RNGseed = 123)))
library(cytomapper) library(BiocParallel) data(pancreasSCE) pancreasSCE <- buildSpatialGraph(pancreasSCE, img_id = "ImageNb", type = "knn", k = 3) # Classic style calculation - setting the seed inside SerialParam for reproducibility (out <- testInteractions(pancreasSCE, group_by = "ImageNb", label = "CellType", method = "classic", colPairName = "knn_interaction_graph", iter = 1000, BPPARAM = SerialParam(RNGseed = 123))) # Histocat style calculation (out <- testInteractions(pancreasSCE, group_by = "ImageNb", label = "CellType", method = "histocat", colPairName = "knn_interaction_graph", iter = 1000, BPPARAM = SerialParam(RNGseed = 123))) # Patch style calculation (out <- testInteractions(pancreasSCE, group_by = "ImageNb", label = "CellType", method = "patch", patch_size = 3, colPairName = "knn_interaction_graph", iter = 1000, BPPARAM = SerialParam(RNGseed = 123)))