R/aggregateNeighbors.R
aggregateNeighbors.Rd
Function to summarize categorical or expression values of all neighbors of each cell.
a SingleCellExperiment
or SpatialExperiment
object
single character indicating the colPair(object)
entry containing the neighbor information.
character specifying whether the neighborhood should be
summarized by cellular features stored in colData(object)
(aggregate_by = "metdata"
) or by marker expression of the
neighboring cells (aggregate_by = "expression"
).
for summarize_by = "metadata"
, a single character
specifying the colData(object)
entry containing the cellular
metadata that should be summarized across each cell's neighborhood.
single logical indicating whether aggregated metadata should be returned in form of proportions instead of absolute counts.
for summarize_by = "expression"
, single character
indicating the assay slot to use.
for summarize_by = "expression"
, an integer, logical
or character vector specifying the features to use. If NULL, defaults to
all features.
for summarize_by = "expression"
, a single character
specifying the statistic to be used for summarizing the expression values
across all neighboring cells. Supported entries are "mean", "median", "sd",
"var". Defaults to "mean" if not specified.
single character specifying the name of the data frame to be
saved in the colData(object)
. Defaults to "aggregatedNeighbors" when
summarize_by = "metadata"
or "statistic_aggregatedExpression" when
summarize_by = "expression"
.
returns an object of class(object)
containing the aggregated
values in form of a DataFrame
object in
colData(object)[[name]]
.
library(cytomapper)
#> Loading required package: EBImage
#>
#> Attaching package: ‘EBImage’
#> The following object is masked from ‘package:SummarizedExperiment’:
#>
#> resize
#> The following object is masked from ‘package:Biobase’:
#>
#> channel
#> The following objects are masked from ‘package:GenomicRanges’:
#>
#> resize, tile
#> The following objects are masked from ‘package:IRanges’:
#>
#> resize, tile
#>
#> Attaching package: ‘cytomapper’
#> The following objects are masked from ‘package:Biobase’:
#>
#> channelNames, channelNames<-
data(pancreasSCE)
sce <- buildSpatialGraph(pancreasSCE, img_id = "ImageNb",
type = "knn", k = 3)
#> The returned object is ordered by the 'ImageNb' entry.
# Aggregating neighboring cell-types
sce <- aggregateNeighbors(sce, colPairName = "knn_interaction_graph",
aggregate_by = "metadata",
count_by = "CellType")
sce$aggregatedNeighbors
#> DataFrame with 362 rows and 3 columns
#> celltype_A celltype_B celltype_C
#> <numeric> <numeric> <numeric>
#> 1 0 0 1
#> 2 0 0 1
#> 3 0 0 1
#> 4 0 0 1
#> 5 0 0 1
#> ... ... ... ...
#> 358 0 0.666667 0.333333
#> 359 0 1.000000 0.000000
#> 360 0 0.000000 1.000000
#> 361 0 0.333333 0.666667
#> 362 0 0.000000 1.000000
# Aggregating neighboring expression values
sce <- aggregateNeighbors(sce, colPairName = "knn_interaction_graph",
aggregate_by = "expression",
assay_type = "exprs",
statistic = "mean")
sce$mean_aggregatedExpression
#> DataFrame with 362 rows and 5 columns
#> H3 CD99 PIN CD8a CDH
#> <numeric> <numeric> <numeric> <numeric> <numeric>
#> 1 2.13779 0.926813 0.101286 1.047837 2.11764
#> 2 2.46185 0.978440 0.228181 0.227142 2.31938
#> 3 3.17958 0.720834 0.153769 0.272481 1.41583
#> 4 2.61744 0.704024 0.134673 0.183676 2.37763
#> 5 2.07004 1.195031 0.273259 0.163717 3.00868
#> ... ... ... ... ... ...
#> 358 3.33523 2.464446 0.2588993 0.337226 2.23754
#> 359 3.15450 2.516081 0.2071171 0.294574 1.55607
#> 360 3.44273 0.622654 0.1657428 0.126924 2.69050
#> 361 3.03402 1.932238 0.1586025 0.192782 3.15072
#> 362 2.83251 0.822080 0.0583222 0.329829 1.54897