Introduction

Why use SETA?

SETA (Single Cell Ecological Taxonomic Analysis) performs compositional analysis of single-cell data by combining established compositional data analysis (CoDA) methods with ecological principles. While other excellent Bioconductor packages address cell-type compositions, SETA is uniquely designed to teach and apply a variety of CoDA methods.

Compositional Analysis of Single-Cell RNA-seq Data

Single-cell RNA-seq experiments produce data where cell-type abundances are inherently compositional - the proportions of different cell types sum to 1 (or 100%). This compositional constraint means that changes in one cell type affect all others, and standard statistical methods can produce misleading results. SETA applies well-established compositional data analysis (CoDA) methods from ecology to properly handle these constraints.

What SETA Does Well

Several Bioconductor packages already address single-cell compositional analysis, each with strong statistical background, including DCATS, speckle, muscat and sccomp. While these packages, especially sccomp, provides the statistical rigor for complex analyses, SETA aims to open up the entire compositional analysis workflow to users who want to understand and apply these methods step-by-step. SETA focuses on simpler models that respect the compositional constraint.

SETA’s unique strengths: - Explicit distance calculations: Aitchison distances and compositional geometry - Ecological metrics: Balance transforms, diversity indices, and ecological framing - Comprehensive transforms: CLR, ALR, ILR, and user-defined balance transforms - Educational focus: separating each step of the workflow into separate functions, SETA aims to make the workflow transparent and easier to explain.

Why These Steps Should Be Executed

The workflow presented in this vignette is essential because:

  1. Proper Data Handling: Cell-type counts must be treated as compositional data to avoid spurious correlations and misleading statistical results
  2. Appropriate Transforms: Log-ratio transforms (CLR, ALR, ILR) maintain compositional properties while enabling standard statistical analysis
  3. Ecological Insights: Applying ecological methods provides novel perspectives on cell-type community structure and diversity
  4. Educational Value: SETA opens up the entire compositional analysis workflow to users, teaching the mathematical foundations rather than hiding them in black-box functions
  5. Unique Capabilities: SETA provides explicit distance calculations, ecological metrics, and comprehensive compositional transforms in a single, educational framework

Package Overview

SETA provides a set of functions for compositional analysis of single-cell RNA-seq data.

In this vignette we show how to:

  • Extract a taxonomic counts matrix from various objects (e.g. a Seurat object)
  • Apply compositional transforms such as CLR and ALR
  • Compute a latent space using PCA (with options for PCoA or NMDS)

This example works with the Tabula Muris Senis lung dataset which can be downloaded directly with bioconductor It contains lung single-cell profiles from 14 donor mice (mouse_id). Each mouse is assigned to one of five age groups, from young adult to old, recorded in the age column. Standard cell-type labels (free_annotation) and additional metadata (e.g. sex) are also included.

This vignette is modular; users can update or replace the dataset-specific settings (e.g., metadata column names, reference cell type) as appropriate for their own data.

Installation

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

Loading the Data

library(SingleCellExperiment)
library(SETA)
library(ggplot2)
library(dplyr)
library(TabulaMurisSenisData)
library(reshape2)
# library(tidytext)

Load and prepare data

sce <- TabulaMurisSenisDroplet(tissues = "Lung")$Lung

# Two samples in this dataset have had certain celltype lineages depleted. 
# For compositional comparisons, we remove these for simplicity.
sce <- sce[, colData(sce)$subtissue != "immune-endo-depleted"]

Quick data exploration:

table(sce$free_annotation, sce$age)
##                                    
##                                       1m   3m  18m  21m  30m
##   Adventitial Fibroblast              29   17   56   43   23
##   Airway Smooth Muscle                 2    0   10    0    0
##   Alveolar Epithelial Type 2           7    2   35   18   11
##   Alveolar Fibroblast                 87   37  166  101   26
##   Alveolar Macrophage                517  116  193   91  276
##   Artery                              13   16   32   26   11
##   B                                   98  135  459  198  180
##   Basophil                             9    8   19   46   40
##   CD4+ T                             143  122   72  107  106
##   CD8+ T                              38   70  234  335  193
##   Capillary                          194  119  617  522  343
##   Capillary Aerocyte                 125   35  129  148   93
##   Ccr7+ Dendritic                      2    2    2    5    7
##   Ciliated                            10    0    9    4    3
##   Classical Monocyte                 485  111  968  242 3615
##   Club                                 2    0    4    5    3
##   Intermediate Monocyte              111   29   95   19 1495
##   Interstitial Macrophage              7   13   16   21 1098
##   Ly6g5b+ T                            0    0   53   65   42
##   Lympatic                             6    2   11   12   10
##   Myeloid Dendritic Type 1            38   20   31   26   44
##   Myeloid Dendritic Type 2             8   10   17   16   26
##   Myofibroblast                       23   13    9    9    0
##   Natural Killer                     140  201  153  224  120
##   Natural Killer T                    20    7  152   50  191
##   Neuroendocrine                       0    0    0    0    0
##   Neutrophil                          64   22  413   29   22
##   Nonclassical Monocyte              157  160  187  194  287
##   Pericyte                            11    7    8    1    8
##   Plasma                               3    1    5   16   23
##   Plasmacytoid Dendritic              21    8   15   10   14
##   Proliferating Alveolar Macrophage   60    9    7    0   14
##   Proliferating Classical Monocyte     1    0   13    0 2473
##   Proliferating Dendritic              4    3    8    4    8
##   Proliferating NK                     5    2    5    3    9
##   Proliferating T                     15    9   11   18   36
##   Regulatory T                         8   11   28   20   55
##   Vein                                19   28   93  107   73
##   Zbtb32+ B                           26   18   41  216  115
table(sce$free_annotation, sce$sex)
##                                    
##                                     female male
##   Adventitial Fibroblast                85   83
##   Airway Smooth Muscle                   9    3
##   Alveolar Epithelial Type 2            24   49
##   Alveolar Fibroblast                  217  200
##   Alveolar Macrophage                  310  883
##   Artery                                48   50
##   B                                    583  487
##   Basophil                              57   65
##   CD4+ T                               267  283
##   CD8+ T                               474  396
##   Capillary                            851  944
##   Capillary Aerocyte                   204  326
##   Ccr7+ Dendritic                        7   11
##   Ciliated                               7   19
##   Classical Monocyte                   514 4907
##   Club                                   6    8
##   Intermediate Monocyte                 90 1659
##   Interstitial Macrophage               34 1121
##   Ly6g5b+ T                             95   65
##   Lympatic                              16   25
##   Myeloid Dendritic Type 1              61   98
##   Myeloid Dendritic Type 2              33   44
##   Myofibroblast                         25   29
##   Natural Killer                       506  332
##   Natural Killer T                     165  255
##   Neuroendocrine                         0    0
##   Neutrophil                            61  489
##   Nonclassical Monocyte                457  528
##   Pericyte                              12   23
##   Plasma                                20   28
##   Plasmacytoid Dendritic                23   45
##   Proliferating Alveolar Macrophage     12   78
##   Proliferating Classical Monocyte       1 2486
##   Proliferating Dendritic               10   17
##   Proliferating NK                       6   18
##   Proliferating T                       33   56
##   Regulatory T                          57   65
##   Vein                                 168  152
##   Zbtb32+ B                            256  160
# Set up a color palette for plots
# 3-group distinct categoricals
age_palette <- c(
    "1m" = "#90EE90",
    "3m" = "#4CBB17",
    "18m" = "#228B22",
    "21m" = "#355E3B",
    "30m" = "#023020")

Extracting the Taxonomic Counts Matrix

The first step in compositional analysis is to extract a cell-type counts matrix from single-cell objects. The setaCounts() function converts cell-level metadata into a sample-by-cell-type count matrix, where each row represents a sample and each column represents a cell type. This transformation is crucial because compositional analysis requires sample-level data, not individual cell data.

For this dataset, we want the cell-type annotations to come from the free_annotation column and the sample identifiers from mouse.id. The function aggregates individual cells by sample and cell type, creating a counts matrix suitable for compositional analysis. If necessary, we reassign these columns before extraction to match the expected column names.

# Extract the taxonomic counts matrix using custom metadata column names
# (Users can specify these column names according to their dataset.)
df <- data.frame(colData(sce))

# Fix special characters in mouse.id before working with the data
df$mouse.id <- gsub("/", "", df$mouse.id)

taxa_counts <- setaCounts(
    df,
    cell_type_col = "free_annotation",
    sample_col = "mouse.id",
    bc_col = "cell")
head(taxa_counts)
##          
##           Adventitial Fibroblast Airway Smooth Muscle
##   1-M-62                      20                    1
##   1-M-63                       9                    1
##   18-F-50                      2                    1
##   18-F-51                     23                    8
##   18-M-52                     16                    1
##   18-M-53                     15                    0
##          
##           Alveolar Epithelial Type 2 Alveolar Fibroblast Alveolar Macrophage
##   1-M-62                           7                  35                 148
##   1-M-63                           0                  52                 369
##   18-F-50                          0                   4                  15
##   18-F-51                          4                  75                  88
##   18-M-52                         20                  36                  49
##   18-M-53                         11                  51                  41
##          
##           Artery   B Basophil CD4+ T CD8+ T Capillary Capillary Aerocyte
##   1-M-62       5  72        6     48     16       119                 90
##   1-M-63       8  26        3     95     22        75                 35
##   18-F-50      1  83        2      5     12        38                  7
##   18-F-51      5 167        1     33     57       172                 14
##   18-M-52     12 100       10     11     79       229                 59
##   18-M-53     14 109        6     23     86       178                 49
##          
##           Ccr7+ Dendritic Ciliated Classical Monocyte Club
##   1-M-62                1        5                201    1
##   1-M-63                1        5                284    1
##   18-F-50               0        0                 39    1
##   18-F-51               0        3                122    0
##   18-M-52               1        2                462    1
##   18-M-53               1        4                345    2
##          
##           Intermediate Monocyte Interstitial Macrophage Ly6g5b+ T Lympatic
##   1-M-62                     50                       3         0        2
##   1-M-63                     61                       4         0        4
##   18-F-50                     5                       0         7        0
##   18-F-51                    37                       0        23        2
##   18-M-52                    23                       6        15        6
##   18-M-53                    30                      10         8        3
##          
##           Myeloid Dendritic Type 1 Myeloid Dendritic Type 2 Myofibroblast
##   1-M-62                        17                        4            10
##   1-M-63                        21                        4            13
##   18-F-50                        4                        1             1
##   18-F-51                       11                        6             2
##   18-M-52                        3                        7             2
##   18-M-53                       13                        3             4
##          
##           Natural Killer Natural Killer T Neutrophil Nonclassical Monocyte
##   1-M-62              53               10         13                    80
##   1-M-63              87               10         51                    77
##   18-F-50             27               11          3                    41
##   18-F-51             54               97          7                    62
##   18-M-52             34                3        307                    45
##   18-M-53             38               41         96                    39
##          
##           Pericyte Plasma Plasmacytoid Dendritic
##   1-M-62         8      1                      6
##   1-M-63         3      2                     15
##   18-F-50        2      0                      0
##   18-F-51        2      3                      5
##   18-M-52        3      0                      3
##   18-M-53        1      2                      7
##          
##           Proliferating Alveolar Macrophage Proliferating Classical Monocyte
##   1-M-62                                 10                                1
##   1-M-63                                 50                                0
##   18-F-50                                 0                                0
##   18-F-51                                 3                                1
##   18-M-52                                 1                               12
##   18-M-53                                 3                                0
##          
##           Proliferating Dendritic Proliferating NK Proliferating T Regulatory T
##   1-M-62                        2                1               4            1
##   1-M-63                        2                4              11            7
##   18-F-50                       1                0               2            4
##   18-F-51                       2                1               4           22
##   18-M-52                       1                1               1            1
##   18-M-53                       4                3               4            1
##          
##           Vein Zbtb32+ B
##   1-M-62    11        16
##   1-M-63     8        10
##   18-F-50    3         9
##   18-F-51   30        13
##   18-M-52   30        11
##   18-M-53   30         8

Applying Compositional Transforms

Once we have our cell-type counts matrix, we must apply appropriate compositional transforms before statistical analysis. Raw count data cannot be used directly because it violates the compositional constraint (proportions must sum to 1). Log-ratio transforms solve this problem by converting the data to a space where standard statistical methods can be applied.

Below we demonstrate four different compositional transforms, each with specific advantages:

CLR Transformation

# CLR Transformation with a pseudocount of 1 (default)
clr_out <- setaCLR(taxa_counts, pseudocount = 1)
print(clr_out$counts[1:5, 1:5])
##          
##           Adventitial Fibroblast Airway Smooth Muscle
##   1-M-62               0.6493376           -1.7020377
##   1-M-63              -0.2008121           -1.8102500
##   18-F-50             -0.2928298           -0.6982949
##   18-F-51              0.7495280           -0.2313013
##   18-M-52              0.4036570           -1.7364092
##          
##           Alveolar Epithelial Type 2 Alveolar Fibroblast Alveolar Macrophage
##   1-M-62                  -0.3157433           1.1883341           2.6087615
##   1-M-63                  -2.5033972           1.4668947           3.4101058
##   18-F-50                 -1.3914421           0.2179958           1.3811467
##   18-F-51                 -0.8190880           1.9022075           2.0601105
##   18-M-52                  0.6149661           1.1813615           1.4824666

For the ALR transform we need to choose a reference cell type whose numbers remain stable across conditions pre-sampling. In an ideal scenario this could a spike-in of a known number of cells. In this example, we assume that one of the reference cell types are “NK cells”. If it does not exist in your dataset, change the reference accordingly or use CLR.

# ALR Transformation: using "NK cell" as the reference cell type
if ("Natural Killer" %in% colnames(taxa_counts)) {
    alr_out <- setaALR(taxa_counts, ref = "Natural Killer", pseudocount = 1)
    print(alr_out$counts[1:5, 1:5])
} else {
    message("Reference 'NK cell' not found in taxa_counts. 
            Please choose an available cell type.")
}
##          
##           Adventitial Fibroblast Airway Smooth Muscle
##   1-M-62             -0.94446161          -3.29583687
##   1-M-63             -2.17475172          -3.78418963
##   18-F-50            -2.23359222          -2.63905733
##   18-F-51            -0.82927935          -1.81010861
##   18-M-52            -0.72213472          -2.86220088
##          
##           Alveolar Epithelial Type 2 Alveolar Fibroblast Alveolar Macrophage
##   1-M-62                 -1.90954250         -0.40546511          1.01496226
##   1-M-63                 -4.47733681         -0.50704490          1.43616619
##   18-F-50                -3.33220451         -1.72276660         -0.55961579
##   18-F-51                -2.39789527          0.32340016          0.48130318
##   18-M-52                -0.51082562          0.05556985          0.35667494

ILR Transform using Helmert basis (boxcox_p = 0 for standard ILR) This is the method used by Cacoa (Viktor Petukhov & Kharchenko Lab)

ilr_out <- setaILR(taxa_counts, boxcox_p = 0, pseudocount = 1)
print(ilr_out$counts[1:5, 1:5])
##          
##                 [,1]       [,2]      [,3]     [,4]       [,5]
##   1-M-62  -1.6626734  0.1719597 1.4241631 2.373621 -0.9942587
##   1-M-63  -1.1380445 -1.2230026 2.5735802 3.731548 -0.3456867
##   18-F-50 -0.2867071 -0.7314827 0.8765777 1.719348 -0.4944201
##   18-F-51 -0.6935510 -0.8803477 1.7342112 1.484547 -1.2497730
##   18-M-52 -1.5132553  1.0462115 1.2302961 1.222300 -0.2317007

Simple Percentage Transform

pct_out <- setaPercent(taxa_counts)
print(pct_out$counts[1:5, 1:5])
##          
##           Adventitial Fibroblast Airway Smooth Muscle
##   1-M-62              1.85528757           0.09276438
##   1-M-63              0.62937063           0.06993007
##   18-F-50             0.60422961           0.30211480
##   18-F-51             1.98446937           0.69025022
##   18-M-52             0.99812851           0.06238303
##          
##           Alveolar Epithelial Type 2 Alveolar Fibroblast Alveolar Macrophage
##   1-M-62                  0.64935065          3.24675325         13.72912801
##   1-M-63                  0.00000000          3.63636364         25.80419580
##   18-F-50                 0.00000000          1.20845921          4.53172205
##   18-F-51                 0.34512511          6.47109577          7.59275237
##   18-M-52                 1.24766064          2.24578915          3.05676856

LogCPM (counts per 10k) Transform

lcpm_out <- setaLogCPM(taxa_counts, pseudocount = 1)
print(lcpm_out$counts[1:5, 1:5])
##          
##           Adventitial Fibroblast Airway Smooth Muscle
##   1-M-62               14.248407            10.856089
##   1-M-63               12.770689            10.448761
##   18-F-50              13.141492            12.556529
##   18-F-51              14.336622            12.921584
##   18-M-52              13.371573            10.284110
##          
##           Alveolar Epithelial Type 2 Alveolar Fibroblast Alveolar Macrophage
##   1-M-62                   12.856089           15.026014           17.075258
##   1-M-63                    9.448761           15.176681           17.980142
##   18-F-50                  11.556529           13.878457           15.556529
##   18-F-51                  12.073588           15.999587           16.227393
##   18-M-52                  13.676428           14.493564           14.927966

Latent Space Analysis

After applying compositional transforms, we can perform dimensionality reduction to visualize and analyze the relationships between samples. The transformed data can be analyzed using standard multivariate methods, but it’s important to remember that we’re now working in log-ratio space, not the original count space.

We demonstrate three different latent space methods, each with specific advantages:

latent_pca <- setaLatent(clr_out, method = "PCA", dims = 5)
# PCA Latent Space Coordinates:
head(latent_pca$latentSpace)
##                PC1        PC2        PC3         PC4         PC5
## 1-M-62  -0.3350223  3.1938325 -0.4848707  0.04635990 -0.76225078
## 1-M-63  -0.5097368  4.8874855  0.9040002  0.04858906  0.36921601
## 18-F-50 -1.5751990 -0.8085481  1.1262639 -1.31190158 -0.78190046
## 18-F-51 -2.3212913 -0.6116102  1.3477482 -2.63828824 -2.03933177
## 18-M-52  0.3601475 -0.5633859 -4.5736259 -1.47565863 -0.01092873
## 18-M-53 -1.1101814  0.2754001 -2.4996341 -0.80132618  0.66184605
# Variance Explained:
latent_pca$varExplained
## [1] 0.49890109 0.15856671 0.10540212 0.06273258 0.05150588

Visualization

Below we show how to create basic visualizations of SETA latent spaces

Variance Explained Plot

Plot the variance explained by the first few principal components.

ve_df <- data.frame(
    PC = factor(seq_along(latent_pca$varExplained), 
                levels = seq_along(latent_pca$varExplained)),
    VarianceExplained = cumsum(latent_pca$varExplained)
)
ggplot(ve_df, aes(x = PC, y = VarianceExplained)) +
    geom_line(color = "#2ca1db", size = 1.5, group = 1) +
    geom_point(color = "#f44323", size = 3) +
    labs(title = "Cumulative Variance Explained by Principal Components",
        x = "Principal Component", y = "Cumulative Variance Explained (%)") +
    theme_minimal() +
    ylim(c(0, 1)) +
    theme(text = element_text(size = 12))
## Warning: Using `size` aesthetic for lines was deprecated in ggplot2 3.4.0.
## ℹ Please use `linewidth` instead.
## This warning is displayed once every 8 hours.
## Call `lifecycle::last_lifecycle_warnings()` to see where this warning was
## generated.

PCA Scatter plot

We overlay the PCA latent space coordinates (using the CLR transform) with metadata from the Seurat object. In this example, points are colored by age.

# Extract latent space coordinates from PCA result.
pca_coords <- latent_pca$latentSpace

# Extract metadata
meta_df <- setaMetadata(
    df,
    sample_col = "mouse.id",
    meta_cols = c("age", "sex"))

# Merge latent coordinates with metadata
pca_plot_df <- cbind(pca_coords, meta_df)
# Create the scatter plot
ggplot(pca_plot_df, aes(x = PC1, y = PC2, color = age)) +
    geom_text(aes(label = sample_id), size = 5) +
    scale_color_manual(
        values = age_palette
    ) +
    labs(title = "PCA Scatter Plot",
        x = "PC1", y = "PC2", color = "Age") +
    theme_linedraw() +
    xlab(sprintf("PC1 (%s%%)", signif(latent_pca$varExplained[1], 4) * 100)) +
    ylab(sprintf("PC2 (%s%%)", signif(latent_pca$varExplained[2], 4) * 100))

Loadings Plot

Plot the variance explained by interesting principal components.

loadings_df <- as.data.frame(latent_pca$loadings)
loadings_df$Celltype <- rownames(loadings_df)

# Melt the loadings into long format for `ggplot2`
loadings_long <- melt(loadings_df, id.vars = "Celltype",
                    variable.name = "PC", value.name = "Loading")

Visualize the PC loadings for each cell

# Subset rows where PC is "PC1" or "PC2"
loading_df <- loadings_long[loadings_long$PC %in% c("PC1", "PC2"), ]

# Append PC to Celltype to mimic tidytext::reorder_within()
loading_df$Celltype_PC <- paste(loading_df$Celltype, loading_df$PC, sep = "___")

# Reorder by Loading *within each PC*
loading_df$Celltype_PC <- with(loading_df, reorder(Celltype_PC, Loading))

ggplot(loading_df,
        aes(x = Loading, y = Celltype, label = Celltype, fill = Loading)) +
    geom_col() +
    facet_wrap(~PC, scales = "free") +
    # scale_y_reordered() +
    labs(title = "PC Loadings",
        x = "Loading",
        y = "Cell Type") +
    theme_linedraw() +
    scale_fill_viridis_c() +
    expand_limits(x = c(-0.2, 0.4))

Conclusion

This vignette illustrated a basic workflow using the SETA package:

  1. Extracting a taxonomic counts matrix from a Seurat Object with setaCounts()

  2. Applying compositional transforms (CLR, ALR, ILR, and simple transforms)

  3. Performing latent space analysis (PCA)

  4. Visualizing sample-level latent spaces

Each step is modular, so you can easily swap in your own data, update metadata column names, or choose different parameters (e.g., pseudocount, reference cell type) to suit your analysis needs.

For further details, see the package documentation and relevant literature on compositional data analysis (e.g., Aitchison, 1982).

Session Info

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] reshape2_1.4.4              rhdf5_2.53.6               
##  [3] TabulaMurisSenisData_1.15.1 caret_7.0-1                
##  [5] lattice_0.22-7              corrplot_0.95              
##  [7] tidyr_1.3.1                 dplyr_1.1.4                
##  [9] ggplot2_4.0.0               SETA_0.99.7                
## [11] SingleCellExperiment_1.31.1 SummarizedExperiment_1.39.2
## [13] Biobase_2.69.1              GenomicRanges_1.61.6       
## [15] Seqinfo_0.99.3              IRanges_2.43.5             
## [17] S4Vectors_0.47.4            BiocGenerics_0.55.4        
## [19] generics_0.1.4              MatrixGenerics_1.21.0      
## [21] matrixStats_1.5.0          
## 
## loaded via a namespace (and not attached):
##   [1] DBI_1.2.3            pROC_1.19.0.1        httr2_1.2.1         
##   [4] rlang_1.1.6          magrittr_2.0.4       e1071_1.7-16        
##   [7] compiler_4.5.1       RSQLite_2.4.3        gdata_3.0.1         
##  [10] png_0.1-8            vctrs_0.6.5          stringr_1.5.2       
##  [13] pkgconfig_2.0.3      shape_1.4.6.1        crayon_1.5.3        
##  [16] fastmap_1.2.0        dbplyr_2.5.1         XVector_0.49.1      
##  [19] labeling_0.4.3       rmarkdown_2.30       prodlim_2025.04.28  
##  [22] bit_4.6.0            purrr_1.1.0          xfun_0.53           
##  [25] glmnet_4.1-10        cachem_1.1.0         jsonlite_2.0.0      
##  [28] blob_1.2.4           recipes_1.3.1        rhdf5filters_1.21.4 
##  [31] DelayedArray_0.35.3  Rhdf5lib_1.31.1      parallel_4.5.1      
##  [34] R6_2.6.1             bslib_0.9.0          stringi_1.8.7       
##  [37] RColorBrewer_1.1-3   parallelly_1.45.1    rpart_4.1.24        
##  [40] lubridate_1.9.4      jquerylib_0.1.4      Rcpp_1.1.0          
##  [43] iterators_1.0.14     knitr_1.50           future.apply_1.20.0 
##  [46] Matrix_1.7-4         splines_4.5.1        nnet_7.3-20         
##  [49] igraph_2.2.0         timechange_0.3.0     tidyselect_1.2.1    
##  [52] dichromat_2.0-0.1    abind_1.4-8          yaml_2.3.10         
##  [55] timeDate_4051.111    codetools_0.2-20     curl_7.0.0          
##  [58] listenv_0.9.1        tibble_3.3.0         plyr_1.8.9          
##  [61] KEGGREST_1.49.2      withr_3.0.2          S7_0.2.0            
##  [64] evaluate_1.0.5       future_1.67.0        survival_3.8-3      
##  [67] proxy_0.4-27         BiocFileCache_2.99.6 ExperimentHub_2.99.6
##  [70] Biostrings_2.77.2    filelock_1.0.3       pillar_1.11.1       
##  [73] BiocManager_1.30.26  foreach_1.5.2        BiocVersion_3.22.0  
##  [76] scales_1.4.0         gtools_3.9.5         BiocStyle_2.37.1    
##  [79] globals_0.18.0       class_7.3-23         glue_1.8.0          
##  [82] tools_4.5.1          AnnotationHub_3.99.6 data.table_1.17.8   
##  [85] ModelMetrics_1.2.2.2 gower_1.0.2          tidygraph_1.3.1     
##  [88] grid_4.5.1           AnnotationDbi_1.71.2 ipred_0.9-15        
##  [91] nlme_3.1-168         HDF5Array_1.37.0     cli_3.6.5           
##  [94] rappdirs_0.3.3       viridisLite_0.4.2    S4Arrays_1.9.1      
##  [97] lava_1.8.1           gtable_0.3.6         sass_0.4.10         
## [100] digest_0.6.37        SparseArray_1.9.1    farver_2.1.2        
## [103] memoise_2.0.1        htmltools_0.5.8.1    lifecycle_1.0.4     
## [106] h5mread_1.1.1        httr_1.4.7           hardhat_1.4.2       
## [109] bit64_4.6.0-1        MASS_7.3-65