scGraphVerse is a modular and extensible R package for constructing, comparing, and visualizing gene regulatory networks (GRNs) from single-cell RNAseq data. It includes several inference algorithms, utilities, and visualization functions.
Single-cell data presents unique challenges: high sparsity, batch variation, and the need to distinguish shared versus condition-specific regulation. scGraphVerse helps address these through multi-method inference consensus. It starts with SingleCellExperiment, Seurat or matrix-based objects, enabling flexible pipelines across diverse experimental setups.
While several GRN inference packages exist in Bioconductor (e.g. SCENIC), most are limited to single-method inference or lack support for multi-condition, multi-replicate, and multi-method comparisons. scGraphVerse adds value by:
The typical scGraphVerse workflow begins with one or more count matrix, or objects like SingleCellExperiment or Seurat. These are passed to the core function infer_networks(), which runs the selected inference methods and returns an adjacency matrix or list of matrices.
For benchmarking or teaching purposes, synthetic datasets can be generated using zinb_simdata() based on a known adjacency matrix. This supports method comparison across early, late, and joint integration strategies using earlyj() and compare_consensus().
An end-to-end workflow, from network inference to community detection and consensus visualization, is provided in the rest of vignettes.
As intended for both high-throughput analysis and intuitive usage, scGraphVerse supports custom workflows using its modular functions.
In this simulation study, we use scGraphVerse to:
We use a pre-computed toy adjacency matrix (toy_adj_matrix
) as our ground
truth network. This matrix represents known regulatory relationships between
genes and serves as the reference for evaluating network inference performance.
# Load the toy ground truth adjacency matrix
data(toy_adj_matrix)
adj_truth <- toy_adj_matrix
# Visualize network
if (requireNamespace("igraph", quietly = TRUE) &&
requireNamespace("ggraph", quietly = TRUE)) {
gtruth <- igraph::graph_from_adjacency_matrix(adj_truth,mode="undirected")
ggraph::ggraph(gtruth, layout = "fr") +
ggraph::geom_edge_link(color = "gray") +
ggraph::geom_node_point(color = "steelblue") +
ggplot2::ggtitle(paste0(
"Ground Truth: ",
igraph::vcount(gtruth),
" nodes, ",
igraph::ecount(gtruth),
" edges"
)) +
ggplot2::theme_minimal()
}
Now we’ll generate synthetic single-cell count data that respects our
ground-truth network structure. We simulate three experimental conditions
(batches) with 50 cells each, using different parameters to model realistic
batch effects and biological variability. The zinb_simdata()
function creates
zero-inflated negative binomial count data where gene expression levels are
influenced by the regulatory relationships defined in our ground-truth network.
Key simulation parameters:
- mu_range
: Different mean expression levels across conditions
- theta
: Overdispersion parameter (lower = more variable)
- pi
: Zero-inflation probability (dropout rate)
- depth_range
: Sequencing depth variation per cell
# Simulation parameters
nodes <- nrow(adj_truth)
sims <- zinb_simdata(
n = 50,
p = nodes,
B = adj_truth,
mu_range = list(c(1, 4), c(1, 7)),
mu_noise = c(1, 3),
theta = c(1, 0.7),
pi = c(0.2, 0.2),
kmat = 2,
depth_range = c(0.8 * nodes * 3, 1.2 * nodes * 3)
)
# Transpose to cells × genes
count_matrices <- lapply(sims, t)
Next, we apply the GENIE3 algorithm to infer gene regulatory networks from our simulated count data. GENIE3 uses random forests to predict each gene’s expression based on all other genes, then ranks the importance of potential regulatory relationships. We process all three batches to capture condition-specific and shared regulatory patterns.
Note: The package now uses Bioconductor data structures:
- count_matrices
is now a MultiAssayExperiment
(MAE) object
- Network outputs are stored in SummarizedExperiment
(SE) objects
GENIE3 workflow: 1. Create a MAE object from the count matrices list 2. For each target gene, build random forest model using all other genes as predictors 3. Extract feature importance scores as edge weights 4. Generate weighted adjacency matrices stored in a SummarizedExperiment 5. Symmetrize the matrices using mean weights to create undirected networks
# Create MultiAssayExperiment from count matrices
mae <- create_mae(count_matrices)
# Infer networks
networks_joint <- infer_networks(
count_matrices_list = mae,
method = "GENIE3",
nCores = 1
)
# Generate weighted adjacency matrices (returns SummarizedExperiment)
wadj_se <- generate_adjacency(networks_joint)
# Symmetrize weights (returns SummarizedExperiment)
swadj_se <- symmetrize(wadj_se, weight_function = "mean")
Now we evaluate how well our inferred networks recover the ground-truth regulatory relationships using ROC curve analysis. The ROC curve plots True Positive Rate (sensitivity) against False Positive Rate (1-specificity) across different edge weight thresholds. A perfect classifier would achieve AUC = 1.0, while random guessing gives AUC = 0.5.
ROC Analysis Steps: 1. Use continuous edge weights from symmetrized adjacency matrices 2. Compare against binary ground-truth network 3. Calculate Area Under Curve (AUC) as overall performance metric 4. Higher AUC indicates better recovery of true regulatory edges
if (requireNamespace("pROC", quietly = TRUE)) {
roc_res <- plotROC(
swadj_se,
adj_truth,
plot_title = "ROC Curve: GENIE3 Network Inference",
is_binary = FALSE
)
roc_res$plot
auc_joint <- roc_res$auc
}
To convert our continuous edge weights into binary networks, we need to
determine appropriate cutoff thresholds. The cutoff_adjacency()
function uses
a data-driven approach: it creates shuffled versions of our count data, infers
networks from these null datasets, and sets the threshold at the 95th
percentile of shuffled edge weights. This ensures we keep only edges that are
much stronger than expected by random chance.
Cutoff Method: 1. Generate shuffled count matrices (randomize gene expression profiles) 2. Run GENIE3 on shuffled data to get null distribution of edge weights 3. Set threshold at 95th percentile of shuffled weights 4. Apply threshold to original networks to get binary adjacency matrices 5. Calculate precision metrics: TPR (sensitivity), FPR, Precision, F1-score
# Binary cutoff at 95th percentile (returns SummarizedExperiment)
binary_se <- cutoff_adjacency(
count_matrices = mae,
weighted_adjm_list = swadj_se,
n = 1,
method = "GENIE3",
quantile_threshold = 0.95,
nCores = 1,
debug = TRUE
)
#> [Method: GENIE3] Matrix 1 → Cutoff = 0.06060
#> [Method: GENIE3] Matrix 2 → Cutoff = 0.05617
# Precision scores
pscores_joint <- pscores(adj_truth, binary_se)
#> Registered S3 methods overwritten by 'fmsb':
#> method from
#> print.roc pROC
#> plot.roc pROC
head(pscores_joint)
#> $Statistics
#> Predicted_Matrix TP TN FP FN TPR FPR Precision F1
#> 1 Matrix 1 14 537 27 17 0.4516129 0.04787234 0.3414634 0.3888889
#> 2 Matrix 2 6 530 34 25 0.1935484 0.06028369 0.1500000 0.1690141
#> MCC
#> 1 0.3542224
#> 2 0.1182658
#>
#> $Radar
#> $Radar$data
#> TPR Specificity Precision F1 MCC
#> Max 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000
#> Min 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000
#> Matrix 1 0.4516129 0.9521277 0.3414634 0.3888889 0.3542224
#> Matrix 2 0.1935484 0.9397163 0.1500000 0.1690141 0.1182658
#>
#> $Radar$plot
# Network plot
if (requireNamespace("ggraph", quietly = TRUE)) {
plotg(binary_se)
}
Since we have three binary networks (one per experimental condition), we can create a consensus network that captures edges consistently detected across conditions. The “vote” method includes an edge in the consensus only if it appears in at least 2 out of 3 individual networks, making it more robust to condition-specific noise.
Consensus Building: - Vote method: Edge included if present in majority of networks (≥2/3) - Union method: Edge included if present in any network (≥1/3) - INet method: Uses weighted evidence combination (more sophisticated)
# Consensus matrix (returns SummarizedExperiment)
consensus <- create_consensus(binary_se, method = "vote")
# Plot consensus network
if (requireNamespace("ggraph", quietly = TRUE)) {
plotg(consensus)
}
# Compare consensus to truth using classify_edges and plot_network_comparison
if (requireNamespace("ggplot2", quietly = TRUE) &&
requireNamespace("ggraph", quietly = TRUE)) {
# Wrap toy_adj_matrix in a SummarizedExperiment for classify_edges
adj_truth_se <- build_network_se(list(ground_truth = adj_truth))
# classify_edges expects SummarizedExperiment objects
edge_classification <- classify_edges(
consensus_matrix = consensus,
reference_matrix = adj_truth_se
)
print(edge_classification)
# Visualize network comparison (show_fp = FALSE to hide false positives)
comparison_plot <- plot_network_comparison(
edge_classification = edge_classification,
show_fp = FALSE
)
print(comparison_plot)
}
#> $tp_edges
#> [1] "ACTG1-PFN1" "ACTG1-TMSB4X" "BTF3-NACA" "CD3D-CD3E" "CD3E-HLA-A"
#> [6] "CD74-CXCR4" "EEF2-GNB2L1" "EEF2-UBA52" "FOS-JUN" "FOS-JUNB"
#> [11] "FTH1-FTL" "HLA-A-HLA-B" "HLA-A-HLA-E" "HLA-B-HLA-C" "HLA-C-HLA-E"
#> [16] "JUN-JUNB" "UBA52-UBC"
#>
#> $fp_edges
#> [1] "ACTG1-BTF3" "ACTG1-COX4I1" "ACTG1-EEF1A1" "ACTG1-HLA-E"
#> [5] "ARPC2-CD3D" "ARPC2-TMSB4X" "ARPC2-UBC" "ARPC3-CXCR4"
#> [9] "ARPC3-EIF4A2" "ARPC3-HLA-A" "ARPC3-JUN" "BTF3-CFL1"
#> [13] "BTF3-FOS" "CD3D-CD74" "CD3D-MYL12B" "CD3D-PABPC1"
#> [17] "CD3D-TMSB4X" "CD3D-UBC" "CD3E-COX7C" "CD3E-EIF4A2"
#> [21] "CD3E-TMSB4X" "CD74-EEF1D" "CD74-UBC" "CFL1-COX4I1"
#> [25] "CFL1-EEF1D" "COX4I1-HLA-E" "COX7C-EEF1A1" "COX7C-TMSB4X"
#> [29] "CXCR4-FTL" "EEF1A1-EEF2" "EEF1A1-FTH1" "EEF1A1-JUN"
#> [33] "EEF1A1-NACA" "EEF1A1-PFN1" "EEF1D-HLA-A" "EEF1D-UBC"
#> [37] "EEF2-EIF3K" "EIF1-MYL12B" "EIF1-PFN1" "EIF4A2-JUN"
#> [41] "EIF4A2-MYL12B" "EIF4A2-NACA" "FOS-HLA-B" "FOS-HLA-C"
#> [45] "FOS-UBA52" "FTL-HLA-B" "HLA-A-JUN" "HLA-A-MYL12B"
#> [49] "HLA-A-TMSB4X" "HLA-B-PFN1" "HLA-C-JUNB" "HLA-C-UBA52"
#> [53] "HLA-E-MYL12B" "HLA-E-MYL6" "HLA-E-NACA" "JUNB-PABPC1"
#> [57] "MYL12B-TMSB4X" "PABPC1-UBA52"
#>
#> $fn_edges
#> [1] "ACTG1-ARPC2" "ACTG1-ARPC3" "ACTG1-CFL1" "ACTG1-MYL6"
#> [5] "ARPC2-ARPC3" "COX4I1-COX7C" "EEF1A1-EEF1D" "EIF1-EIF3K"
#> [9] "EIF1-GNB2L1" "EIF4A2-PABPC1" "GNB2L1-UBA52" "HLA-A-HLA-C"
#> [13] "HLA-B-HLA-E" "MYL12B-MYL6"
#>
#> $consensus_graph
#> IGRAPH ae74c40 UN-- 35 75 --
#> + attr: name (v/c)
#> + edges from ae74c40 (vertex names):
#> [1] ACTG1 --BTF3 ACTG1 --COX4I1 ACTG1 --EEF1A1 ACTG1 --HLA-E ACTG1 --PFN1
#> [6] ACTG1 --TMSB4X ARPC2 --CD3D ARPC2 --TMSB4X ARPC2 --UBC ARPC3 --CXCR4
#> [11] ARPC3 --EIF4A2 ARPC3 --HLA-A ARPC3 --JUN BTF3 --CFL1 BTF3 --FOS
#> [16] BTF3 --NACA CD3D --CD3E CD3D --CD74 CD3D --MYL12B CD3D --PABPC1
#> [21] CD3D --TMSB4X CD3D --UBC CD3E --COX7C CD3E --EIF4A2 CD3E --HLA-A
#> [26] CD3E --TMSB4X CD74 --CXCR4 CD74 --EEF1D CD74 --UBC CFL1 --COX4I1
#> [31] CFL1 --EEF1D COX4I1--HLA-E COX7C --EEF1A1 COX7C --TMSB4X CXCR4 --FTL
#> [36] EEF1A1--EEF2 EEF1A1--FTH1 EEF1A1--JUN EEF1A1--NACA EEF1A1--PFN1
#> + ... omitted several edges
#>
#> $reference_graph
#> IGRAPH 6d1d2d6 UN-- 35 31 --
#> + attr: name (v/c)
#> + edges from 6d1d2d6 (vertex names):
#> [1] ACTG1 --ARPC2 ACTG1 --ARPC3 ACTG1 --CFL1 ACTG1 --MYL6 ACTG1 --PFN1
#> [6] ACTG1 --TMSB4X ARPC2 --ARPC3 BTF3 --NACA CD3D --CD3E CD3E --HLA-A
#> [11] CD74 --CXCR4 COX4I1--COX7C EEF1A1--EEF1D EEF2 --GNB2L1 EEF2 --UBA52
#> [16] EIF1 --EIF3K EIF1 --GNB2L1 EIF4A2--PABPC1 FOS --JUN FOS --JUNB
#> [21] FTH1 --FTL GNB2L1--UBA52 HLA-A --HLA-B HLA-A --HLA-C HLA-A --HLA-E
#> [26] HLA-B --HLA-C HLA-B --HLA-E HLA-C --HLA-E JUN --JUNB MYL12B--MYL6
#> [31] UBA52 --UBC
#>
#> $edge_colors
#> [1] "blue" "blue" "blue" "blue" "red" "red" "blue" "red" "red" "red"
#> [11] "red" "blue" "blue" "red" "red" "blue" "blue" "blue" "red" "red"
#> [21] "red" "blue" "red" "blue" "red" "red" "blue" "red" "red" "blue"
#> [31] "red"
#>
#> $use_stringdb
#> [1] FALSE
Community detection identifies groups of highly interconnected genes that likely share biological functions or regulatory mechanisms. We’ll detect communities in both our inferred consensus network and the ground-truth network, then compare their similarity using several metrics.
Community Detection Steps: 1. Apply Louvain algorithm to identify network modules 2. Compare community assignments using multiple similarity measures: - Variation of Information (VI): Lower values = more similar - Normalized Mutual Information (NMI): Higher values = more similar - Adjusted Rand Index (ARI): Higher values = more similar
Now we’ll detect communities in our ground-truth network to establish the reference community structure:
if (requireNamespace("igraph", quietly = TRUE)) {
comm_truth <- community_path(adj_truth)
}
#> Detecting communities...
#> Running pathway enrichment...
#> 'select()' returned 1:1 mapping between keys and columns
#> 'select()' returned 1:1 mapping between keys and columns
#> 'select()' returned 1:1 mapping between keys and columns
# Detect communities in consensus network
if (requireNamespace("igraph", quietly = TRUE)) {
comm_cons <- community_path(consensus)
}
#> Detecting communities...
#> Running pathway enrichment...
#> 'select()' returned 1:1 mapping between keys and columns
#> 'select()' returned 1:1 mapping between keys and columns
#> 'select()' returned 1:1 mapping between keys and columns
#> 'select()' returned 1:1 mapping between keys and columns
#> 'select()' returned 1:1 mapping between keys and columns
Finally, we’ll quantify how similar the community structures are between our inferred consensus network and the ground-truth network using the broken-down functions:
# Calculate community metrics (VI, NMI, ARI)
if (requireNamespace("igraph", quietly = TRUE)) {
comm_metrics <- compute_community_metrics(comm_truth, list(comm_cons))
print(comm_metrics)
# Calculate topology metrics (modularity, density, etc.)
# Note: compute_topology_metrics expects community_path outputs
# Returns a list with $topology_measures and $control_topology
topo_metrics <- compute_topology_metrics(comm_truth, list(comm_cons))
print(topo_metrics)
# Visualize community comparison
if (requireNamespace("fmsb", quietly = TRUE)) {
plot_community_comparison(
community_metrics = comm_metrics,
topology_measures = topo_metrics$topology_measures,
control_topology = topo_metrics$control_topology
)
}
}
#> VI NMI ARI
#> Predicted_1 2.134704 0.4507208 0.09604142
#> $topology_measures
#> Modularity Communities Density Transitivity
#> Predicted_1 0.4152 6 0.1260504 0.1782178
#>
#> $control_topology
#> Modularity Communities Density Transitivity
#> 0.81373569 10.00000000 0.05210084 0.46666667
Edge mining provides a detailed analysis of network reconstruction performance by categorizing each predicted edge as True Positive (TP), False Positive (FP), True Negative (TN), or False Negative (FN). This gives us insight into which specific regulatory relationships were successfully recovered. The function now returns a list of dataframes matching the documentation.
Edge Mining Analysis: - True Positives (TP): Correctly predicted edges that exist in ground truth - False Positives (FP): Incorrectly predicted edges not in ground truth - False Negatives (FN): Missed edges that exist in ground truth - True Negatives (TN): Correctly identified non-edges
em <- edge_mining(
consensus,
ground_truth = adj_truth,
query_edge_types = "TP"
)
head(em[[1]])
#> gene1 gene2 edge_type pubmed_hits
#> 6 CD3D CD3E TP 60
#> 15 CD74 CXCR4 TP 136
#> 29 FTH1 FTL TP 100
#> 30 EEF2 GNB2L1 TP 0
#> 33 CD3E HLA-A TP 0
#> 37 HLA-A HLA-B TP 2922
#> PMIDs
#> 6 40969714,40831791,40750862,40725440,40702236,40577225,40394466,39958562,39780208,39493321,39220810,38816839,37942472,37841886,37627776,37571885,37478215,37274867,37113759,37031292,36806614,36291815,36276947,36146529,35720370,35303369,35198572,35090306,34966587,34692491,34108992,33901225,33748188,33688252,33628209,33604380,32509861,32424701,31921117,31842801,30885360,29977931,29920501,29653965,28368009,26459776,25946140,25557485,24748432,21883749,21856934,20478055,16264327,15546002,11261926,7757067,1386345,1981047,2331673,3248386
#> 15 40952055,40951943,40948762,40932351,40913267,40831537,40796616,40766141,40740771,40736341,40682854,40575951,40421217,40356050,40342960,40255236,40186663,40066443,40055309,40018995,39915484,39914814,39834330,39712016,39691708,39633216,39502000,39351536,39345457,39323959,39323779,39111632,39084404,39079349,39079345,39030653,39010084,38992165,38874510,38786064,38712609,38674069,38663915,38650025,38455420,38386050,38319415,38237304,38169915,38106024,38061122,37961378,37946732,37725104,37671160,37648811,37642473,37616250,37529341,37508563,37335089,36895934,36816799,36738160,36709495,36096455,35885004,35784460,35563296,35509079,35276027,35252348,35011861,34424927,34098338,34041839,33447370,33239628,33173092,33008022,32104230,31811089,31745887,31745866,31581595,31414920,31105032,31009094,30737274,30716779,30682543,30680604,30567353,30506423,30439447,30371153,30160778,30127875,29935880,29476963
#> 29 41068355,41044642,40948800,40679566,40409664,40358793,40257585,40139473,40048981,39988734,39954839,39906141,39849491,39800242,39534872,39475272,39419454,39334701,39294687,39294443,39148171,39147356,39053339,39041134,38923019,38765702,38740757,38685226,38683121,38561847,38556118,38505909,38393256,38262246,38221933,38160863,38056574,38007415,37976599,37964337,37867341,37851191,37839303,37788777,37664922,37660254,37605582,37555866,37283515,37270049,37198940,37143164,37101205,37094356,37086630,36922863,36909561,36831385,36778397,36728677,36647288,36608444,36562685,36547083,36477858,36066504,35918659,35661779,35598199,35543324,35449524,35249107,35100981,35008695,34938499,34864093,34825691,34787054,34732689,34707088,34660571,34547407,34258619,33864445,33402128,32377595,31920471,31401526,31320750,30325535,30076742,26765579,24496804,24448401,24007662,23766848,21555518,21029774,16760464,15099026
#> 30 <NA>
#> 33 <NA>
#> 37 41069639,41050417,41046214,41035653,41034087,41011272,41003814,41001029,40962674,40944463,40941618,40910546,40896631,40868188,40873549,40866583,40835749,40820324,40791414,40770645,40739281,40721532,40668377,40637775,40617889,40577315,40575193,40574847,40565849,40535491,40517559,40510810,40508101,40508057,40504263,28520367,40475155,40469303,40451247,40405507,40390041,40387827,40366340,40360437,40342296,40314039,40303409,40284905,40239304,40230599,40200142,40171627,40170583,40149438,40127639,40087529,40074802,40071333,40052889,40020215,40010968,39979805,39937552,39933824,39928675,39914254,39901240,39898909,39892763,39890193,39836838,39833815,39824805,39804685,39768617,39738923,39728756,39727970,39726387,39723380,39716692,39715961,39713404,39694865,39684614,39672954,39642776,39605664,39586371,39558667,39530176,39469985,39459542,39456988,39447013,39435845,39374941,39337366,39336927,39306605
#> query_status
#> 6 hits_found
#> 15 hits_found
#> 29 hits_found
#> 30 no_hits
#> 33 no_hits
#> 37 hits_found
sessionInfo()
#> R version 4.5.1 Patched (2025-08-23 r88802)
#> Platform: x86_64-pc-linux-gnu
#> Running under: Ubuntu 24.04.3 LTS
#>
#> Matrix products: default
#> BLAS: /home/biocbuild/bbs-3.22-bioc/R/lib/libRblas.so
#> LAPACK: /usr/lib/x86_64-linux-gnu/lapack/liblapack.so.3.12.0 LAPACK version 3.12.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] TENxPBMCData_1.27.0 HDF5Array_1.37.0
#> [3] h5mread_1.1.1 rhdf5_2.53.6
#> [5] DelayedArray_0.35.3 SparseArray_1.9.1
#> [7] S4Arrays_1.9.1 abind_1.4-8
#> [9] Matrix_1.7-4 SingleCellExperiment_1.31.1
#> [11] SummarizedExperiment_1.39.2 Biobase_2.69.1
#> [13] GenomicRanges_1.61.5 Seqinfo_0.99.2
#> [15] IRanges_2.43.5 S4Vectors_0.47.4
#> [17] BiocGenerics_0.55.4 generics_0.1.4
#> [19] MatrixGenerics_1.21.0 matrixStats_1.5.0
#> [21] scGraphVerse_0.99.20 BiocStyle_2.37.1
#>
#> loaded via a namespace (and not attached):
#> [1] R.methodsS3_1.8.2 dichromat_2.0-0.1
#> [3] gld_2.6.8 Biostrings_2.77.2
#> [5] vctrs_0.6.5 ggtangle_0.0.7
#> [7] perturbR_0.1.3 digest_0.6.37
#> [9] png_0.1-8 shape_1.4.6.1
#> [11] proxy_0.4-27 Exact_3.3
#> [13] pcaPP_2.0-5 BiocBaseUtils_1.11.2
#> [15] gypsum_1.5.0 ggrepel_0.9.6
#> [17] bst_0.3-24 magick_2.9.0
#> [19] hdrcde_3.4 MASS_7.3-65
#> [21] fontLiberation_0.1.0 reshape2_1.4.4
#> [23] foreach_1.5.2 qvalue_2.41.0
#> [25] withr_3.0.2 xfun_0.53
#> [27] ggfun_0.2.0 survival_3.8-3
#> [29] doRNG_1.8.6.2 memoise_2.0.1
#> [31] ggbeeswarm_0.7.2 clusterProfiler_4.17.0
#> [33] gson_0.1.0 systemfonts_1.3.1
#> [35] tidytree_0.4.6 networkD3_0.4.1
#> [37] gtools_3.9.5 R.oo_1.27.1
#> [39] WeightSVM_1.7-16 KEGGREST_1.49.2
#> [41] httr_1.4.7 GENIE3_1.31.0
#> [43] fmsb_0.7.6 hash_2.2.6.3
#> [45] rhdf5filters_1.21.4 rstudioapi_0.17.1
#> [47] DOSE_4.3.0 curl_7.0.0
#> [49] ScaledMatrix_1.17.0 ggraph_2.2.2
#> [51] polyclip_1.10-7 ExperimentHub_2.99.5
#> [53] stringr_1.5.2 pracma_2.4.4
#> [55] doParallel_1.0.17 evaluate_1.0.5
#> [57] BiocFileCache_2.99.6 hms_1.1.4
#> [59] glmnet_4.1-10 bookdown_0.45
#> [61] irlba_2.3.5.1 colorspace_2.1-2
#> [63] filelock_1.0.3 reticulate_1.43.0
#> [65] readxl_1.4.5 magrittr_2.0.4
#> [67] readr_2.1.5 viridis_0.6.5
#> [69] ggtree_3.99.1 lattice_0.22-7
#> [71] XML_3.99-0.19 scuttle_1.19.0
#> [73] cowplot_1.2.0 class_7.3-23
#> [75] pillar_1.11.1 nlme_3.1-168
#> [77] iterators_1.0.14 caTools_1.18.3
#> [79] compiler_4.5.1 beachmat_2.25.5
#> [81] stringi_1.8.7 DescTools_0.99.60
#> [83] plyr_1.8.9 mpath_0.4-2.26
#> [85] fda_6.3.0 crayon_1.5.3
#> [87] scater_1.37.0 gbm_2.2.2
#> [89] gridGraphics_0.5-1 chron_2.3-62
#> [91] haven_2.5.5 graphlayouts_1.2.2
#> [93] org.Hs.eg.db_3.22.0 bit_4.6.0
#> [95] rootSolve_1.8.2.4 dplyr_1.1.4
#> [97] fastmatch_1.1-6 codetools_0.2-20
#> [99] BiocSingular_1.25.0 bslib_0.9.0
#> [101] e1071_1.7-16 lmom_3.2
#> [103] alabaster.ranges_1.9.1 fds_1.8
#> [105] MultiAssayExperiment_1.35.9 splines_4.5.1
#> [107] Rcpp_1.1.0 dbplyr_2.5.1
#> [109] sparseMatrixStats_1.21.0 cellranger_1.1.0
#> [111] knitr_1.50 blob_1.2.4
#> [113] BiocVersion_3.22.0 robin_2.0.0
#> [115] fs_1.6.6 DelayedMatrixStats_1.31.0
#> [117] pscl_1.5.9 expm_1.0-0
#> [119] ggplotify_0.1.3 sqldf_0.4-11
#> [121] tibble_3.3.0 tzdb_0.5.0
#> [123] tweenr_2.0.3 pkgconfig_2.0.3
#> [125] tools_4.5.1 cachem_1.1.0
#> [127] RSQLite_2.4.3 viridisLite_0.4.2
#> [129] DBI_1.2.3 numDeriv_2016.8-1.1
#> [131] distributions3_0.2.3 celldex_1.19.0
#> [133] fastmap_1.2.0 rmarkdown_2.30
#> [135] scales_1.4.0 grid_4.5.1
#> [137] AnnotationHub_3.99.6 sass_0.4.10
#> [139] patchwork_1.3.2 BiocManager_1.30.26
#> [141] graph_1.87.0 alabaster.schemas_1.9.0
#> [143] SingleR_2.11.4 rpart_4.1.24
#> [145] farver_2.1.2 tidygraph_1.3.1
#> [147] gsubfn_0.7 yaml_2.3.10
#> [149] deSolve_1.40 cli_3.6.5
#> [151] purrr_1.1.0 lifecycle_1.0.4
#> [153] askpass_1.2.1 rainbow_3.8
#> [155] mvtnorm_1.3-3 BiocParallel_1.43.4
#> [157] gtable_0.3.6 pROC_1.19.0.1
#> [159] parallel_4.5.1 ape_5.8-1
#> [161] jsonlite_2.0.0 bitops_1.0-9
#> [163] ggplot2_4.0.0 bit64_4.6.0-1
#> [165] yulab.utils_0.2.1 alabaster.matrix_1.9.0
#> [167] BiocNeighbors_2.3.1 proto_1.0.0
#> [169] jquerylib_0.1.4 alabaster.se_1.9.0
#> [171] GOSemSim_2.35.2 R.utils_2.13.0
#> [173] lazyeval_0.2.2 alabaster.base_1.9.5
#> [175] htmltools_0.5.8.1 enrichplot_1.29.2
#> [177] GO.db_3.22.0 rappdirs_0.3.3
#> [179] data.tree_1.2.0 tinytex_0.57
#> [181] glue_1.8.0 STRINGdb_2.21.0
#> [183] httr2_1.2.1 XVector_0.49.1
#> [185] gdtools_0.4.4 RCurl_1.98-1.17
#> [187] qpdf_1.4.1 treeio_1.33.0
#> [189] mclust_6.1.1 ks_1.15.1
#> [191] gridExtra_2.3 boot_1.3-32
#> [193] igraph_2.2.0 R6_2.6.1
#> [195] tidyr_1.3.1 fdatest_2.1.1
#> [197] ggiraph_0.9.2 gplots_3.2.0
#> [199] forcats_1.0.1 labeling_0.4.3
#> [201] cluster_2.1.8.1 rngtools_1.5.2
#> [203] Rhdf5lib_1.31.1 aplot_0.2.9
#> [205] plotrix_3.8-4 tidyselect_1.2.1
#> [207] vipor_0.4.7 ggforce_0.5.0
#> [209] fontBitstreamVera_0.1.1 AnnotationDbi_1.71.2
#> [211] rsvd_1.0.5 KernSmooth_2.23-26
#> [213] S7_0.2.0 fontquiver_0.2.1
#> [215] data.table_1.17.8 htmlwidgets_1.6.4
#> [217] fgsea_1.35.8 RColorBrewer_1.1-3
#> [219] rlang_1.1.6 rentrez_1.2.4
#> [221] beeswarm_0.4.0