vignettes/DatasetMarkerGeneAlignment.Rmd
DatasetMarkerGeneAlignment.Rmd
In the realm of single-cell genomics, the ability to compare and integrate data across different conditions, datasets, or methodologies is crucial for deriving meaningful biological insights. This vignette introduces several functions designed to facilitate such comparisons and analyses by providing robust tools for evaluating and visualizing similarities and differences in high-dimensional data.
compareCCA()
: This function enables the comparison
of datasets by applying Canonical Correlation Analysis (CCA). It helps
assess how well two datasets align with each other, providing insights
into the relationship between different single-cell experiments or
conditions.
comparePCA()
: This function allows you to compare
datasets using Principal Component Analysis (PCA). It evaluates how
similar or different the principal components are between two datasets,
offering a way to understand the underlying structure and variance in
your data.
comparePCASubspace()
: Extending the comparison to
specific subspaces, this function focuses on subsets of principal
components. It provides a detailed analysis of how subspace structures
differ or align between datasets, which is valuable for fine-grained
comparative studies.
plotPairwiseDistancesDensity()
: To visualize the
distribution of distances between pairs of samples, this function
generates density plots. It helps in understanding the variation and
relationships between samples in high-dimensional spaces.
plotWassersteinDistance()
: This function visualizes
the Wasserstein distance, a metric for comparing distributions, across
datasets. It provides an intuitive view of how distributions differ
between datasets, aiding in the evaluation of alignment and
discrepancies.
calculateCramerPValue()
: This function computes the
Cramér’s V statistic, which measures the strength of association between
categorical variables. Cramér’s V is particularly useful for quantifying
the association strength in contingency tables, providing insight into
how strongly different cell types or conditions are related.
calculateHotellingPValue()
: This function performs
Hotelling’s T-squared test, a multivariate statistical test used to
assess the differences between two groups of observations. It is useful
for comparing the mean vectors of multiple variables and is commonly
employed in single-cell data to evaluate group differences in principal
component space.
calculateNearestNeighborProbabilities()
: This
function calculates the probability of nearest neighbor assignments in a
high-dimensional space. It is used to evaluate the likelihood of a cell
belonging to a specific cluster based on its proximity to neighboring
cells. This function helps in assessing clustering quality and the
robustness of cell type assignments.
calculateAveragePairwiseCorrelation()
: This function
computes the average pairwise correlation between variables or features
across cells. It provides an overview of the relationships between
different markers or genes, helping to identify correlated features and
understand their interactions within the dataset.
regressPC()
: This function performs linear
regression of a covariate of interest onto one or more principal
components derived from single-cell data. This analysis helps quantify
the variance explained by the covariate and is useful for tasks such as
assessing batch effects, clustering homogeneity, and alignment of query
and reference datasets.
calculateHVGOverlap()
: To assess the similarity
between datasets based on highly variable genes, this function computes
the overlap coefficient. It measures how well the sets of highly
variable genes from different datasets correspond to each
other.
calculateVarImpOverlap()
: Using Random Forest, this
function identifies and compares the importance of genes for
differentiating cell types between datasets. It highlights which genes
are most critical in each dataset and compares their importance,
providing insights into shared and unique markers.
These functions collectively offer a comprehensive toolkit for comparing and analyzing single-cell data. Whether you are assessing alignment between datasets, visualizing distance distributions, or identifying key genes, these tools are designed to enhance your ability to derive meaningful insights from complex, high-dimensional data.
In this vignette, we will guide you through the practical use of each function, demonstrate how to interpret their outputs, and show how they can be integrated into your single-cell genomics workflow.
In the context of the scDiagnostics
package, this
vignette demonstrates how to leverage various functions to evaluate and
compare single-cell data across two distinct datasets:
reference_data
: This dataset features meticulously
curated cell type annotations assigned by experts. It serves as a robust
benchmark for evaluating the accuracy and consistency of cell type
annotations across different datasets, offering a reliable standard
against which other annotations can be assessed.query_data
: This dataset contains cell type annotations
from both expert assessments and those generated using the SingleR
package. By comparing these annotations with those from the reference
dataset, you can identify discrepancies between manual and automated
results, highlighting potential inconsistencies or areas requiring
further investigation.
# Load library
library(scDiagnostics)
# Load datasets
data("reference_data")
data("query_data")
# Set seed for reproducibility
set.seed(0)
Some functions in the vignette are designed to work with SingleCellExperiment objects that contain data from only one cell type. We will create separate SingleCellExperiment objects that only CD4 cells, to ensure compatibility with these functions.
# Load library
library(scran)
library(scater)
# Subset to CD4 cells
ref_data_cd4 <- reference_data[, which(
reference_data$expert_annotation == "CD4")]
query_data_cd4 <- query_data_cd4 <- query_data[, which(
query_data$expert_annotation == "CD4")]
# Select highly variable genes
ref_top_genes <- getTopHVGs(ref_data_cd4, n = 500)
query_top_genes <- getTopHVGs(query_data_cd4, n = 500)
common_genes <- intersect(ref_top_genes, query_top_genes)
# Subset data by common genes
ref_data_cd4 <- ref_data_cd4[common_genes,]
query_data_cd4 <- query_data_cd4[common_genes,]
# Run PCA on both datasets
ref_data_cd4 <- runPCA(ref_data_cd4)
query_data_cd4 <- runPCA(query_data_cd4)
compareCCA()
In single-cell genomics, datasets from different sources or
experimental conditions often need to be compared to understand how well
they align. The compareCCA
function facilitates this
comparison by performing Canonical Correlation Analysis (CCA) between a
query dataset and a reference dataset. This function is particularly
useful when datasets have different sample sizes or distributions, and
it helps in assessing the similarity between them after projecting them
onto a common principal component space.
compareCCA()
performs the following steps:
# Perform CCA
cca_comparison <- compareCCA(query_data = query_data_cd4,
reference_data = ref_data_cd4,
query_cell_type_col = "expert_annotation",
ref_cell_type_col = "expert_annotation",
pc_subset = 1:5)
plot(cca_comparison)
Canonical Correlation Analysis (CCA) produces several key outputs:
In summary, the compareCCA
function allows you to
compare how well two datasets are aligned by projecting them into a
shared PCA space, performing CCA, and then evaluating the similarity of
the canonical variables. This approach is valuable for integrative
analyses and understanding the relationships between different datasets
in single-cell studies.
comparePCA()
The comparePCA()
function compares the PCA subspaces
between the query and reference datasets. It calculates the principal
angles between the PCA subspaces to assess the alignment and similarity
between them. This is useful for understanding how well the
dimensionality reduction spaces of your datasets match.
comparePCA()
operates as follows:
# Perform PCA
pca_comparison <- comparePCA(query_data = query_data_cd4,
reference_data = ref_data_cd4,
query_cell_type_col = "expert_annotation",
ref_cell_type_col = "expert_annotation",
pc_subset = 1:5,
metric = "cosine")
plot(pca_comparison)
comparePCASubspace()
In single-cell RNA-seq analysis, it is essential to assess the
similarity between the subspaces spanned by the top principal components
(PCs) of query and reference datasets. This is particularly important
when comparing the structure and variation captured by each dataset. The
comparePCASubspace()
function is designed to provide
insights into how well the subspaces align by computing the cosine
similarity between the loadings of the top variables for each PC. This
analysis helps in determining the degree of correspondence between the
datasets, which is critical for accurate cell type annotation and data
integration.
comparePCASubspace()
performs the following
operations:
# Compare PCA subspaces between query and reference data
subspace_comparison <- comparePCASubspace(
query_data = query_data_cd4,
reference_data = ref_data_cd4,
query_cell_type_col = "expert_annotation",
ref_cell_type_col = "expert_annotation",
pc_subset = 1:5
)
# View weighted cosine similarity score
subspace_comparison$weighted_cosine_similarity
#> [1] 0.2578375
# Plot output for PCA subspace comparison (if a plot method is available)
plot(subspace_comparison)
In the results:
By using comparePCASubspace()
, you can quantify the
alignment of PCA subspaces between query and reference datasets, aiding
in the evaluation of dataset integration and the reliability of cell
type annotations.
plotPairwiseDistancesDensity()
The plotPairwiseDistancesDensity()
function is designed
to calculate and visualize the pairwise distances or correlations
between cell types in query and reference datasets. This function is
particularly useful in single-cell RNA sequencing (scRNA-seq) analysis,
where it is essential to evaluate the consistency and reliability of
cell type annotations by comparing their expression profiles.
The function operates on SingleCellExperiment objects, which are commonly used to store single-cell data, including expression matrices and associated metadata. Users specify the cell types of interest in both the query and reference datasets, and the function computes either the distances or correlation coefficients between these cells.
When principal component analysis (PCA) is applied, the function projects the expression data into a lower-dimensional PCA space, which can be specified by the user. This allows for a more focused analysis of the major sources of variation in the data. Alternatively, if no dimensionality reduction is desired, the function can directly use the expression data for computation.
Depending on the user’s choice, the function can calculate pairwise Euclidean distances or correlation coefficients. The resulting values are used to compare the relationships between cells within the same dataset (either query or reference) and between cells across the two datasets.
The output of the function is a density plot generated by
ggplot2
, which displays the distribution of pairwise
distances or correlations. The plot provides three key comparisons:
By examining these density curves, users can assess the similarity of cells within each dataset and across datasets. For example, a higher density of lower distances in the “Query vs. Reference” comparison would suggest that the query and reference cells are similar in their expression profiles, indicating consistency in the annotation of the cell type across the datasets.
This visual approach offers an intuitive way to diagnose potential discrepancies in cell type annotations, identify outliers, or confirm the reliability of the cell type assignments.
# Example usage of the function
plotPairwiseDistancesDensity(query_data = query_data,
reference_data = reference_data,
query_cell_type_col = "expert_annotation",
ref_cell_type_col = "expert_annotation",
cell_type_query = "CD4",
cell_type_ref = "CD4",
pc_subset = 1:5,
distance_metric = "euclidean")
This example demonstrates how to compare CD4 cells between a query and reference dataset, with PCA applied to the first five principal components and pairwise Euclidean distances calculated. The output is a density plot that helps visualize the distribution of these distances, aiding in the interpretation of the similarity between the two datasets.
calculateWassersteinDistance()
the calculateWassersteinDistance()
function uses the
concept of Wasserstein distances, a measure rooted in optimal transport
theory, to quantitatively compare these datasets. By projecting both
datasets into a common PCA space, the function enables a meaningful
comparison even when they originate from different experimental
conditions.
Before diving into a code example, it’s important to understand the output of this function. The result is a list containing three key components: the null distribution of Wasserstein distances computed from the reference dataset alone, the mean Wasserstein distance between the query and reference datasets, and a vector of the unique cell types present in the reference dataset. This comparison provides insights into how similar or different the query dataset is relative to the reference dataset.
What the function reveals is particularly valuable: if the mean Wasserstein distance between the query and reference is significantly different from the null distribution, it indicates a notable variation between the datasets, reflecting differences in cell-type compositions, technical artifacts, or biological questions of interest.
# Generate the Wasserstein distance density plot
wasserstein_data <- calculateWassersteinDistance(
query_data = query_data_cd4,
reference_data = ref_data_cd4,
query_cell_type_col = "expert_annotation",
ref_cell_type_col = "expert_annotation",
pc_subset = 1:10,
)
plot(wasserstein_data)
The plotting function generates a density plot that provides a clear graphical representation of the Wasserstein distances. On the plot, you’ll see the distribution of Wasserstein distances calculated within the reference dataset, referred to as the null distribution.
Within this density plot, two significant vertical lines are overlayed: one representing the significance threshold and another indicating the mean Wasserstein distance between the reference and query datasets. The significance threshold is determined based on your specified alpha level, which is the probability threshold for considering differences between the datasets as statistically significant. This threshold line helps to establish a boundary; if the observed reference-query distance exceeds this threshold, it suggests that the differences between your query and reference datasets are unlikely to occur by random chance.
The line for the reference-query distance helps identify where this calculated distance lies in comparison to the null distribution. If this dashed line is noticeably separate from the bulk of the null distribution and beyond the significance threshold, it implies a significant deviation between the query dataset and the reference dataset. In practical terms, this could indicate distinct biological differences, batch effects, or other analytical artefacts that differentiate the datasets.
calculateCramerPValue()
The calculateCramerPValue function is designed to perform the Cramer test for comparing multivariate empirical cumulative distribution functions (ECDFs) between two samples in single-cell data. This test is particularly useful for assessing whether the distributions of principal components (PCs) differ significantly between the reference and query datasets for specific cell types.
To use this function, you first need to provide two key inputs:
reference_data
and query_data
, both of which
should be SingleCellExperiment
objects containing numeric expression matrices. If
query_data
is not supplied, the function will only use the
reference_data.
You should also specify the column names
for cell type annotations in both datasets via
ref_cell_type_col
and query_cell_type_col.
If
cell_types
is not provided, the function will automatically
include all unique cell types found in the datasets. The
pc_subset
parameter allows you to define which principal
components to include in the analysis, with the default being the first
five PCs.
The function performs the following steps: it first projects the data into PCA space, subsets the data according to the specified cell types and principal components, and then applies the Cramer test to compare the ECDFs between the reference and query datasets. The result is a named vector of p-values from the Cramer test for each cell type, which indicates whether there is a significant difference in the distributions of PCs between the two datasets.
Here’s an example of how to use the
calculateCramerPValue()
function:
# Perform Cramer test for the specified cell types and principal components
cramer_test <- calculateCramerPValue(
reference_data = reference_data,
query_data = query_data,
ref_cell_type_col = "expert_annotation",
query_cell_type_col = "SingleR_annotation",
cell_types = c("CD4", "CD8", "B_and_plasma", "Myeloid"),
pc_subset = 1:5
)
# Display the Cramer test p-values
print(cramer_test)
#> CD4 CD8 B_and_plasma Myeloid
#> 0.0000000 0.0000000 0.1378621 0.4575425
In this example, the function compares the distributions of the first five principal components between the reference and query datasets for cell types such as “CD4”, “CD8”, “B_and_plasma”, and “Myeloid”. The output is a vector of p-values indicating whether the differences in PC distributions are statistically significant.
calculateHotellingPValue()
The calculateHotellingPValue()
function is designed to
compute Hotelling’s T-squared test statistic and corresponding p-values
for comparing multivariate means between reference and query datasets in
the context of single-cell RNA-seq data. This statistical test is
particularly useful for assessing whether the mean vectors of principal
components (PCs) differ significantly between the two datasets, which
can be indicative of differences in the cell type distributions.
To use this function, you need to provide two SingleCellExperiment
objects: query_data
and reference_data
, each
containing numeric expression matrices. You also need to specify the
column names for cell type annotations in both datasets with
query_cell_type_col
and ref_cell_type_col.
The
cell_types
parameter allows you to choose which cell types
to include in the analysis, and if not specified, the function will
automatically include all cell types present in the datasets. The
pc_subset
parameter determines which principal components
to consider, with the default being the first five PCs. Additionally,
n_permutation
specifies the number of permutations for
calculating p-values, with a default of 500.
The function works by first projecting the data into PCA space and then performing Hotelling’s T-squared test for each specified cell type to compare the means between the reference and query datasets. It uses permutation testing to determine the p-values, indicating whether the observed differences in means are statistically significant. The result is a named numeric vector of p-values for each cell type.
Here is an example of how to use the
calculateHotellingPValue()
function:
# Calculate Hotelling's T-squared test p-values
p_values <- calculateHotellingPValue(
query_data = query_data,
reference_data = reference_data,
query_cell_type_col = "SingleR_annotation",
ref_cell_type_col = "expert_annotation",
pc_subset = 1:10
)
# Display the p-values rounded to five decimal places
print(round(p_values, 5))
#> CD4 CD8 B_and_plasma Myeloid
#> 0.00 0.00 0.09 0.00
In this example, the function calculates p-values for Hotelling’s T-squared test, focusing on the first ten principal components to assess whether the multivariate means differ significantly between the reference and query datasets for each cell type. The p-values indicate the likelihood of observing the differences by chance and help identify significant disparities in cell type distributions.
calculateAveragePairwiseCorrelation()
The calculateAveragePairwiseCorrelation()
function is
designed to compute the average pairwise correlations between specified
cell types in single-cell gene expression data. This function operates
on SingleCellExperiment
objects and is ideal for single-cell analysis workflows. It calculates
pairwise correlations between query and reference cells using a
specified correlation method, and then averages these correlations for
each cell type pair. This helps in assessing the similarity between
cells in the reference and query datasets and provides insights into the
reliability of cell type annotations.
To use the calculateAveragePairwiseCorrelation()
function, you need to supply it with two SingleCellExperiment
objects: one for the query cells and one for the reference cells. The
function also requires column names specifying cell type annotations in
both datasets, and optionally a vector of cell types to include in the
analysis. Additionally, you can specify a subset of principal components
to use, or use the raw data directly if pc_subset
is set to
NULL
. The correlation method can be either “spearman” or
“pearson”.
Here’s an example of how to use this function:
# Compute pairwise correlations between specified cell types
cor_matrix_avg <- calculateAveragePairwiseCorrelation(
query_data = query_data,
reference_data = reference_data,
query_cell_type_col = "SingleR_annotation",
ref_cell_type_col = "expert_annotation",
cell_types = c("CD4", "CD8", "B_and_plasma"),
pc_subset = 1:5,
correlation_method = "spearman"
)
# Visualize the average pairwise correlation matrix
plot(cor_matrix_avg)
In this example, calculateAveragePairwiseCorrelation()
computes the average pairwise correlations for the cell types “CD4”,
“CD8”, and “B_and_plasma”. This function is particularly useful for
understanding the relationships between different cell types in your
single-cell datasets and evaluating how well cell types in the query
data correspond to those in the reference data.
regressPC()
The regressPC()
function performs linear regression of a
covariate of interest onto one or more principal components using data
from a SingleCellExperiment
object. This method helps quantify the variance explained by a
covariate, which can be useful in applications such as quantifying batch
effects, assessing clustering homogeneity, and evaluating alignment
between query and reference datasets in cell type annotation
settings.
The function calculates the R-squared value from the linear regression of the covariate onto each principal component. The variance contribution of the covariate effect per principal component is computed as the product of the variance explained by the principal component and the R-squared value. The total variance explained by the covariate is obtained by summing these contributions across all principal components.
Here is how you can use the regressPC()
function:
# Perform regression analysis using only reference data
regress_res <- regressPC(
reference_data = reference_data,
ref_cell_type_col = "expert_annotation",
cell_types = c("CD4", "CD8", "B_and_plasma", "Myeloid"),
pc_subset = 1:15
)
# Plot results showing R-squared values
plot(regress_res, plot_type = "r_squared")
# Plot results showing p-values
plot(regress_res, plot_type = "p-value")
In this example, regressPC()
is used to perform
regression analysis on principal components 1 to 15 from the reference
dataset. The results are then visualized using plots showing the
R-squared values and p-values.
If you also have a query dataset and want to compare it with the reference dataset, you can include it in the analysis:
# Perform regression analysis using both reference and query data
regress_res <- regressPC(
reference_data = reference_data,
query_data = query_data,
ref_cell_type_col = "expert_annotation",
query_cell_type_col = "SingleR_annotation",
cell_types = c("CD4", "CD8", "B_and_plasma", "Myeloid"),
pc_subset = 1:15
)
# Plot results showing R-squared values
plot(regress_res, plot_type = "r_squared")
# Plot results showing p-values
plot(regress_res, plot_type = "p-value")
In this case, the function projects the query data onto the principal components of the reference data and performs regression analysis. The results are visualized in the same way, providing insights into how well the principal components from the query data align with those from the reference data.
This function is useful for understanding the impact of different covariates on principal components and for comparing how cell types in different datasets are represented in the principal component space.
calculateHVGOverlap()
The calculateHVGOverlap()
function computes the overlap
coefficient between two sets of highly variable genes (HVGs) from a
reference dataset and a query dataset. The overlap coefficient is a
measure of similarity between the two sets, reflecting how much the HVGs
in one dataset overlap with those in the other.
The function begins by ensuring that the input vectors
reference_genes
and query_genes
are character
vectors and that neither of them is empty. Once these checks are
complete, the function identifies the common genes between the two sets
using the intersect function, which finds the intersection of the two
gene sets.
Next, the function calculates the size of this intersection, representing the number of genes common to both sets. The overlap coefficient is then computed by dividing the size of the intersection by the size of the smaller set of genes. This ensures that the coefficient is a value between 0 and 1, where 0 indicates no overlap and 1 indicates complete overlap.
Finally, the function rounds the overlap coefficient to two decimal places before returning it as the output.
The overlap coefficient quantifies the extent to which the HVGs in the reference dataset align with those in the query dataset. A higher coefficient indicates a stronger similarity between the two datasets in terms of their HVGs, which can suggest that the datasets are more comparable or that they capture similar biological variability. Conversely, a lower coefficient indicates less overlap, suggesting that the datasets may be capturing different biological signals or that they are less comparable.
# Load library to get top HVGs
library(scran)
# Select the top 500 highly variable genes from both datasets
ref_var_genes <- getTopHVGs(reference_data, n = 500)
query_var_genes <- getTopHVGs(query_data, n = 500)
# Calculate the overlap coefficient between the reference and query HVGs
overlap_coefficient <- calculateHVGOverlap(reference_genes = ref_var_genes,
query_genes = query_var_genes)
# Display the overlap coefficient
overlap_coefficient
#> [1] 0.93
calculateVarImpOverlap()
The calculateVarImpOverlap()
function helps you identify
and compare the most important genes for distinguishing cell types in
both a reference dataset and a query dataset. It does this using the
Random Forest algorithm, which calculates how important each gene is in
differentiating between cell types.
To use the function, you need to provide a reference dataset containing expression data and cell type annotations. Optionally, you can also provide a query dataset if you want to compare gene importance across both datasets. The function allows you to specify which cell types to analyze and how many trees to use in the Random Forest model. Additionally, you can decide how many top genes you want to compare between the datasets.
Let’s say you have a reference dataset (reference_data
)
and a query dataset (query_data
). Both datasets contain
expression data and cell type annotations, stored in columns named
“expert_annotation” and SingleR_annotation
, respectively.
You want to calculate the importance of genes using 500 trees and
compare the top 50 genes between the datasets.
Here’s how you would use the function:
# RF function to compare which genes are best at differentiating cell types
rf_output <- calculateVarImpOverlap(reference_data = reference_data,
query_data = query_data,
query_cell_type_col = "expert_annotation",
ref_cell_type_col = "expert_annotation",
n_tree = 500,
n_top = 50)
# Comparison table
rf_output$var_imp_comparison
#> CD4-CD8 CD4-B_and_plasma CD4-Myeloid
#> 0.84 0.84 0.88
#> CD8-B_and_plasma CD8-Myeloid B_and_plasma-Myeloid
#> 0.84 0.80 0.76
After running the function, you’ll receive the importance scores of genes for each pair of cell types in your reference dataset. If you provided a query dataset, you’ll also get the importance scores for those cell types. The function will tell you how much the top genes in the reference and query datasets overlap, which helps you understand if the same genes are important for distinguishing cell types across different datasets.
For example, if there’s a high overlap, it suggests that similar genes are crucial in both datasets for differentiating the cell types, which could be important for validating your findings or identifying robust markers.
R version 4.4.1 (2024-06-14)
Platform: x86_64-pc-linux-gnu
Running under: Ubuntu 22.04.5 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.20.so; LAPACK version 3.10.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] scater_1.32.1 ggplot2_3.5.1
[3] scran_1.32.0 scuttle_1.14.0
[5] SingleCellExperiment_1.26.0 SummarizedExperiment_1.34.0
[7] Biobase_2.64.0 GenomicRanges_1.56.2
[9] GenomeInfoDb_1.40.1 IRanges_2.38.1
[11] S4Vectors_0.42.1 BiocGenerics_0.50.0
[13] MatrixGenerics_1.16.0 matrixStats_1.4.1
[15] scDiagnostics_1.1.0 BiocStyle_2.32.1
loaded via a namespace (and not attached):
[1] DBI_1.2.3 gridExtra_2.3
[3] cramer_0.9-4 rlang_1.1.4
[5] magrittr_2.0.3 ggridges_0.5.6
[7] compiler_4.4.1 DelayedMatrixStats_1.26.0
[9] systemfonts_1.1.0 vctrs_0.6.5
[11] pkgconfig_2.0.3 crayon_1.5.3
[13] fastmap_1.2.0 XVector_0.44.0
[15] labeling_0.4.3 utf8_1.2.4
[17] rmarkdown_2.28 UCSC.utils_1.0.0
[19] ggbeeswarm_0.7.2 ragg_1.3.3
[21] xfun_0.48 bluster_1.14.0
[23] zlibbioc_1.50.0 cachem_1.1.0
[25] beachmat_2.20.0 jsonlite_1.8.9
[27] DelayedArray_0.30.1 BiocParallel_1.38.0
[29] irlba_2.3.5.1 parallel_4.4.1
[31] cluster_2.1.6 biglm_0.9-3
[33] R6_2.5.1 bslib_0.8.0
[35] ranger_0.16.0 limma_3.60.6
[37] boot_1.3-30 jquerylib_0.1.4
[39] Rcpp_1.0.13 bookdown_0.41
[41] knitr_1.48 Matrix_1.7-0
[43] igraph_2.1.1 tidyselect_1.2.1
[45] abind_1.4-8 yaml_2.3.10
[47] viridis_0.6.5 codetools_0.2-20
[49] lattice_0.22-6 tibble_3.2.1
[51] withr_3.0.2 evaluate_1.0.1
[53] desc_1.4.3 pillar_1.9.0
[55] BiocManager_1.30.25 generics_0.1.3
[57] sparseMatrixStats_1.16.0 munsell_0.5.1
[59] scales_1.3.0 glue_1.8.0
[61] metapod_1.12.0 tools_4.4.1
[63] speedglm_0.3-5 data.table_1.16.2
[65] BiocNeighbors_1.22.0 ScaledMatrix_1.12.0
[67] locfit_1.5-9.10 fs_1.6.4
[69] grid_4.4.1 edgeR_4.2.2
[71] colorspace_2.1-1 GenomeInfoDbData_1.2.12
[73] beeswarm_0.4.0 BiocSingular_1.20.0
[75] vipor_0.4.7 cli_3.6.3
[77] rsvd_1.0.5 textshaping_0.4.0
[79] fansi_1.0.6 S4Arrays_1.4.1
[81] viridisLite_0.4.2 dplyr_1.1.4
[83] gtable_0.3.6 sass_0.4.9
[85] digest_0.6.37 SparseArray_1.4.8
[87] ggrepel_0.9.6 dqrng_0.4.1
[89] htmlwidgets_1.6.4 farver_2.1.2
[91] htmltools_0.5.8.1 pkgdown_2.1.1
[93] lifecycle_1.0.4 httr_1.4.7
[95] transport_0.15-4 statmod_1.5.0
[97] MASS_7.3-60.2