This function detects local outliers in spatial transcriptomics data based on standard quality control metrics, such as library size, unique genes, and mitochondrial ratio. Local outliers are defined as spots with low/high quality metrics compared to their surrounding neighbors, based on a modified z-score statistic.
Usage
localOutliers(
spe,
metric = "detected",
direction = "lower",
n_neighbors = 36,
samples = "sample_id",
log = TRUE,
cutoff = 3,
workers = 1
)
Arguments
- spe
SpatialExperiment or SingleCellExperiment object
- metric
colData QC metric to use for outlier detection
- direction
Direction of outlier detection (higher, lower, or both)
- n_neighbors
Number of nearest neighbors to use for outlier detection
- samples
Column name in colData to use for sample IDs
- log
Logical indicating whether to log1p transform the features (default is TRUE)
- cutoff
Cutoff for outlier detection (default is 3)
- workers
Number of workers for parallel processing (default is 1)
Examples
library(SpotSweeper)
library(SpatialExperiment)
# load example data
spe <- STexampleData::Visium_humanDLPFC()
#> see ?STexampleData and browseVignettes('STexampleData') for documentation
#> loading from cache
# change from gene id to gene names
rownames(spe) <- rowData(spe)$gene_name
# drop out-of-tissue spots
spe <- spe[, spe$in_tissue == 1]
spe <- spe[, !is.na(spe$ground_truth)]
# Identifying the mitochondrial transcripts in our SpatialExperiment.
is.mito <- rownames(spe)[grepl("^MT-", rownames(spe))]
# Calculating QC metrics for each spot using scuttle
spe <- scuttle::addPerCellQCMetrics(spe, subsets = list(Mito = is.mito))
colnames(colData(spe))
#> [1] "barcode_id" "sample_id" "in_tissue"
#> [4] "array_row" "array_col" "ground_truth"
#> [7] "reference" "cell_count" "sum"
#> [10] "detected" "subsets_Mito_sum" "subsets_Mito_detected"
#> [13] "subsets_Mito_percent" "total"
# Identifying local outliers using SpotSweeper
spe <- localOutliers(spe,
metric = "sum",
direction = "lower",
log = TRUE
)