Skip to contents

Cluster the correlograms to find patterns in length scales of spatial autocorrelation. All the correlograms clustered must be computed with the same method and have the same number of lags. Correlograms are clustered jointly across samples.

Usage

clusterCorrelograms(
  sfe,
  features,
  BLUSPARAM,
  sample_id = "all",
  method = "I",
  colGeometryName = NULL,
  annotGeometryName = NULL,
  reducedDimName = NULL,
  swap_rownames = NULL,
  name = "sp.correlogram"
)

Arguments

sfe

A SpatialFeatureExperiment object with correlograms computed for features of interest.

features

Features whose correlograms to cluster.

BLUSPARAM

A BlusterParam object specifying the algorithm to use.

sample_id

Sample(s) in the SFE object whose cells/spots to use. Can be "all" to compute metric for all samples; the metric is computed separately for each sample.

method

"corr" for correlation, "I" for Moran's I, "C" for Geary's C

colGeometryName

Name of colGeometry from which to look for features.

annotGeometryName

Name of annotGeometry from which to look for features.

reducedDimName

Name of a dimension reduction, can be seen in reducedDimNames. colGeometryName and annotGeometryName have precedence over reducedDimName.

swap_rownames

Column name of rowData(object) to be used to identify features instead of rownames(object) when labeling plot elements. If not found in rowData, then rownames of the gene count matrix will be used.

name

Name under which the correlogram results are stored, which is by default "sp.correlogram".

Value

A data frame with 3 columns: feature for the features, cluster a factor for cluster membership of the features within each sample, and sample_id for the sample.

Examples

library(SpatialFeatureExperiment)
library(SFEData)
library(bluster)
sfe <- McKellarMuscleData("small")
#> see ?SFEData and browseVignettes('SFEData') for documentation
#> loading from cache
colGraph(sfe, "visium") <- findVisiumGraph(sfe)
inds <- c(1, 3, 4, 5)
sfe <- runUnivariate(sfe,
    type = "sp.correlogram",
    features = rownames(sfe)[inds],
    exprs_values = "counts", order = 5
)
clust <- clusterCorrelograms(sfe,
    features = rownames(sfe)[inds],
    BLUSPARAM = KmeansParam(2)
)