1 Installation

The package can be installed using bioconductor install manager:

if (!require("BiocManager", quietly = TRUE))
    install.packages("BiocManager")

BiocManager::install("ClusterFoldSimilarity")
library(ClusterFoldSimilarity)

2 Introduction

Comparing single-cell data across different datasets, samples and batches has demonstrated to be challenging. ClusterFoldSimilarity aims to solve the complexity of comparing different single-cell datasets by computing similarity scores between clusters (or user-defined groups) from any number of independent single-cell experiments, including different species and sequencing technologies. It accomplishes this by identifying analogous fold-change patterns across cell groups that share a common set of features (such as genes). Additionally, it selects and reports the top important features that have contributed to the observed similarity, serving as a tool for feature selection.

The output is a table that contains the similarity values for all the combinations of cluster-pairs from the independent datasets. ClusterFoldSimilarity also includes various plotting utilities to enhance the interpretability of the similarity scores.

2.0.1 Cross-species analysis and sequencing technologies (e.g.: Human vs Mouse, ATAC-Seq vs RNA-Seq)

ClusterFoldSimilarity is able to compare any number of independent experiments, including different organisms, making it useful for matching cell populations across different organisms, and thus, useful for inter-species analysis. Additionally, it can be used with single-cell RNA-Seq data, single-cell ATAC-Seq data, or more broadly, with continuous numerical data that shows changes in feature abundance across a set of common features between different groups.

2.0.2 Compatibility

It can be easily integrated on any existing single-cell analysis pipeline, and it is compatible with the most used single-cell objects: Seurat and SingleCellExperiment.

Parallel computing is available through the option parallel=TRUE which make use of BiocParallel.

3 Using ClusterFoldSimilarity to find similar clusters/cell-groups across datasets

Typically, ClusterFoldSimilarity will receive as input either a list of two or more Seurat or SingleCellExperiment objects.

ClusterFoldSimilarity will obtain the raw count data from these objects ( GetAssayData(assay, slot = "counts") in the case of Seurat, or counts() for SingleCellExperiment object), and group or cluster label information (using Idents() function from Seurat, or colLabels() for SingleCellExperiment ).

For the sake of illustration, we will employ the scRNAseq package, which contains numerous individual-cell datasets ready for download and encompassing samples from both human and mouse origins. In this example, we specifically utilize 2 human single-cell datasets obtained from the pancreas.

library(Seurat)
library(scRNAseq)
library(dplyr)
# Human pancreatic single cell data 1
pancreasMuraro <- scRNAseq::MuraroPancreasData(ensembl=FALSE)
pancreasMuraro <- pancreasMuraro[,rownames(colData(pancreasMuraro)[!is.na(colData(pancreasMuraro)$label),])]
colData(pancreasMuraro)$cell.type <- colData(pancreasMuraro)$label
rownames(pancreasMuraro) <- make.names(unlist(lapply(strsplit(rownames(pancreasMuraro), split="__"), function(x)x[[1]])), unique = TRUE)
singlecell1Seurat <- CreateSeuratObject(counts=counts(pancreasMuraro), meta.data=as.data.frame(colData(pancreasMuraro)))

Table 1: Cell-types on pancreas dataset from Muraro et al.
Var1 Freq
acinar 219
alpha 812
beta 448
delta 193
duct 245
endothelial 21
epsilon 3
mesenchymal 80
pp 101
unclear 4
# Human pancreatic single cell data 2
pancreasBaron <- scRNAseq::BaronPancreasData(which="human", ensembl=FALSE)
colData(pancreasBaron)$cell.type <- colData(pancreasBaron)$label
rownames(pancreasBaron) <- make.names(rownames(pancreasBaron), unique = TRUE)

singlecell2Seurat <- CreateSeuratObject(counts=counts(pancreasBaron), meta.data=as.data.frame(colData(pancreasBaron)))

Table 2: Cell-types on pancreas dataset from Baron et al.
Var1 Freq
acinar 958
activated_stellate 284
alpha 2326
beta 2525
delta 601
ductal 1077
endothelial 252
epsilon 18
gamma 255
macrophage 55
mast 25
quiescent_stellate 173
schwann 13
t_cell 7

As we want to perform clustering analysis for later comparison of these cluster groups using ClusterFoldSimilarity, we first need to normalize and identify variable features for each dataset independently.

Note: these steps should be done tailored to each independent dataset, here we apply the same parameters for the sake of simplicity:

# Create a list with the unprocessed single-cell datasets
singlecellObjectList <- list(singlecell1Seurat, singlecell2Seurat)
# Apply the same processing to each dataset and return a list of single-cell analysis
singlecellObjectList <- lapply(X=singlecellObjectList, FUN=function(scObject){
scObject <- NormalizeData(scObject)
scObject <- FindVariableFeatures(scObject, selection.method="vst", nfeatures=2000)
scObject <- ScaleData(scObject, features=VariableFeatures(scObject))
scObject <- RunPCA(scObject, features=VariableFeatures(object=scObject))
scObject <- FindNeighbors(scObject, dims=seq(16))
scObject <- FindClusters(scObject, resolution=0.4)
})

Once we have all of our single-cell datasets analyzed independently, we can compute the similarity values. clusterFoldSimilarity() takes as arguments:

  • scList: a list of single-cell objects (mandatory) either of class Seurat or of class SingleCellExperiment.
  • sampleNames: vector with names for each of the datasets. If not set the datasets will be named in the given order as: 1, 2, …, N.
  • topN: the top n most similar clusters/groups to report for each cluster/group (default: 1, the top most similar cluster). If set to Inf it will return the values from all the possible cluster-pairs.
  • topNFeatures: the top n features (e.g.: genes) that contribute to the observed similarity between the pair of clusters (default: 1, the top contributing gene). If a negative number, the tool will report the n most dissimilar features.
  • nSubsampling: number of subsamplings (1/3 of cells on each iteration) at group level for calculating the fold-changes (default: 15). At start, the tool will report a message with the recommended number of subsamplings for the given data (average n of subsamplings needed to observe all cells).
  • parallel: whether to use parallel computing with multiple threads or not (default: FALSE). If we want to use a specific single-cell experiment for annotation (from which we know a ground-truth label, e.g. cell type, cell cycle, treatment… etc.), we can use that label to directly compare the single-cell datasets.

Here we will use the annotated pancreas cell-type labels from the dataset 1 to illustrate how to match clusters to cell-types using a reference dataset:

# Assign cell-type annotated from the original study to the cell labels:
Idents(singlecellObjectList[[1]]) <- factor(singlecellObjectList[[1]][[]][, "cell.type"])

library(ClusterFoldSimilarity)
similarityTable <- clusterFoldSimilarity(scList=singlecellObjectList, 
                                         sampleNames=c("human", "humanNA"),
                                         topN=1, 
                                         nSubsampling=24)

Table 3: A table of the first 10 rows of the similarity results.
similarityValue w datasetL clusterL datasetR clusterR topFeatureConserved featureScore
acinar 50.38169 6.82 human acinar humanNA 2 REG1A 50.427213
alpha 34.65354 6.51 human alpha humanNA 3 GCG 34.122804
beta 33.92137 6.59 human beta humanNA 4 INS 30.337947
delta 29.38112 6.29 human delta humanNA 4 PCSK1 7.573988
duct 33.42910 6.42 human duct humanNA 5 SPP1 13.135512
endothelial 45.73913 7.01 human endothelial humanNA 9 PLVAP 18.286168

A data.frame with the results is returned containing:

  • similarityValue: The top similarity value calculated between datasetL:clusterL and datasetR.
  • w: Weight associated with the similarity score value.
  • datasetL: Dataset left, the dataset/sample which has been used to be compared.
  • clusterL: Cluster left, the cluster source from datasetL which has been compared.
  • datasetR: Dataset right, the dataset/sample used for comparison against datasetL.
  • clusterR: Cluster right, the cluster target from datasetR which is being compared with the clusterL from datasetL.
  • topFeatureConserved: The features (e.g.: genes, peaks…) that most contributed to the similarity between clusterL & clusterR.
  • featureScore: The similarity score contribution for the specific topFeatureConserved (e.g.: genes, peaks…).

By default, clusterFoldSimilarity() will plot a graph network that visualizes the connections between the clusters from the different datasets using the similarity table that has been obtained. The arrows point in the direction of the similarity (datasetL:clusterL -> datasetR:clusterR); it can be useful for identifying relationships between groups of clusters and cell-populations that tend to be more similar. The graph plot can also be obtained by using the function plotClustersGraph() from this package, using as input the similarity table.

In this example, as we have information regarding cell-type labels, we can check how the cell types match by calculating the most abundant cell type on each of the similar clusters:

typeCount <- singlecellObjectList[[2]][[]] %>% 
  group_by(seurat_clusters) %>% 
  count(cell.type) %>%
  arrange(desc(n), .by_group = TRUE) %>% 
  filter(row_number()==1)

Table 4: Cell-type label matching on pancreas single-cell data.
seurat_clusters cell.type n matched.type
0 alpha 1308 alpha
1 beta 1022 beta
2 acinar 927 acinar
3 alpha 967 alpha
4 beta 958 beta
5 ductal 817 duct
6 delta 586 delta
7 beta 494 beta
8 activated_stellate 283 mesenchymal
9 endothelial 250 endothelial
10 gamma 211 epsilon
11 ductal 174 acinar
12 macrophage 55 unclear
13 ductal 77 unclear

3.1 Retrieving the top-n similarities

If we suspect that clusters could be related with more than one cluster of other datasets, we can retrieve the top n similarities for each cluster:

# Retrieve the top 3 similar cluster for each of the clusters:
similarityTable3Top <- clusterFoldSimilarity(scList=singlecellObjectList, 
                                             topN=3,
                                             sampleNames=c("human", "humanNA"), 
                                             nSubsampling=24)

Table 5: Similarity results showing the top 3 most similar clusters.
similarityValue w datasetL clusterL datasetR clusterR topFeatureConserved featureScore
acinar.3 50.26528 6.84 human acinar humanNA 2 REG1A 49.727196
acinar.12 42.79969 6.65 human acinar humanNA 11 REG3A 40.431410
acinar.6 26.84455 5.95 human acinar humanNA 5 SERPINA3 14.144226
alpha.4 34.54703 6.51 human alpha humanNA 3 GCG 33.834157
alpha.1 31.30748 6.34 human alpha humanNA 0 GCG 29.992194
alpha.5 29.91542 6.37 human alpha humanNA 4 INS 8.058919

3.2 Obtaining the top-n feature markers

If we are interested on the features that contribute the most to the similarity, we can retrieve the top n features:

# Retrieve the top 5 features that contribute the most to the similarity between each pair of clusters:
similarityTable5TopFeatures <- clusterFoldSimilarity(scList=singlecellObjectList, 
                                                     topNFeatures=5, 
                                                     nSubsampling=24)

Table 6: Similarity results showing the top 5 features that most contributed to the similarity.
similarityValue w datasetL clusterL datasetR clusterR topFeatureConserved featureScore
acinar.11 50.17496 6.82 1 acinar 2 2 REG1A 50.56400
acinar.12 50.17496 6.82 1 acinar 2 2 CTRB2 46.65210
acinar.13 50.17496 6.82 1 acinar 2 2 REG1B 46.30476
acinar.14 50.17496 6.82 1 acinar 2 2 PRSS1 45.14041
acinar.15 50.17496 6.82 1 acinar 2 2 CELA3A 40.36778
alpha.16 34.56059 6.52 1 alpha 2 3 GCG 33.88686
alpha.17 34.56059 6.52 1 alpha 2 3 TTR 22.29361
alpha.18 34.56059 6.52 1 alpha 2 3 CHGB 13.28128
alpha.19 34.56059 6.52 1 alpha 2 3 PCSK2 11.74380
alpha.20 34.56059 6.52 1 alpha 2 3 TM4SF4 11.52553

3.3 Retrieving all the similarity values and plotting a similarity heatmap

Sometimes it is useful to retrieve all the similarity values for downstream analysis (e.g. identify more than one cluster that is similar to a cluster of interest, finding the most dissimilar clusters, etc). To obtain all the values, we need to specify topN=Inf.

By default, clusterFoldSimilarity creates a heatmap plot with the computed similarity values (from the perspective of the first dataset found on scList; to modify this plot see the following section). The top 2 similarities for each group within dataset 1 (heatmap row-wise) are highlighted with colored borders.

similarityTableAllValues <- clusterFoldSimilarity(scList=singlecellObjectList, 
                                                  sampleNames=c("human", "humanNA"), 
                                                  topN=Inf)

dim(similarityTableAllValues)
## [1] 280   8

For downstream analysis of the similarities, it can be convenient to create a matrix with all the scores from the comparison of two datasets:

library(dplyr)
dataset1 <- "human"
dataset2 <- "humanNA"
similarityTable2 <- similarityTableAllValues %>% 
  filter(datasetL == dataset1 & datasetR == dataset2) %>% 
  arrange(desc(as.numeric(clusterL)), as.numeric(clusterR))
cls <- unique(similarityTable2$clusterL)
cls2 <- unique(similarityTable2$clusterR)
similarityMatrixAll <- t(matrix(similarityTable2$similarityValue, ncol=length(unique(similarityTable2$clusterL))))
rownames(similarityMatrixAll) <- cls
colnames(similarityMatrixAll) <- cls2
Table 7: A 2 datasets Similarity matrix.
0 1 2 3 4 5 6 7 8 9 10 11 12 13
acinar -10.80 31.40 23.84 24.02 -9.33 -12.06 -1.36 -11.17 23.70 -15.89 -11.12 27.56 31.66 27.28
alpha -10.71 -11.53 -8.71 -9.28 22.99 -17.53 50.36 -10.54 -11.14 -12.06 28.71 -11.45 -26.14 17.14
beta 6.59 -19.93 -16.49 34.59 25.98 26.02 -13.56 -15.64 22.63 -16.84 26.66 -14.32 -17.47 29.83
delta 33.98 29.61 -15.33 -14.77 15.82 -14.87 25.13 -14.71 26.90 -14.98 -14.93 -13.58 33.63 16.49
duct -8.69 16.00 -13.98 26.65 -15.48 25.88 25.63 29.02 -12.34 -13.29 19.79 -13.85 24.40 -11.77
endothelial -19.19 27.78 33.01 28.57 -15.60 -14.60 21.96 -15.94 23.96 -9.38 -6.09 -17.48 -15.73 -14.73
epsilon 12.20 32.73 -13.66 48.53 -15.65 -9.01 -9.08 -17.79 -15.95 -15.11 8.82 45.50 -9.37 32.36
mesenchymal -15.53 20.84 -18.31 25.19 21.79 23.88 -14.26 -13.68 29.69 -16.25 27.21 -10.93 42.79 -18.99
pp -19.21 -17.72 29.76 -1.80 -15.37 -5.46 -12.68 30.29 -11.34 -14.40 -13.70 -12.21 11.23 24.74
unclear 26.92 20.35 -13.82 34.81 22.24 -20.49 -20.01 -17.86 30.24 17.47 25.49 2.81 -17.92 41.86

4 Using ClusterFoldSimilarity across species and numerous datasets:

ClusterFoldSimilarity can compare any number of independent studies, including different organisms, making it useful for inter-species analysis. Also, it can be used on different sequencing data technologies: e.g.: compare single-cell ATAC-Seq VS RNA-seq.

In this example, we are going to add a pancreas single-cell dataset from Mouse to the 2 existing ones from Human that we have processed in the previous steps.

# Mouse pancreatic single cell data
pancreasBaronMM <- scRNAseq::BaronPancreasData(which="mouse", ensembl=FALSE)

Table 8: Cell-types on pancreas dataset from Baron et al.
Var1 Freq
B_cell 10
T_cell 7
activated_stellate 14
alpha 191
beta 894
delta 218
ductal 275
endothelial 139
gamma 41
immune_other 8
macrophage 36
quiescent_stellate 47
schwann 6
colData(pancreasBaronMM)$cell.type <- colData(pancreasBaronMM)$label
# Translate mouse gene ids to human ids
# *for the sake of simplicity we are going to transform to uppercase all mouse gene names
rownames(pancreasBaronMM) <- make.names(toupper(rownames(pancreasBaronMM)), unique=TRUE)
# Create seurat object
singlecell3Seurat <- CreateSeuratObject(counts=counts(pancreasBaronMM), meta.data=as.data.frame(colData(pancreasBaronMM)))

# We append the single-cell object to our list
singlecellObjectList[[3]] <- singlecell3Seurat

Now, we process the new single-cell dataset from mouse, and we calculate the similarity scores between the 3 independent datasets.

scObject <- singlecellObjectList[[3]]
scObject <- NormalizeData(scObject)
scObject <- FindVariableFeatures(scObject, selection.method="vst", nfeatures=2000)
scObject <- ScaleData(scObject, features=VariableFeatures(scObject))
scObject <- RunPCA(scObject, features=VariableFeatures(object=scObject))
scObject <- FindNeighbors(scObject, dims=seq(16))
scObject <- FindClusters(scObject, resolution=0.4)
singlecellObjectList[[3]] <- scObject

This time we will make use of the option parallel=TRUE. We can set the specific number of CPUs to use using BiocParallel::register()

# We use the cell labels as a second reference, but we can also use the cluster labels if our interest is to match clusters
Idents(singlecellObjectList[[3]]) <- factor(singlecellObjectList[[3]][[]][,"cell.type"])

# We subset the most variable genes in each experiment
singlecellObjectListVariable <- lapply(singlecellObjectList, function(x){x[VariableFeatures(x),]})

# Setting the number of CPUs with BiocParallel:
BiocParallel::register(BPPARAM =  BiocParallel::MulticoreParam(workers = 6))

similarityTableHumanMouse <- clusterFoldSimilarity(scList=singlecellObjectListVariable,
                                                        sampleNames=c("human", "humanNA", "mouse"),
                                                        topN=1, 
                                                        nSubsampling=24,
                                                        parallel=TRUE)

We can compute and visualize with a heatmap all the similarities for each cluster/group of cells from the 3 datasets using topN=Inf. Additionally, we can use the function similarityHeatmap() from this package to plot the heatmap with the datasets in a different order, or just plot the 2 datasets we are interested in. The top 2 similarities are highlighted to help visualizing the best matching groups.

similarityTableHumanMouseAll <- clusterFoldSimilarity(scList=singlecellObjectListVariable,
                                                          sampleNames=c("human", "humanNA", "mouse"),
                                                          topN=Inf,
                                                          nSubsampling=24,
                                                          parallel=TRUE)

As the similarity values might not be symmetric (e.g. a cluster A from D1 showing the top similarity to B from D2, might not be the top similar cluster to B from D2), we can select which dataset to plot in the Y-axis:

ClusterFoldSimilarity::similarityHeatmap(similarityTable=similarityTableHumanMouseAll, mainDataset="humanNA")

Additionally, we can turn-off the highlight using highlightTop=FALSE

# Turn-off the highlighting:
ClusterFoldSimilarity::similarityHeatmap(similarityTable=similarityTableHumanMouseAll, mainDataset="humanNA", highlightTop=FALSE)

5 Similarity score calculation

ClusterFoldSimilarity does not need to integrate the data, or apply any batch correction techniques across the datasets that we aim to analyze, which makes it less prone to data-loss or noise. The similarity value is based on the fold-changes between clusters/groups of cells defined by the user. These fold-changes from different independent datasets are first computed using a Bayesian approach, we calculate this fold-change distribution using a permutation analysis that shrink the fold-changes with no biological meaning. These differences in abundance are then combined using a pairwise dot product approach, after adding these feature contributions and applying a fold-change concordance weight, a similarity value is obtained for each of the clusters of each of the datasets present.

Session information

sessionInfo()
## R version 4.4.0 RC (2024-04-16 r86468)
## Platform: x86_64-pc-linux-gnu
## Running under: Ubuntu 22.04.4 LTS
## 
## Matrix products: default
## BLAS:   /home/biocbuild/bbs-3.20-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] dplyr_1.1.4                 scRNAseq_2.17.10           
##  [3] SingleCellExperiment_1.27.0 SummarizedExperiment_1.35.0
##  [5] Biobase_2.65.0              GenomicRanges_1.57.0       
##  [7] GenomeInfoDb_1.41.0         IRanges_2.39.0             
##  [9] S4Vectors_0.43.0            BiocGenerics_0.51.0        
## [11] MatrixGenerics_1.17.0       matrixStats_1.3.0          
## [13] Seurat_5.0.3                SeuratObject_5.0.1         
## [15] sp_2.1-4                    ClusterFoldSimilarity_1.1.0
## [17] BiocStyle_2.33.0           
## 
## loaded via a namespace (and not attached):
##   [1] ProtGenerics_1.37.0      spatstat.sparse_3.0-3    bitops_1.0-7            
##   [4] httr_1.4.7               RColorBrewer_1.1-3       tools_4.4.0             
##   [7] sctransform_0.4.1        alabaster.base_1.5.0     utf8_1.2.4              
##  [10] R6_2.5.1                 HDF5Array_1.33.0         lazyeval_0.2.2          
##  [13] uwot_0.2.2               rhdf5filters_1.17.0      withr_3.0.0             
##  [16] gridExtra_2.3            progressr_0.14.0         cli_3.6.2               
##  [19] spatstat.explore_3.2-7   fastDummies_1.7.3        labeling_0.4.3          
##  [22] alabaster.se_1.5.0       sass_0.4.9               spatstat.data_3.0-4     
##  [25] ggridges_0.5.6           pbapply_1.7-2            Rsamtools_2.21.0        
##  [28] systemfonts_1.0.6        svglite_2.1.3            parallelly_1.37.1       
##  [31] rstudioapi_0.16.0        RSQLite_2.3.6            generics_0.1.3          
##  [34] BiocIO_1.15.0            ica_1.0-3                spatstat.random_3.2-3   
##  [37] Matrix_1.7-0             fansi_1.0.6              abind_1.4-5             
##  [40] lifecycle_1.0.4          yaml_2.3.8               rhdf5_2.49.0            
##  [43] SparseArray_1.5.0        BiocFileCache_2.13.0     Rtsne_0.17              
##  [46] grid_4.4.0               blob_1.2.4               promises_1.3.0          
##  [49] ExperimentHub_2.13.0     crayon_1.5.2             miniUI_0.1.1.1          
##  [52] lattice_0.22-6           cowplot_1.1.3            GenomicFeatures_1.57.0  
##  [55] KEGGREST_1.45.0          magick_2.8.3             pillar_1.9.0            
##  [58] knitr_1.46               rjson_0.2.21             future.apply_1.11.2     
##  [61] codetools_0.2-20         leiden_0.4.3.1           glue_1.7.0              
##  [64] data.table_1.15.4        vctrs_0.6.5              png_0.1-8               
##  [67] gypsum_1.1.0             spam_2.10-0              gtable_0.3.5            
##  [70] cachem_1.0.8             xfun_0.43                S4Arrays_1.5.0          
##  [73] mime_0.12                survival_3.6-4           tinytex_0.50            
##  [76] fitdistrplus_1.1-11      ROCR_1.0-11              nlme_3.1-164            
##  [79] bit64_4.0.5              alabaster.ranges_1.5.0   filelock_1.0.3          
##  [82] RcppAnnoy_0.0.22         bslib_0.7.0              irlba_2.3.5.1           
##  [85] KernSmooth_2.23-22       colorspace_2.1-0         DBI_1.2.2               
##  [88] tidyselect_1.2.1         bit_4.0.5                compiler_4.4.0          
##  [91] curl_5.2.1               httr2_1.0.1              xml2_1.3.6              
##  [94] ggdendro_0.2.0           DelayedArray_0.31.0      plotly_4.10.4           
##  [97] bookdown_0.39            rtracklayer_1.65.0       scales_1.3.0            
## [100] lmtest_0.9-40            rappdirs_0.3.3           stringr_1.5.1           
## [103] digest_0.6.35            goftest_1.2-3            spatstat.utils_3.0-4    
## [106] alabaster.matrix_1.5.0   rmarkdown_2.26           XVector_0.45.0          
## [109] htmltools_0.5.8.1        pkgconfig_2.0.3          highr_0.10              
## [112] dbplyr_2.5.0             fastmap_1.1.1            ensembldb_2.29.0        
## [115] rlang_1.1.3              htmlwidgets_1.6.4        UCSC.utils_1.1.0        
## [118] shiny_1.8.1.1            farver_2.1.1             jquerylib_0.1.4         
## [121] zoo_1.8-12               jsonlite_1.8.8           BiocParallel_1.39.0     
## [124] RCurl_1.98-1.14          magrittr_2.0.3           kableExtra_1.4.0        
## [127] GenomeInfoDbData_1.2.12  dotCall64_1.1-1          patchwork_1.2.0         
## [130] Rhdf5lib_1.27.0          munsell_0.5.1            Rcpp_1.0.12             
## [133] reticulate_1.36.1        stringi_1.8.3            alabaster.schemas_1.5.0 
## [136] zlibbioc_1.51.0          MASS_7.3-60.2            AnnotationHub_3.13.0    
## [139] plyr_1.8.9               parallel_4.4.0           listenv_0.9.1           
## [142] ggrepel_0.9.5            deldir_2.0-4             Biostrings_2.73.0       
## [145] splines_4.4.0            tensor_1.5               igraph_2.0.3            
## [148] spatstat.geom_3.2-9      RcppHNSW_0.6.0           paws.common_0.7.2       
## [151] reshape2_1.4.4           BiocVersion_3.20.0       XML_3.99-0.16.1         
## [154] evaluate_0.23            BiocManager_1.30.22      httpuv_1.6.15           
## [157] RANN_2.6.1               tidyr_1.3.1              purrr_1.0.2             
## [160] polyclip_1.10-6          future_1.33.2            scattermore_1.2         
## [163] alabaster.sce_1.5.0      ggplot2_3.5.1            paws.storage_0.5.0      
## [166] xtable_1.8-4             restfulr_0.0.15          AnnotationFilter_1.29.0 
## [169] RSpectra_0.16-1          later_1.3.2              viridisLite_0.4.2       
## [172] tibble_3.2.1             memoise_2.0.1            AnnotationDbi_1.67.0    
## [175] GenomicAlignments_1.41.0 cluster_2.1.6            globals_0.16.3