Getting started with 'SpotSweeper'
Michael Totty
Johns Hopkins Bloomberg School of Public Health, Baltimore, MD, USABoyi Guo
Johns Hopkins Bloomberg School of Public Health, Baltimore, MD, USAStephanie Hicks
Johns Hopkins Bloomberg School of Public Health, Baltimore, MD, USA2025-10-22
Source:vignettes/getting_started.Rmd
      getting_started.RmdIntroduction
SpotSweeper is an R package for spatial transcriptomics
data quality control (QC). It provides functions for detecting and
visualizing spot-level local outliers and artifacts using
spatially-aware methods. The package is designed to work with SpatialExperiment
objects, and is compatible with data from 10X Genomics Visium and other
spatial transcriptomics platforms.
Installation
Currently, the only way to install SpotSweeper is by
downloading the development version which can be installed from GitHub using the
following:
if (!require("devtools")) install.packages("devtools")
remotes::install_github("MicTott/SpotSweeper")Once accepted in Bioconductor,
SpotSweeper will be installable using:
if (!requireNamespace("BiocManager", quietly = TRUE)) {
    install.packages("BiocManager")
}
BiocManager::install("SpotSweeper")Spot-level local outlier detection
Loading example data
Here we’ll walk you through the standard workflow for using
‘SpotSweeper’ to detect and visualize local outliers in spatial
transcriptomics data. We’ll use the Visium_humanDLPFC
dataset from the STexampleData package, which is a
SpatialExperiment object.
Because local outliers will be saved in the colData of
the SpatialExperiment object, we’ll first view the
colData and drop out-of-tissue spots before calculating
quality control (QC) metrics and running SpotSweeper.
library(SpotSweeper)
# load  Maynard et al DLPFC daatset
spe <- STexampleData::Visium_humanDLPFC()## see ?STexampleData and browseVignettes('STexampleData') for documentation## loading from cache
# show column data before SpotSweeper
colnames(colData(spe))## [1] "barcode_id"   "sample_id"    "in_tissue"    "array_row"    "array_col"   
## [6] "ground_truth" "reference"    "cell_count"
# drop out-of-tissue spots
spe <- spe[, spe$in_tissue == 1]Calculating QC metrics using scuttle
We’ll use the scuttle package to calculate QC metrics.
To do this, we’ll need to first change the rownames from
gene id to gene names. We’ll then get the mitochondrial transcripts and
calculate QC metrics for each spot using
scuttle::addPerCellQCMetrics.
# change from gene id to gene names
rownames(spe) <- rowData(spe)$gene_name
# identifying the mitochondrial transcripts
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
We can now use SpotSweeper to identify local outliers in
the spatial transcriptomics data. We’ll use the
localOutliers function to detect local outliers based on
the unique detected genes, total library size, and percent of the total
reads that are mitochondrial. These methods assume a normal
distribution, so we’ll use the log-transformed sum of the counts and the
log-transformed number of detected genes. For mitochondrial percent,
we’ll use the raw mitochondrial percentage.
# library size
spe <- localOutliers(spe,
    metric = "sum",
    direction = "lower",
    log = TRUE
)
# unique genes
spe <- localOutliers(spe,
    metric = "detected",
    direction = "lower",
    log = TRUE
)
# mitochondrial percent
spe <- localOutliers(spe,
    metric = "subsets_Mito_percent",
    direction = "higher",
    log = FALSE
)The localOutlier function automatically outputs the
results to the colData with the naming convention
X_outliers, where X is the name of the input
colData. We can then combine all outliers into a single
column called local_outliers in the colData of
the SpatialExperiment object.
# combine all outliers into "local_outliers" column
spe$local_outliers <- as.logical(spe$sum_outliers) |
    as.logical(spe$detected_outliers) |
    as.logical(spe$subsets_Mito_percent_outliers)Visualizing local outliers
We can visualize the local outliers using the
plotQCmetrics function. This function creates a scatter
plot of the specified metric and highlights the local outliers in red
using the escheR package. Here, we’ll visualize local
outliers of library size, unique genes, mitochondrial percent, and
finally, all local outliers. We’ll then arrange these plots in a grid
using ggpubr::arrange.
## Loading required package: ggplot2
# all local outliers
plotQCmetrics(spe, metric = "sum_log", outliers = "local_outliers", point_size = 1.1, 
       stroke = 0.75) +
      ggtitle("All Local Outliers")
Removing technical artifacts using SpotSweeper
Loading example data
# load in DLPFC sample with hangnail artifact
data(DLPFC_artifact)
spe <- DLPFC_artifact
# inspect colData before artifact detection
colnames(colData(spe))##  [1] "sample_id"          "in_tissue"          "array_row"         
##  [4] "array_col"          "key"                "sum_umi"           
##  [7] "sum_gene"           "expr_chrM"          "expr_chrM_ratio"   
## [10] "ManualAnnotation"   "subject"            "region"            
## [13] "sex"                "age"                "diagnosis"         
## [16] "sample_id_complete" "count"              "sizeFactor"Visualizing technical artifacts
Technical artifacts can commonly be visualized by standard QC
metrics, including library size, unique genes, or mitochondrial
percentage. We can first visualize the technical artifacts using the
plotQCmetrics function. In this sample, we can clearly see
a hangnail artifact on the right side of the tissue section in the
mitochondrial ratio plot.
plotQCmetrics(spe,
    metric = "expr_chrM_ratio",
    outliers = NULL, point_size = 1.1
) +
    ggtitle("Mitochondrial Percent")
Identifying artifacts using SpotSweeper
We can then use the findArtifacts function to identify
artifacts in the spatial transcriptomics (data. This function identifies
technical artifacts based on the first principle component of the local
variance of the specified QC metric (mito_percent) at
numerous neighorhood sizes (n_order=5). Currently,
kmeans clustering is used to cluster the technical artifact
vs high-quality Visium spots. Similar to localOutliers, the
findArtifacts function then outputs the results to the
colData.
# find artifacts using SpotSweeper
spe <- findArtifacts(spe,
    mito_percent = "expr_chrM_ratio",
    mito_sum = "expr_chrM",
    n_order = 5,
    name = "artifact"
)
# check that "artifact" is now in colData
colnames(colData(spe))##  [1] "sample_id"           "in_tissue"           "array_row"          
##  [4] "array_col"           "key"                 "sum_umi"            
##  [7] "sum_gene"            "expr_chrM"           "expr_chrM_ratio"    
## [10] "ManualAnnotation"    "subject"             "region"             
## [13] "sex"                 "age"                 "diagnosis"          
## [16] "sample_id_complete"  "count"               "sizeFactor"         
## [19] "expr_chrM_ratio_log" "coords"              "k6"                 
## [22] "k18"                 "k36"                 "k60"                
## [25] "k90"                 "artifact"Visualizing artifacts
We can visualize the artifacts using the escheR package.
Here, we’ll visualize the artifacts using the plotQCmetrics
function and arrange these plots using ggpubr::arrange.
plotQCmetrics(spe,
    metric = "expr_chrM_ratio",
    outliers = "artifact", point_size = 1.1
) +
    ggtitle("Hangnail artifact") # Session information
# Session information
utils::sessionInfo()## R version 4.5.1 (2025-06-13)
## Platform: x86_64-pc-linux-gnu
## Running under: Ubuntu 24.04.3 LTS
## 
## Matrix products: default
## BLAS:   /usr/lib/x86_64-linux-gnu/openblas-pthread/libblas.so.3 
## LAPACK: /usr/lib/x86_64-linux-gnu/openblas-pthread/libopenblasp-r0.3.26.so;  LAPACK version 3.12.0
## 
## locale:
##  [1] LC_CTYPE=C.UTF-8       LC_NUMERIC=C           LC_TIME=C.UTF-8       
##  [4] LC_COLLATE=C.UTF-8     LC_MONETARY=C.UTF-8    LC_MESSAGES=C.UTF-8   
##  [7] LC_PAPER=C.UTF-8       LC_NAME=C              LC_ADDRESS=C          
## [10] LC_TELEPHONE=C         LC_MEASUREMENT=C.UTF-8 LC_IDENTIFICATION=C   
## 
## time zone: UTC
## tzcode source: system (glibc)
## 
## attached base packages:
## [1] stats4    stats     graphics  grDevices utils     datasets  methods  
## [8] base     
## 
## other attached packages:
##  [1] escheR_1.8.0                ggplot2_4.0.0              
##  [3] STexampleData_1.16.0        SpatialExperiment_1.18.1   
##  [5] SingleCellExperiment_1.30.1 SummarizedExperiment_1.38.1
##  [7] Biobase_2.68.0              GenomicRanges_1.60.0       
##  [9] GenomeInfoDb_1.44.3         IRanges_2.42.0             
## [11] S4Vectors_0.46.0            MatrixGenerics_1.20.0      
## [13] matrixStats_1.5.0           ExperimentHub_2.16.1       
## [15] AnnotationHub_3.16.1        BiocFileCache_2.16.2       
## [17] dbplyr_2.5.1                BiocGenerics_0.54.1        
## [19] generics_0.1.4              SpotSweeper_1.5.0          
## [21] BiocStyle_2.36.0           
## 
## loaded via a namespace (and not attached):
##  [1] DBI_1.2.3               rlang_1.1.6             magrittr_2.0.4         
##  [4] spatialEco_2.0-3        compiler_4.5.1          RSQLite_2.4.3          
##  [7] png_0.1-8               systemfonts_1.3.1       vctrs_0.6.5            
## [10] pkgconfig_2.0.3         crayon_1.5.3            fastmap_1.2.0          
## [13] magick_2.9.0            XVector_0.48.0          labeling_0.4.3         
## [16] scuttle_1.18.0          rmarkdown_2.30          UCSC.utils_1.4.0       
## [19] ragg_1.5.0              purrr_1.1.0             bit_4.6.0              
## [22] xfun_0.53               cachem_1.1.0            beachmat_2.24.0        
## [25] jsonlite_2.0.0          blob_1.2.4              DelayedArray_0.34.1    
## [28] BiocParallel_1.42.2     terra_1.8-70            parallel_4.5.1         
## [31] R6_2.6.1                bslib_0.9.0             RColorBrewer_1.1-3     
## [34] jquerylib_0.1.4         Rcpp_1.1.0              bookdown_0.45          
## [37] knitr_1.50              Matrix_1.7-3            tidyselect_1.2.1       
## [40] abind_1.4-8             yaml_2.3.10             codetools_0.2-20       
## [43] curl_7.0.0              lattice_0.22-7          tibble_3.3.0           
## [46] withr_3.0.2             KEGGREST_1.48.1         S7_0.2.0               
## [49] evaluate_1.0.5          desc_1.4.3              Biostrings_2.76.0      
## [52] pillar_1.11.1           BiocManager_1.30.26     filelock_1.0.3         
## [55] BiocVersion_3.21.1      scales_1.4.0            glue_1.8.0             
## [58] tools_4.5.1             BiocNeighbors_2.2.0     fs_1.6.6               
## [61] grid_4.5.1              AnnotationDbi_1.70.0    GenomeInfoDbData_1.2.14
## [64] cli_3.6.5               rappdirs_0.3.3          textshaping_1.0.4      
## [67] S4Arrays_1.8.1          viridisLite_0.4.2       dplyr_1.1.4            
## [70] gtable_0.3.6            sass_0.4.10             digest_0.6.37          
## [73] SparseArray_1.8.1       rjson_0.2.23            htmlwidgets_1.6.4      
## [76] farver_2.1.2            memoise_2.0.1           htmltools_0.5.8.1      
## [79] pkgdown_2.1.3           lifecycle_1.0.4         httr_1.4.7             
## [82] mime_0.13               bit64_4.6.0-1           MASS_7.3-65