Contents

0.1 Introduction

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.

0.2 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")

0.3 Spot-level local outlier detection

0.3.1 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]

0.3.2 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"

0.3.3 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)

0.3.4 Visualizing local outliers

We can visualize the local outliers using the plotQC 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.

library(escheR)
## Loading required package: ggplot2
library(ggpubr)

# library size
p1 <- plotQC(spe,
    metric = "sum_log",
    outliers = "sum_outliers", point_size = 1.1
) +
    ggtitle("Library Size")

# unique genes
p2 <- plotQC(spe,
    metric = "detected_log",
    outliers = "detected_outliers", point_size = 1.1
) +
    ggtitle("Unique Genes")

# mitochondrial percent
p3 <- plotQC(spe,
    metric = "subsets_Mito_percent",
    outliers = "subsets_Mito_percent_outliers", point_size = 1.1
) +
    ggtitle("Mitochondrial Percent")

# all local outliers
p4 <- plotQC(spe,
    metric = "sum_log",
    outliers = "local_outliers", point_size = 1.1, stroke = 0.75
) +
    ggtitle("All Local Outliers")

# plot
plot_list <- list(p1, p2, p3, p4)
ggarrange(
    plotlist = plot_list,
    ncol = 2, nrow = 2,
    common.legend = FALSE
)

0.4 Removing technical artifacts using SpotSweeper

0.4.1 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"

0.4.2 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 plotQC function. This function plots the Visium spots with the specified QC metric.We’ll then again arrange these plots using ggpubr::arrange.

# library size
p1 <- plotQC(spe,
    metric = "sum_umi",
    outliers = NULL, point_size = 1.1
) +
    ggtitle("Library Size")

# unique genes
p2 <- plotQC(spe,
    metric = "sum_gene",
    outliers = NULL, point_size = 1.1
) +
    ggtitle("Unique Genes")

# mitochondrial percent
p3 <- plotQC(spe,
    metric = "expr_chrM_ratio",
    outliers = NULL, point_size = 1.1
) +
    ggtitle("Mitochondrial Percent")

# plot
plot_list <- list(p1, p2, p3)
ggarrange(
    plotlist = plot_list,
    ncol = 3, nrow = 1,
    common.legend = FALSE
)

0.4.3 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_rings=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_rings = 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"

0.4.4 Visualizing artifacts

We can visualize the artifacts using the escheR package. Here, we’ll visualize the artifacts using the plotQC function and arrange these plots using ggpubr::arrange.

plotQC(spe,
    metric = "expr_chrM_ratio",
    outliers = "artifact", point_size = 1.1
) +
    ggtitle("Hangnail artifact")

# Session information

utils::sessionInfo()
## R version 4.4.0 (2024-04-24)
## Platform: x86_64-pc-linux-gnu
## Running under: Ubuntu 22.04.4 LTS
## 
## Matrix products: default
## BLAS:   /home/biocbuild/bbs-3.19-bioc/R/lib/libRblas.so 
## LAPACK: /usr/lib/x86_64-linux-gnu/lapack/liblapack.so.3.10.0
## 
## locale:
##  [1] LC_CTYPE=en_US.UTF-8       LC_NUMERIC=C              
##  [3] LC_TIME=en_GB              LC_COLLATE=C              
##  [5] LC_MONETARY=en_US.UTF-8    LC_MESSAGES=en_US.UTF-8   
##  [7] LC_PAPER=en_US.UTF-8       LC_NAME=C                 
##  [9] LC_ADDRESS=C               LC_TELEPHONE=C            
## [11] LC_MEASUREMENT=en_US.UTF-8 LC_IDENTIFICATION=C       
## 
## time zone: America/New_York
## tzcode source: system (glibc)
## 
## attached base packages:
## [1] stats4    stats     graphics  grDevices utils     datasets  methods  
## [8] base     
## 
## other attached packages:
##  [1] ggpubr_0.6.0                escheR_1.4.0               
##  [3] ggplot2_3.5.1               STexampleData_1.12.0       
##  [5] SpatialExperiment_1.14.0    SingleCellExperiment_1.26.0
##  [7] SummarizedExperiment_1.34.0 Biobase_2.64.0             
##  [9] GenomicRanges_1.56.0        GenomeInfoDb_1.40.0        
## [11] IRanges_2.38.0              S4Vectors_0.42.0           
## [13] MatrixGenerics_1.16.0       matrixStats_1.3.0          
## [15] ExperimentHub_2.12.0        AnnotationHub_3.12.0       
## [17] BiocFileCache_2.12.0        dbplyr_2.5.0               
## [19] BiocGenerics_0.50.0         SpotSweeper_1.0.0          
## [21] BiocStyle_2.32.0           
## 
## loaded via a namespace (and not attached):
##  [1] DBI_1.2.2                 rlang_1.1.3              
##  [3] magrittr_2.0.3            spatialEco_2.0-2         
##  [5] compiler_4.4.0            RSQLite_2.3.6            
##  [7] DelayedMatrixStats_1.26.0 png_0.1-8                
##  [9] vctrs_0.6.5               pkgconfig_2.0.3          
## [11] crayon_1.5.2              fastmap_1.2.0            
## [13] backports_1.4.1           magick_2.8.3             
## [15] XVector_0.44.0            labeling_0.4.3           
## [17] scuttle_1.14.0            utf8_1.2.4               
## [19] rmarkdown_2.26            UCSC.utils_1.0.0         
## [21] tinytex_0.51              purrr_1.0.2              
## [23] bit_4.0.5                 xfun_0.44                
## [25] zlibbioc_1.50.0           cachem_1.1.0             
## [27] beachmat_2.20.0           jsonlite_1.8.8           
## [29] blob_1.2.4                highr_0.10               
## [31] DelayedArray_0.30.1       BiocParallel_1.38.0      
## [33] terra_1.7-71              broom_1.0.5              
## [35] parallel_4.4.0            R6_2.5.1                 
## [37] bslib_0.7.0               car_3.1-2                
## [39] jquerylib_0.1.4           Rcpp_1.0.12              
## [41] bookdown_0.39             knitr_1.46               
## [43] Matrix_1.7-0              tidyselect_1.2.1         
## [45] abind_1.4-5               yaml_2.3.8               
## [47] codetools_0.2-20          curl_5.2.1               
## [49] lattice_0.22-6            tibble_3.2.1             
## [51] withr_3.0.0               KEGGREST_1.44.0          
## [53] evaluate_0.23             Biostrings_2.72.0        
## [55] pillar_1.9.0              BiocManager_1.30.23      
## [57] filelock_1.0.3            carData_3.0-5            
## [59] generics_0.1.3            BiocVersion_3.19.1       
## [61] sparseMatrixStats_1.16.0  munsell_0.5.1            
## [63] scales_1.3.0              glue_1.7.0               
## [65] tools_4.4.0               BiocNeighbors_1.22.0     
## [67] ggsignif_0.6.4            cowplot_1.1.3            
## [69] grid_4.4.0                tidyr_1.3.1              
## [71] AnnotationDbi_1.66.0      colorspace_2.1-0         
## [73] GenomeInfoDbData_1.2.12   cli_3.6.2                
## [75] rappdirs_0.3.3            fansi_1.0.6              
## [77] viridisLite_0.4.2         S4Arrays_1.4.0           
## [79] dplyr_1.1.4               gtable_0.3.5             
## [81] rstatix_0.7.2             sass_0.4.9               
## [83] digest_0.6.35             SparseArray_1.4.4        
## [85] farver_2.1.2              rjson_0.2.21             
## [87] memoise_2.0.1             htmltools_0.5.8.1        
## [89] lifecycle_1.0.4           httr_1.4.7               
## [91] mime_0.12                 bit64_4.0.5              
## [93] MASS_7.3-60.2