Over the last decade, we have seen exponential growth of the scale of scRNA-seq datasets to millions of cells sequenced in a single study. This has enabled researchers to characterize the gene expression profiles of various cell types across tissues. The rapid growth of scRNA-seq data has also created an unique set of challenges, for instance, there is a pressing need for scalable approaches for scRNA-seq data visualization.
This vignette introduces scBubbletree, a transparent workflow for quantitative exploration of single cell RNA-seq data.
In short, the algorithm of scBubbletree performs clustering to identify clusters (“bubbles”) of transcriptionally similar cells, and then visualizes these clusters as leafs in a hierarchical dendrogram (“bubbletree”) which describes their natural relationships. The workflow comprises four steps: 1. determining the clustering resolution, 2. clustering, 3. hierarchical cluster grouping and 4. visualization. We explain each step in the following using real scRNA-seq dataset of five cancer cell lines.
To run this vignette we need to load a few R-packages:
library(scBubbletree)
library(ggplot2)
library(ggtree)
library(patchwork)
Here we will analyze a scRNA-seq dataset containing a mixture of 3,918 cells from five human lung adenocarcinoma cell lines (HCC827, H1975, A549, H838 and H2228). The dataset is available here2 https://github.com/LuyiTian/sc_mixology/blob/master/data/ sincell_with_class_5cl.RData.
The library has been prepared with 10x Chromium platform and sequenced with Illumina NextSeq 500 platform. Raw data has been processed with Cellranger. The tool demuxlet has been used to predict the identity of each cell based on known genetic differences between the different cell lines.
Data processing was performed with R-package Seurat. Gene expressions were normalized with the function SCTransform using default parameters, and principal component analysis (PCA) was performed with function RunPCA based on the 5,000 most variable genes in the dataset identified with the function FindVariableFeatures.
In both datasets we saw that the first 15 principal components capture most of the variance in the data, and the proportion of variance explained by each subsequent principal component was negligible. Thus, we used the single cell projections (embeddings) in 15-dimensional feature space, \(A^{3,918\times 15}\).
# # This script can be used to generate data("d_ccl", package = "scBubbletree")
#
# # create directory
# dir.create(path = "case_study/")
#
# # download the data from:
# https://github.com/LuyiTian/sc_mixology/raw/master/data/
# sincell_with_class_5cl.RData
#
# # load the data
# load(file = "case_study/sincell_with_class_5cl.RData")
#
# # we are only interested in the 10x data object 'sce_sc_10x_5cl_qc'
# d <- sce_sc_10x_5cl_qc
#
# # remove the remaining objects (cleanup)
# rm(sc_Celseq2_5cl_p1, sc_Celseq2_5cl_p2, sc_Celseq2_5cl_p3, sce_sc_10x_5cl_qc)
#
# # get the meta data for each cell
# meta <- colData(d)[,c("cell_line_demuxlet","non_mt_percent","total_features")]
#
# # create Seurat object from the raw counts and append the meta data to it
# d <- Seurat::CreateSeuratObject(counts = d@assays$data$counts,
# project = '')
#
# # check if all cells are matched between d and meta
# # table(rownames(d@meta.data) == meta@rownames)
# d@meta.data <- cbind(d@meta.data, meta@listData)
#
# # cell type predictions are provided as part of the meta data
# table(d@meta.data$cell_line)
#
# # select 5,000 most variable genes
# d <- Seurat::FindVariableFeatures(object = d,
# selection.method = "vst",
# nfeatures = 5000)
#
# # Preprocessing with Seurat: SCT transformation + PCA
# d <- SCTransform(object = d,
# variable.features.n = 5000)
# d <- RunPCA(object = d,
# npcs = 50,
# features = VariableFeatures(object = d))
#
# # perform UMAP + t-SNE
# d <- RunUMAP(d, dims = 1:15)
# d <- RunTSNE(d, dims = 1:15)
#
# # save the preprocessed data
# save(d, file = "case_study/d.RData")
#
# # save the PCA matrix 'A', meta data 'm' and
# # marker genes matrix 'e'
# d <- get(load(file ="case_study/d.RData"))
# A <- d@reductions$pca@cell.embeddings[, 1:15]
# m <- d@meta.data
# e <- t(as.matrix(d@assays$SCT@data[
# rownames(d@assays$SCT@data) %in%
# c("ALDH1A1",
# "PIP4K2C",
# "SLPI",
# "CT45A2",
# "CD74"), ]))
#
# d_ccl <- list(A = A, m = m, e = e)
# save(d_ccl, file = "data/d_ccl.RData")
Load the processed PCA matrix and the meta data
data("d_ccl", package = "scBubbletree")
A <- d_ccl$A
m <- d_ccl$m
e <- d_ccl$e
We will analyze this data with scBubbletree.
As first input scBubbletree uses matrix \(A^{n\times f}\) which represents a low-dimensional projection of the original scRNA-seq data, with \(n\) rows as cells and \(f\) columns as low-dimension features.
We will use the PCA data generated by Seurat as \(A\). In particular, we will use the first 15 principal components (PCs) as every additional PC explains negligible amount of variance in the data.
Important remark about \(A\): the scBubbletree workflow works directly with the numeric matrix \(A^{n\times f}\) and is agnostic to the initial data processing protocol. This enables seamless integration of scBubbletree with computational pipelines using objects generated by the R-packages Seurat and SingleCellExperiment. The users simply have to extract \(A\) from the corresponding Seurat or SingleCellExperiment objects.
# A has n=cells as rows, f=features as columns (e.g. from PCA)
dim(A)
FALSE [1] 3918 15
The scBubbletree workflow performs the following steps:
If we use graph-based community detection (GCD, recommended for scRNA-seq) with e.g. the Louvain or Leiden method, then we need to find appropriate value for the resolution parameter \(r\). Otherwise, we can use the simpler k-means clustering algorithm in which case we need to find an appropriate value of the number of clusters \(k\). In the next we will try both, GCD and k-means, clustering.
How many clusters (cell types) are there are in the data? Can we guess a reasonable value of \(k\)?
To find a reasonable value of \(k\) we can study the literature or databases
such as the human protein atlas database (HPA). We can also use the function
get_k
for data-driven inference of \(k\) based on the Gap statistic and the
within-cluster sum of squares (WCSS).
As this is a toy dataset, we will skip the first approach and perform a
data-driven search for \(k\) using get_k
. As input we need to provide the
matrix \(A\) as input and a vector of \(k\)s. The output will be the Gap
statistic and WCSS estimates for each \(k\).
Lets run get_k
now:
b_k <- get_k(B_gap = 5,
ks = 1:10,
x = A,
n_start = 50,
iter_max = 200,
kmeans_algorithm = "MacQueen",
cores = 1)
The Gap statistic and WCSS curves have a noticeable knee (elbow) at \(k=5\).
Hence, \(k\)=5 appears to be reasonable first choice of \(k\). Means (points) and
95% confidence intervals are shown for the Gap statistic at each \(k\) computed
using B_gap
=5 MCMC simulations.
g0 <- ggplot(data = b_k$gap_stats_summary)+
geom_line(aes(x = k, y = gap_mean))+
geom_point(aes(x = k, y = gap_mean), size = 1)+
geom_errorbar(aes(x = k, y = gap_mean, ymin = L95, ymax = H95), width = 0.1)+
ylab(label = "Gap")|
ggplot(data = b_k$wcss_stats_summary)+
geom_line(aes(x = k, y = wcss_mean))+
geom_point(aes(x = k, y = wcss_mean), size = 1)+
ylab(label = "WCSS")+
scale_y_log10()+
annotation_logticks(base = 10, sides = "l")
g0
For Louvain clustering we need to select a clustering resolution \(r\). Higher resolutions lead to more communities and lower resolutions lead to fewer communities. We can use the same strategy as before to find a reasonable reasonable value of \(r\).
Lets use the function get_r
for data-driven estimation of \(r\) based on
the Gap statistic and WCSS. As input we need to provide the matrix \(A\) and
a vector of \(r\)s. The output will be the Gap statistic and WCSS estimate
for each \(r\) (or the number of communities \(k'\) detected at resolution \(r\)).
b_r <- get_r(B_gap = 5,
rs = 10^seq(from = -4, to = 0, by = 0.5),
x = A,
n_start = 10,
iter_max = 50,
algorithm = "original",
knn_k = 50,
cores = 1)
Once again, the Gap statistic and WCSS curves have noticeable knees (elbows)
at \(k'=5\) (\(r=0.0025\)). Means (points) and 95% confidence intervals are shown
for the Gap statistic at each \(k\) computed using B_gap
=5 MCMC simulations.
g0_r <- (ggplot(data = b_r$gap_stats_summary)+
geom_line(aes(x = k, y = gap_mean))+
geom_point(aes(x = k, y = gap_mean), size = 1)+
geom_errorbar(aes(x = k, y = gap_mean, ymin = L95, ymax = H95), width = 0.1)+
ylab(label = "Gap")+
xlab(label = "k'")|
ggplot(data = b_r$gap_stats_summary)+
geom_line(aes(x = r, y = gap_mean))+
geom_point(aes(x = r, y = gap_mean), size = 1)+
geom_errorbar(aes(x = r, y = gap_mean, ymin = L95, ymax = H95), width = 0.1)+
ylab(label = "Gap")+
xlab(label = "r")+
scale_x_log10()+
annotation_logticks(base = 10, sides = "b"))/
(ggplot(data = b_r$wcss_stats_summary)+
geom_line(aes(x = k, y = wcss_mean))+
geom_point(aes(x = k, y = wcss_mean), size = 1)+
ylab(label = "WCSS")+
xlab(label = "k'")|
ggplot(data = b_r$wcss_stats_summary)+
geom_line(aes(x = r, y = wcss_mean))+
geom_point(aes(x = r, y = wcss_mean), size = 1)+
ylab(label = "WCSS")+
xlab(label = "r")+
scale_x_log10()+
annotation_logticks(base = 10, sides = "b"))
g0_r
A range of resolutions yields \(k=5\) number of communities, i.e. among the tested \(r\)s, we saw \(k=5\) communities for \(r = 0.003\) and \(r = 0.1\). Lets use \(r=0.1\) for clustering.
ggplot(data = b_r$gap_stats_summary)+
geom_point(aes(x = r, y = k), size = 1)+
xlab(label = "r")+
ylab(label = "k'")+
scale_x_log10()+
annotation_logticks(base = 10, sides = "b")
knitr::kable(x = b_r$gap_stats_summary[b_r$gap_stats_summary$k == 5, ],
digits = 4, row.names = FALSE)
gap_mean | r | k | gap_SE | L95 | H95 |
---|---|---|---|---|---|
2.1658 | 0.0032 | 5 | 0.0077 | 2.1507 | 2.1810 |
2.1653 | 0.0100 | 5 | 0.0022 | 2.1609 | 2.1696 |
2.1660 | 0.0316 | 5 | 0.0027 | 2.1606 | 2.1713 |
2.1662 | 0.1000 | 5 | 0.0069 | 2.1527 | 2.1797 |
Now that we found out that \(k=5\) is a reasonable choice based on the data, we will perform k-means clustering with \(k=5\) and \(A\) as inputs. For this we will use the function kmeans (R-package stats) which offers various variants of k-means. Here we will use MacQueen’s k-means variant and perform \(n_\textit{start} = 1000\) (default in scBubbletree) random starts and a maximum number of iterations \(iter_\textit{max}=300\).
Important remark: for smaller datasets (e.g. \(n<50,000\)) \(n_{start}=1000\) and \(n_{iter} = 300\) are unnecessarily high, however for larger datasets this is necessary to make sure that k-means converges.
After the clustering is complete we will organize the bubbles in a natural hierarchy. For this we perform \(B\) bootstrap iterations (default \(B=200\)). In iteration \(b\) the algorithm draws a random subset of \(N_\textit{eff}\) (default \(N_\textit{eff}=200\)) cells with replacement from each cluster and computes the average inter-cluster Euclidean distances. This data is used to populate the distance matrix (\(D^{k\times k}_b\)), which is provided as input for hierarchical clustering with average linkage to generate a hierarchical clustering dendrogram \(H_b\).
The collection of distance matrices that are computed during \(B\) iterations are used to compute a consensus (average) distance matrix (\(\hat{D}^{k\times k}\)) and from this a corresponding consensus hierarchical dendrogram (bubbletree; \(\hat{H}\)) is constructed. The collection of dendrograms are used to quantify the robustness of the bubbletree topology, i.e. to count the number of times each branch in the bubbletree is found among the topologies of the bootstrap dendrograms. Branches can have has variable degrees of support ranging between 0 (no support) and \(B\) (complete support). Distances between bubbles (inter- bubble relationships) are described quantitatively in the bubbletree as sums of branch lengths.
Steps 2.1 and 3. are performed next
k5_kmeans <- get_bubbletree_kmeans(
x = A,
k = 5,
cores = 1,
B = 200,
N_eff = 200,
round_digits = 1,
show_simple_count = FALSE,
kmeans_algorithm = "MacQueen")
… and plot the bubbletree
k5_kmeans$tree
Lets describe the bubbletree:
bubbles: The bubbletree has k=5
bubbles (clusters) shown as leaves. The
absolute and relative cell frequencies in each bubble and the bubble IDs are
shown as labels. Bubble radii scale linearly with absolute cell count in each
bubble, i.e. large bubbles have many cells and small bubbles contain few cells.
Bubble 1 is the largest one in the dendrogram and contains 1,253 cells (\(\approx\) 32% of all cells in the dataset). Bubble 4 is the smallest one and contains only 436 cells (\(\approx\) 11% of all cells in the dataset).
We can access the bubble data shown in the bubbletree
knitr::kable(k5_kmeans$tree_meta,
digits = 2, row.names = FALSE)
label | Cells | n | p | pct | lab_short | lab_long | tree_order |
---|---|---|---|---|---|---|---|
3 | 436 | 3918 | 0.11 | 11.1 | 3 (0.4K, 11.1%) | 3 (436, 11.1%) | 5 |
1 | 593 | 3918 | 0.15 | 15.1 | 1 (0.6K, 15.1%) | 1 (593, 15.1%) | 4 |
5 | 1253 | 3918 | 0.32 | 32.0 | 5 (1.3K, 32%) | 5 (1253, 32%) | 3 |
2 | 760 | 3918 | 0.19 | 19.4 | 2 (0.8K, 19.4%) | 2 (760, 19.4%) | 2 |
4 | 876 | 3918 | 0.22 | 22.4 | 4 (0.9K, 22.4%) | 4 (876, 22.4%) | 1 |
topology: inter-bubble distances are represented by sums of branch
lengths in the dendrogram. Branches of the bubbletree are annotated with
their bootstrap support values (red branch labels). The branch support
value tells us how manytimes a given branch from the bubbletree was found
among the \(B\) bootstrap dendrograms. We ran get_bubbletree_kmeans
with
\(B=200\). All but one branch have complete (200 out of 200) support, and
one branch has lower support of 179 (85%). This tells us that the branch
between bubbles (3, 4) and 1 is not as robust.
Lets also perform clustering with the Louvain algorithm (function FindClusters, R-package Seurat) and resolution parameter \(r=0.1\). There are numerous variants of the Louvain algorithm. Here we will use the original implementation. We will do clustering with \(n_\textit{start} = 20\) random starts and a maximum number of iterations \(iter_\textit{max} = 100\).
Steps 2.2 and 3. (hierarchical clustering) are performed next
k5_louvain <- get_bubbletree_graph(x = A,
r = 0.1,
n_start = 20,
iter_max = 100,
algorithm = "original",
knn_k = 50,
cores = 1,
B = 200,
N_eff = 200,
round_digits = 1,
show_simple_count = FALSE)
… and plot the bubbletree. We see nearly identical dendrogram as the one generated by kmeans clustering. The bubble IDs are different but we see similar bubble sizes, topology and branch robustness values.
k5_louvain$tree
The two dendrograms shown side-by-side:
k5_kmeans$tree|k5_louvain$tree
Given the high degree of similarity between the two clustering solutions we proceed in the next with the k-means results.
To extract biologically useful information from the bubbletree (and also for 2D UMAP or t-SNE plots) we need to adorn it with biologically relevant cell features. This includes both numeric and categorical cell features.
Numeric cell features:
Categorical cell features:
In the next two paragraph we will explain how to ‘attach’ numeric and categorical features to the bubbletree using scBubbletree.
Categorical cell features can be ‘attached’ to the bubbletree using the function
get_cat_tiles
. Here we will show the relative frequency of cell type labels
across the bubbles (parameter integrate_vertical=TRUE
).
Interpretation of the figure below:
w1 <- get_cat_tiles(btd = k5_kmeans,
f = m$cell_line_demuxlet,
integrate_vertical = TRUE,
round_digits = 1,
x_axis_name = 'Cell line',
rotate_x_axis_labels = TRUE,
tile_text_size = 2.75)
(k5_kmeans$tree|w1$plot)+
patchwork::plot_layout(widths = c(1, 1))
We can also show the inter-bubble cell type composition, i.e. the relative
frequencies of different cell types in a specific bubble (with parameter
integrate_vertical=FALSE
).
Interpretation of the figure below:
w2 <- get_cat_tiles(btd = k5_kmeans,
f = m$cell_line_demuxlet,
integrate_vertical = FALSE,
round_digits = 1,
x_axis_name = 'Cell line',
rotate_x_axis_labels = TRUE,
tile_text_size = 2.75)
(k5_kmeans$tree|w2$plot)+
patchwork::plot_layout(widths = c(1, 1))
scBubbletree uses R-package ggtree to visualize the bubbletree, and ggplot2 to visualize annotations. Furthermore, R-package patchwork is used to combine plots.
(k5_kmeans$tree|w1$plot|w2$plot)+
patchwork::plot_layout(widths = c(1, 2, 2))+
patchwork::plot_annotation(tag_levels = "A")
To quantify the purity of a cluster (or bubble) \(i\) with \(n_i\) number of cells, each of which carries one of \(L\) possible labels (e.g. cell lines), we can compute the Gini impurity index:
\(\textit{GI}_i=\sum_{j=1}^{L} \pi_{ij}(1-\pi_{ij})\),
with \(\pi_{ij}\) as the relative frequency of label \(j\) in cluster \(i\). In
homogeneous (pure
) clusters most cells carry a distinct label. Hence, the
\(\pi\)’s are close to either 1 or 0, and GI takes on a small value close to
zero. In impure
clusters cells carry a mixture of different labels. In
this case most \(\pi\) are far from either 1 or 0, and GI diverges from 0
and approaches 1. If the relative frequencies of the different labels in
cluster \(i\) are equal to the (background) relative frequencies of the labels
in the sample, then cluster \(i\) is completely impure
.
To compute the overall Gini impurity of a bubbletree, which represents a clustering solution with \(k\) bubbles, we estimated the weighted Gini impurity (WGI) by computing the weighted (by the cluster size) average of the
\(\textit{WGI}=\sum_{i=1}^{k} \textit{GI}_i \dfrac{n_i}{n}\),
with \(n_i\) as the number of cells in cluster \(i\) and \(n=\sum_i n_i\).
The Gini impurity results are shown below:
# gini
get_gini(labels = m$cell_line_demuxlet,
clusters = k5_kmeans$cluster)$gi
FALSE cluster GI
FALSE 1 1 0.020099588
FALSE 2 4 0.004558391
FALSE 3 3 0.000000000
FALSE 4 2 0.007873961
FALSE 5 5 0.000000000
All cluster-specific GIs are close to 0 and hence also WGI is close to 0.
This indicates nearly perfect mapping of cell lines to bubbles. This analysis
performed for different values of \(k\) with function get_gini_k
, which takes
as main input the output of get_k
or get_r
gini_boot <- get_gini_k(labels = m$cell_line_demuxlet,
obj = b_k)
From the figure we can conclude that WGI drops to 0 at k=5
, and all labels
are nearly perfectly split across the bubbles with each bubble containing cells
exclusively from one cell type.
g1 <- ggplot(data = gini_boot$wgi_summary)+
geom_point(aes(x = k, y = wgi), size = 0.75)+
ylab(label = "WGI")+
ylim(c(0, 1))
g1
We can also “attach” numeric cell features to the bubbletree. We will “attach” the expression of five marker genes, i.e. one marker gene for each of the five cancer cell lines.
We can visualize numeric features in two ways.
First, we can show numeric feature aggregates (e.g. “mean”, “median”, “sum”,
“pct nonzero” or “pct zero”) in the different bubbles with get_num_tiles
w3 <- get_num_tiles(btd = k5_kmeans,
fs = e,
summary_function = "mean",
x_axis_name = 'Gene expression',
rotate_x_axis_labels = TRUE,
round_digits = 1,
tile_text_size = 2.75)
(k5_kmeans$tree|w3$plot)+
patchwork::plot_layout(widths = c(1, 1))
Second, we can visualize the distributions of the numeric cell features in each
bubble as violins with get_num_violins
w4 <- get_num_violins(btd = k5_kmeans,
fs = e,
x_axis_name = 'Gene expression',
rotate_x_axis_labels = TRUE)
(k5_kmeans$tree|w3$plot|w4$plot)+
patchwork::plot_layout(widths = c(1.5, 2, 2.5))+
patchwork::plot_annotation(tag_levels = 'A')
What is the percent of UMIs coming from mitochondrial genes in each bubble?
w_mt_dist <- get_num_violins(btd = k5_kmeans,
fs = 1-m$non_mt_percent,
x_axis_name = 'MT [%]',
rotate_x_axis_labels = TRUE)
w_umi_dist <- get_num_violins(btd = k5_kmeans,
fs = m$nCount_RNA/1000,
x_axis_name = 'RNA count (in thousands)',
rotate_x_axis_labels = TRUE)
w_gene_dist <- get_num_violins(btd = k5_kmeans,
fs = m$nFeature_RNA,
x_axis_name = 'Gene count',
rotate_x_axis_labels = TRUE)
(k5_kmeans$tree|w_mt_dist$plot|w_umi_dist$plot|w_gene_dist$plot)+
patchwork::plot_layout(widths = c(0.7, 1, 1, 1))+
patchwork::plot_annotation(tag_levels = 'A')
The clustering of both algorithms seems to work well, i.e. cells from specific cell lines are clustered together. Two potential issues:
These challenges will become more severe in real datasets (e.g. PBMC samples, see case study B) which contain many clusters of cells that are not clearly separable.
Wide range of clustering approaches are used for clustering of scRNA-seq data.
scBubbletree implements the function get_bubbletree_dummy
to
allow users to incorporate results from various clustering approaches together
with our workflow. With this function we skip the k-means clustering portion of
the workflow and proceed with computing distances between the clusters and
generation of the bubbletree.
Lets try get_bubbletree_dummy
. First, will perform k-medoids clustering with
R-package cluster and then generate the bubbletree:
pam_k5 <- cluster::pam(x = A,
k = 5,
metric = "euclidean")
dummy_k5_pam <- get_bubbletree_dummy(x = A,
cs = pam_k5$clustering,
B = 200,
N_eff = 200,
cores = 2,
round_digits = 1)
dummy_k5_pam$tree|
get_cat_tiles(btd = dummy_k5_pam,
f = m$cell_line_demuxlet,
integrate_vertical = TRUE,
round_digits = 1,
tile_text_size = 2.75,
x_axis_name = 'Cell line',
rotate_x_axis_labels = TRUE)$plot
scBubbletree promotes simple and transparent visual exploration of scRNA-seq. It is not a black-box approach and the user is encouraged to explore the data with different values of \(k\) and \(r\) or custom clustering solutions. Attaching of cell features to the bubbletree is necessary for biological interpretation of the individual bubbles and their relationships which are described by the bubbletree topology.
sessionInfo()
FALSE R version 4.3.1 (2023-06-16)
FALSE Platform: x86_64-pc-linux-gnu (64-bit)
FALSE Running under: Ubuntu 22.04.3 LTS
FALSE
FALSE Matrix products: default
FALSE BLAS: /home/biocbuild/bbs-3.18-bioc/R/lib/libRblas.so
FALSE LAPACK: /usr/lib/x86_64-linux-gnu/lapack/liblapack.so.3.10.0
FALSE
FALSE locale:
FALSE [1] LC_CTYPE=en_US.UTF-8 LC_NUMERIC=C
FALSE [3] LC_TIME=en_GB LC_COLLATE=C
FALSE [5] LC_MONETARY=en_US.UTF-8 LC_MESSAGES=en_US.UTF-8
FALSE [7] LC_PAPER=en_US.UTF-8 LC_NAME=C
FALSE [9] LC_ADDRESS=C LC_TELEPHONE=C
FALSE [11] LC_MEASUREMENT=en_US.UTF-8 LC_IDENTIFICATION=C
FALSE
FALSE time zone: America/New_York
FALSE tzcode source: system (glibc)
FALSE
FALSE attached base packages:
FALSE [1] stats graphics grDevices utils datasets methods base
FALSE
FALSE other attached packages:
FALSE [1] patchwork_1.1.3 ggtree_3.10.0 ggplot2_3.4.4 scBubbletree_1.4.0
FALSE [5] BiocStyle_2.30.0
FALSE
FALSE loaded via a namespace (and not attached):
FALSE [1] deldir_1.0-9 pbapply_1.7-2 gridExtra_2.3
FALSE [4] rlang_1.1.1 magrittr_2.0.3 RcppAnnoy_0.0.21
FALSE [7] spatstat.geom_3.2-7 matrixStats_1.0.0 ggridges_0.5.4
FALSE [10] compiler_4.3.1 png_0.1-8 vctrs_0.6.4
FALSE [13] reshape2_1.4.4 stringr_1.5.0 pkgconfig_2.0.3
FALSE [16] fastmap_1.1.1 magick_2.8.1 ellipsis_0.3.2
FALSE [19] labeling_0.4.3 utf8_1.2.4 promises_1.2.1
FALSE [22] rmarkdown_2.25 purrr_1.0.2 xfun_0.40
FALSE [25] cachem_1.0.8 aplot_0.2.2 jsonlite_1.8.7
FALSE [28] goftest_1.2-3 later_1.3.1 spatstat.utils_3.0-4
FALSE [31] irlba_2.3.5.1 parallel_4.3.1 cluster_2.1.4
FALSE [34] R6_2.5.1 ica_1.0-3 spatstat.data_3.0-3
FALSE [37] bslib_0.5.1 stringi_1.7.12 RColorBrewer_1.1-3
FALSE [40] reticulate_1.34.0 parallelly_1.36.0 scattermore_1.2
FALSE [43] lmtest_0.9-40 jquerylib_0.1.4 Rcpp_1.0.11
FALSE [46] bookdown_0.36 knitr_1.44 tensor_1.5
FALSE [49] future.apply_1.11.0 zoo_1.8-12 sctransform_0.4.1
FALSE [52] httpuv_1.6.12 Matrix_1.6-1.1 splines_4.3.1
FALSE [55] igraph_1.5.1 tidyselect_1.2.0 abind_1.4-5
FALSE [58] yaml_2.3.7 spatstat.random_3.2-1 spatstat.explore_3.2-5
FALSE [61] codetools_0.2-19 miniUI_0.1.1.1 listenv_0.9.0
FALSE [64] lattice_0.22-5 tibble_3.2.1 plyr_1.8.9
FALSE [67] withr_2.5.1 shiny_1.7.5.1 treeio_1.26.0
FALSE [70] ROCR_1.0-11 evaluate_0.22 Rtsne_0.16
FALSE [73] gridGraphics_0.5-1 future_1.33.0 survival_3.5-7
FALSE [76] proxy_0.4-27 polyclip_1.10-6 fitdistrplus_1.1-11
FALSE [79] pillar_1.9.0 BiocManager_1.30.22 Seurat_4.4.0
FALSE [82] KernSmooth_2.23-22 plotly_4.10.3 ggfun_0.1.3
FALSE [85] generics_0.1.3 sp_2.1-1 munsell_0.5.0
FALSE [88] scales_1.2.1 tidytree_0.4.5 globals_0.16.2
FALSE [91] xtable_1.8-4 glue_1.6.2 lazyeval_0.2.2
FALSE [94] tools_4.3.1 data.table_1.14.8 RANN_2.6.1
FALSE [97] fs_1.6.3 leiden_0.4.3 cowplot_1.1.1
FALSE [100] grid_4.3.1 tidyr_1.3.0 ape_5.7-1
FALSE [103] colorspace_2.1-0 nlme_3.1-163 cli_3.6.1
FALSE [106] spatstat.sparse_3.0-3 fansi_1.0.5 viridisLite_0.4.2
FALSE [109] dplyr_1.1.3 uwot_0.1.16 gtable_0.3.4
FALSE [112] yulab.utils_0.1.0 sass_0.4.7 digest_0.6.33
FALSE [115] progressr_0.14.0 ggrepel_0.9.4 ggplotify_0.1.2
FALSE [118] farver_2.1.1 htmlwidgets_1.6.2 SeuratObject_4.1.4
FALSE [121] memoise_2.0.1 htmltools_0.5.6.1 lifecycle_1.0.3
FALSE [124] httr_1.4.7 mime_0.12 MASS_7.3-60