Contents

1 Introduction

BERT (Batch-Effect Removal with Trees) offers flexible and efficient batch effect correction of omics data, while providing maximum tolerance to missing values. Tested on multiple datasets from proteomic analyses, BERT offered a typical 5-10x runtime improvement over existing methods, while retaining more numeric values and preserving batch effect reduction quality.

As such, BERT is a valuable preprocessing tool for data analysis workflows, in particular for proteomic data. By providing BERT via Bioconductor, we make this tool available to a wider research community. An accompanying research paper is currently under preparation and will be made public soon.

BERT addresses the same fundamental data integration challenges than the [HarmonizR][https://github.com/HSU-HPC/HarmonizR] package, which is released on Bioconductor in November 2023. However, various algorithmic modications and optimizations of BERT provide better execution time and better data coverage than HarmonizR. Moreover, BERT offers a more user-friendly design and a less error-prone input format.

Please note that our package BERT is neither affiliated with nor related to Bidirectional Encoder Representations from Transformers as published by Google.

Please report any questions and issues in the GitHub forum, the BioConductor forum or directly contact the authors,

2 Installation

Please download and install a current version of R (Windows binaries). You might want to consider installing a development environment as well, e.g. RStudio. Finally, BERT can be installed via Bioconductor using

if (!require("BiocManager", quietly = TRUE)){
    install.packages("BiocManager")
}
BiocManager::install("BERT")

which will install all required dependencies. To install the development version of BERT, you can use devtools as follows

devtools::install_github("HSU-HPC/BERT")

which may require the manual installation of the dependencies sva and limma.

if (!require("BiocManager", quietly = TRUE)){
    install.packages("BiocManager")
}
BiocManager::install("sva")
BiocManager::install("limma")

3 Data Preparation

As input, BERT requires a dataframe1 Matrices and SummarizedExperiments work as well, but will automatically be converted to dataframes. with samples in rows and features in columns. For each sample, the respective batch should be indicated by an integer or string in a corresponding column labelled Batch. Missing values should be labelled as NA. A valid example dataframe could look like this:

example = data.frame(feature_1 = stats::rnorm(5), feature_2 = stats::rnorm(5), Batch=c(1,1,2,2,2))
example
#>    feature_1  feature_2 Batch
#> 1  1.2554067 -2.3213518     1
#> 2  0.1932058 -0.6996628     1
#> 3  0.2829865 -0.2355690     2
#> 4 -0.9763425  0.5470149     2
#> 5 -0.1558776  0.8077365     2

Note that each batch should contain at least two samples. Optional columns that can be passed are

Note that BERT tries to find all metadata information for a SummarizedExperiment, including the mandatory batch information, using colData. For instance, a valid SummarizedExperiment might be defined as

nrows <- 200
ncols <- 8
expr_values <- matrix(runif(nrows * ncols, 1, 1e4), nrows)
# colData also takes all other metadata information, such as Label, Sample,
# Covariables etc.
colData <- data.frame(Batch=c(1,1,1,1,2,2,2,2), Reference=c(1,1,0,0,1,1,0,0))
dataset_raw = SummarizedExperiment::SummarizedExperiment(assays=list(expr=expr_values), colData=colData)

4 Basic Usage

BERT can be invoked by importing the BERT library and calling the BERT function. The batch effect corrected data is returned as a dataframe that mirrors the input dataframe2 In particular, the row and column names are in the same order and the optional columns are preserved..

library(BERT)
# generate test data with 10% missing values as provided by the BERT library
dataset_raw <- generate_dataset(features=60, batches=10, samplesperbatch=10, mvstmt=0.1, classes=2)
# apply BERT
dataset_adjusted <- BERT(dataset_raw)
#> 2024-10-29 16:20:16.13059 INFO::Formatting Data.
#> 2024-10-29 16:20:16.140123 INFO::Replacing NaNs with NAs.
#> 2024-10-29 16:20:16.149629 INFO::Removing potential empty rows and columns
#> 2024-10-29 16:20:16.434418 INFO::Found  600  missing values.
#> 2024-10-29 16:20:16.448754 INFO::Introduced 0 missing values due to singular proteins at batch/covariate level.
#> 2024-10-29 16:20:16.449506 INFO::Done
#> 2024-10-29 16:20:16.450109 INFO::Acquiring quality metrics before batch effect correction.
#> 2024-10-29 16:20:16.464285 INFO::Starting hierarchical adjustment
#> 2024-10-29 16:20:16.465354 INFO::Found  10  batches.
#> 2024-10-29 16:20:16.465971 INFO::Cores argument is not defined or BPPARAM has been specified. Argument corereduction will not be used.
#> 2024-10-29 16:20:19.652829 INFO::Using default BPPARAM
#> 2024-10-29 16:20:19.653636 INFO::Processing subtree level 1
#> 2024-10-29 16:20:21.569597 INFO::Processing subtree level 2
#> 2024-10-29 16:20:23.437314 INFO::Adjusting the last 1 batches sequentially
#> 2024-10-29 16:20:23.439917 INFO::Done
#> 2024-10-29 16:20:23.440765 INFO::Acquiring quality metrics after batch effect correction.
#> 2024-10-29 16:20:23.447742 INFO::ASW Batch was 0.444124968779769 prior to batch effect correction and is now -0.125734461251687 .
#> 2024-10-29 16:20:23.448641 INFO::ASW Label was 0.368137864957289 prior to batch effect correction and is now 0.841303712305322 .
#> 2024-10-29 16:20:23.450168 INFO::Total function execution time is  7.34107065200806  s and adjustment time is  6.9748649597168 s ( 95.01 )

BERT uses the logging library to convey live information to the user during the adjustment procedure. The algorithm first verifies the shape and suitability of the input dataframe (lines 1-6) before continuing with the actual batch effect correction (lines 8-14). BERT measure batch effects before and after the correction step by means of the average silhouette score (ASW) with respect to batch and labels (lines 7 and 15). The ASW Label should increase in a successful batch effect correction, whereas low values (\(\leq 0\)) are desireable for the ASW Batch3 The optimum of ASW Label is 1, which is typically however not achieved on real-world datasets. Also, the optimum of ASW Batch can vary, depending on the class distributions of the batches.. Finally, BERT prints the total function execution time (including the computation time for the quality metrics).

5 Advanced Options

5.1 Parameters

BERT offers a large number of parameters to customize the batch effect adjustment. The full function call, including all defaults is

BERT(data, cores = NULL, combatmode = 1, corereduction=2, stopParBatches=2, backend="default", method="ComBat", qualitycontrol=TRUE, verify=TRUE, labelname="Label", batchname="Batch", referencename="Reference", samplename="Sample", covariatename=NULL, BPPARAM=NULL, assayname=NULL)

In the following, we list the respective meaning of each parameter: - data: The input dataframe/matrix/SummarizedExperiment to adjust. See Data Preparation for detailed formatting instructions. - data The data for batch-effect correction. Must contain at least two samples per batch and 2 features.

  • cores: BERT uses BiocParallel for parallelization. If the user specifies a value cores, BERT internally creates and uses a new instance of BiocParallelParam, which is however not exhibited to the user. Setting this parameter can speed up the batch effect adjustment considerably, in particular for large datasets and on unix-based operating systems. A value between \(2\) and \(4\) is a reasonable choice for typical commodity hardware. Multi-node computations are not supported as of now. If, however, cores is not specified, BERT will default to BiocParallel::bpparam(), which may have been set by the user or the system. Additionally, the user can directly specify a specific instance of BiocParallelParam to be used via the BPPARAM argument.
  • combatmode An integer that encodes the parameters to use for ComBat.
Value par.prior mean.only
1 TRUE FALSE
2 TRUE TRUE
3 FALSE FALSE
4 FALSE TRUE

The value of this parameter will be ignored, if method!="ComBat".

  • corereduction Positive integer indicating the factor by which the number of processes should be reduced, once no further adjustment is possible for the current number of batches.4 E.g. consider a BERT call with 8 batches and 8 processes. Further adjustment is not possible with this number of processes, since batches are always processed in pairs. With corereduction=2, the number of processes for the following adjustment steps would be set to \(8/2=4\), which is the maximum number of usable processes for this example. This parameter is used only, if the user specified a custom value for parameter cores.

  • stopParBatches Positive integer indicating the minimum number of batches required at a hierarchy level to proceed with parallelized adjustment. If the number of batches is smaller, adjustment will be performed sequentially to avoid communication overheads.

  • backend: The backend to use for inter-process communication. Possible choices are default and file, where the former refers to the default communication backend of the requested parallelization mode and the latter will create temporary .rds files for data communication. ‘default’ is usually faster for small to medium sized datasets.

  • method: The method to use for the underlying batch effect correction steps. Should be either ComBat, limma for limma::removeBatchEffects or ref for adjustment using specified references (cf. Data Preparation). The underlying batch effect adjustment method for ref is a modified version of the limma method.

  • qualitycontrol: A boolean to (de)activate the ASW computation. Deactivating the ASW computations accelerates the computations.

  • verify: A boolean to (de)activate the initial format check of the input data. Deactivating this verification step accelerates the computations.

  • labelname: A string containing the name of the column to use as class labels. The default is “Label”.

  • batchname: A string containing the name of the column to use as batch labels. The default is “Batch”.

  • referencename: A string containing the name of the column to use as reference labels. The default is “Reference”.

  • covariatename: A vector containing the names of columns with categorical covariables.The default is NULL, in which case all column names are matched agains the pattern “Cov”.

  • BPPARAM: An instance of BiocParallelParam that will be used for parallelization. The default is null, in which case the value of cores determines the behaviour of BERT.

  • assayname: If the user chooses to pass a SummarizedExperiment object, they need to specify the name of the assay that they want to apply BERT to here. BERT then returns the input SummarizedExperiment with an additional assay labeled assayname_BERTcorrected.

5.2 Verbosity

BERT utilizes the logging package for output. The user can easily specify the verbosity of BERT by setting the global logging level in the script. For instance

logging::setLevel("WARN") # set level to warn and upwards
result <- BERT(data,cores = 1) # BERT executes silently

5.3 Choosing the Optimal Number of Cores

BERT exhibits a large number of parameters for parallelisation as to provide users with maximum flexibility. For typical scenarios, however, the default parameters are well suited. For very large experiments (\(>15\) batches), we recommend to increase the number of cores (a reasonable value is \(4\) but larger values may be possible on your hardware). Most users should leave all parameters to their respective default.

6 Examples

In the following, we present simple cookbook examples for BERT usage. Note that ASWs (and runtime) will most likely differ on your machine, since the data generating process involves multiple random choices.

6.1 Sequential Adjustment with limma

Here, BERT uses limma as underlying batch effect correction algorithm (method='limma') and performs all computations on a single process (cores parameter is left on default).

# import BERT
library(BERT)
# generate data with 30 batches, 60 features, 15 samples per batch, 15% missing values and 2 classes
dataset_raw <- generate_dataset(features=60, batches=20, samplesperbatch=15, mvstmt=0.15, classes=2)
# BERT
dataset_adjusted <- BERT(dataset_raw, method="limma")
#> 2024-10-29 16:20:23.546412 INFO::Formatting Data.
#> 2024-10-29 16:20:23.547878 INFO::Replacing NaNs with NAs.
#> 2024-10-29 16:20:23.549653 INFO::Removing potential empty rows and columns
#> 2024-10-29 16:20:23.553336 INFO::Found  2700  missing values.
#> 2024-10-29 16:20:23.592385 INFO::Introduced 0 missing values due to singular proteins at batch/covariate level.
#> 2024-10-29 16:20:23.5931 INFO::Done
#> 2024-10-29 16:20:23.593636 INFO::Acquiring quality metrics before batch effect correction.
#> 2024-10-29 16:20:23.607054 INFO::Starting hierarchical adjustment
#> 2024-10-29 16:20:23.607836 INFO::Found  20  batches.
#> 2024-10-29 16:20:23.60834 INFO::Cores argument is not defined or BPPARAM has been specified. Argument corereduction will not be used.
#> 2024-10-29 16:20:23.608948 INFO::Using default BPPARAM
#> 2024-10-29 16:20:23.609441 INFO::Processing subtree level 1
#> 2024-10-29 16:20:24.202383 INFO::Processing subtree level 2
#> 2024-10-29 16:20:24.635275 INFO::Processing subtree level 3
#> 2024-10-29 16:20:25.120103 INFO::Adjusting the last 1 batches sequentially
#> 2024-10-29 16:20:25.121893 INFO::Done
#> 2024-10-29 16:20:25.122475 INFO::Acquiring quality metrics after batch effect correction.
#> 2024-10-29 16:20:25.132688 INFO::ASW Batch was 0.428201566933583 prior to batch effect correction and is now -0.115844958910403 .
#> 2024-10-29 16:20:25.133331 INFO::ASW Label was 0.353473648235071 prior to batch effect correction and is now 0.840396330146731 .
#> 2024-10-29 16:20:25.134171 INFO::Total function execution time is  1.58811616897583  s and adjustment time is  1.51422023773193 s ( 95.35 )

6.2 Parallel Batch Effect Correction with ComBat

Here, BERT uses ComBat as underlying batch effect correction algorithm (method is left on default) and performs all computations on a 2 processes (cores=2).

# import BERT
library(BERT)
# generate data with 30 batches, 60 features, 15 samples per batch, 15% missing values and 2 classes
dataset_raw <- generate_dataset(features=60, batches=20, samplesperbatch=15, mvstmt=0.15, classes=2)
# BERT
dataset_adjusted <- BERT(dataset_raw, cores=2)
#> 2024-10-29 16:20:25.172447 INFO::Formatting Data.
#> 2024-10-29 16:20:25.173232 INFO::Replacing NaNs with NAs.
#> 2024-10-29 16:20:25.17433 INFO::Removing potential empty rows and columns
#> 2024-10-29 16:20:25.176589 INFO::Found  2700  missing values.
#> 2024-10-29 16:20:25.203843 INFO::Introduced 0 missing values due to singular proteins at batch/covariate level.
#> 2024-10-29 16:20:25.20466 INFO::Done
#> 2024-10-29 16:20:25.20527 INFO::Acquiring quality metrics before batch effect correction.
#> 2024-10-29 16:20:25.21967 INFO::Starting hierarchical adjustment
#> 2024-10-29 16:20:25.220696 INFO::Found  20  batches.
#> 2024-10-29 16:20:26.125638 INFO::Set up parallel execution backend with 2 workers
#> 2024-10-29 16:20:26.127246 INFO::Processing subtree level 1 with 20 batches using 2 cores.
#> 2024-10-29 16:20:29.169876 INFO::Adjusting the last 2 batches sequentially
#> 2024-10-29 16:20:29.171456 INFO::Adjusting sequential tree level 1 with 2 batches
#> 2024-10-29 16:20:30.804743 INFO::Done
#> 2024-10-29 16:20:30.805404 INFO::Acquiring quality metrics after batch effect correction.
#> 2024-10-29 16:20:30.814218 INFO::ASW Batch was 0.49499769271194 prior to batch effect correction and is now -0.121505295866779 .
#> 2024-10-29 16:20:30.814796 INFO::ASW Label was 0.290454304285838 prior to batch effect correction and is now 0.808326502231165 .
#> 2024-10-29 16:20:30.815488 INFO::Total function execution time is  5.64313817024231  s and adjustment time is  5.5839102268219 s ( 98.95 )

6.3 Batch Effect Correction Using SummarizedExperiment

Here, BERT takes the input data using a SummarizedExperiment instead. Batch effect correction is then performed using ComBat as underlying algorithm (method is left on default) and all computations are performed on a single process (cores parameter is left on default).

nrows <- 200
ncols <- 8
# SummarizedExperiments store samples in columns and features in rows (in contrast to BERT).
# BERT will automatically account for this.
expr_values <- matrix(runif(nrows * ncols, 1, 1e4), nrows)
# colData also takes further metadata information, such as Label, Sample,
# Reference or Covariables
colData <- data.frame("Batch"=c(1,1,1,1,2,2,2,2), "Label"=c(1,2,1,2,1,2,1,2), "Sample"=c(1,2,3,4,5,6,7,8))
dataset_raw = SummarizedExperiment::SummarizedExperiment(assays=list(expr=expr_values), colData=colData)
dataset_adjusted = BERT(dataset_raw, assayname = "expr")
#> 2024-10-29 16:20:30.871269 INFO::Formatting Data.
#> 2024-10-29 16:20:30.871973 INFO::Recognized SummarizedExperiment
#> 2024-10-29 16:20:30.872465 INFO::Typecasting input to dataframe.
#> 2024-10-29 16:20:30.905877 INFO::Replacing NaNs with NAs.
#> 2024-10-29 16:20:30.907017 INFO::Removing potential empty rows and columns
#> 2024-10-29 16:20:30.910151 INFO::Found  0  missing values.
#> 2024-10-29 16:20:30.916079 INFO::Introduced 0 missing values due to singular proteins at batch/covariate level.
#> 2024-10-29 16:20:30.91662 INFO::Done
#> 2024-10-29 16:20:30.917108 INFO::Acquiring quality metrics before batch effect correction.
#> 2024-10-29 16:20:30.920862 INFO::Starting hierarchical adjustment
#> 2024-10-29 16:20:30.921508 INFO::Found  2  batches.
#> 2024-10-29 16:20:30.921997 INFO::Cores argument is not defined or BPPARAM has been specified. Argument corereduction will not be used.
#> 2024-10-29 16:20:30.922523 INFO::Using default BPPARAM
#> 2024-10-29 16:20:30.923008 INFO::Adjusting the last 2 batches sequentially
#> 2024-10-29 16:20:30.923907 INFO::Adjusting sequential tree level 1 with 2 batches
#> 2024-10-29 16:20:30.976354 INFO::Done
#> 2024-10-29 16:20:30.977184 INFO::Acquiring quality metrics after batch effect correction.
#> 2024-10-29 16:20:30.980988 INFO::ASW Batch was 0.00351969359420185 prior to batch effect correction and is now -0.0986215727352535 .
#> 2024-10-29 16:20:30.981551 INFO::ASW Label was -0.00500673314007958 prior to batch effect correction and is now 0.0118986877594866 .
#> 2024-10-29 16:20:30.98232 INFO::Total function execution time is  0.111027002334595  s and adjustment time is  0.0548160076141357 s ( 49.37 )

6.4 BERT with Covariables

BERT can utilize categorical covariables that are specified in columns Cov_1, Cov_2, .... These columns are automatically detected and integrated into the batch effect correction process.

# import BERT
library(BERT)
# set seed for reproducibility
set.seed(1)
# generate data with 5 batches, 60 features, 30 samples per batch, 15% missing values and 2 classes
dataset_raw <- generate_dataset(features=60, batches=5, samplesperbatch=30, mvstmt=0.15, classes=2)
# create covariable column with 2 possible values, e.g. male/female condition
dataset_raw["Cov_1"] = sample(c(1,2), size=dim(dataset_raw)[1], replace=TRUE)
# BERT
dataset_adjusted <- BERT(dataset_raw)
#> 2024-10-29 16:20:31.020476 INFO::Formatting Data.
#> 2024-10-29 16:20:31.021253 INFO::Replacing NaNs with NAs.
#> 2024-10-29 16:20:31.022267 INFO::Removing potential empty rows and columns
#> 2024-10-29 16:20:31.024028 INFO::Found  1350  missing values.
#> 2024-10-29 16:20:31.024894 INFO::BERT requires at least 2 numeric values per batch/covariate level. This may reduce the number of adjustable features considerably, depending on the quantification technique.
#> 2024-10-29 16:20:31.041695 INFO::Introduced 0 missing values due to singular proteins at batch/covariate level.
#> 2024-10-29 16:20:31.042365 INFO::Done
#> 2024-10-29 16:20:31.042897 INFO::Acquiring quality metrics before batch effect correction.
#> 2024-10-29 16:20:31.047509 INFO::Starting hierarchical adjustment
#> 2024-10-29 16:20:31.048212 INFO::Found  5  batches.
#> 2024-10-29 16:20:31.048691 INFO::Cores argument is not defined or BPPARAM has been specified. Argument corereduction will not be used.
#> 2024-10-29 16:20:31.049278 INFO::Using default BPPARAM
#> 2024-10-29 16:20:31.049776 INFO::Processing subtree level 1
#> 2024-10-29 16:20:31.331149 INFO::Adjusting the last 2 batches sequentially
#> 2024-10-29 16:20:31.332899 INFO::Adjusting sequential tree level 1 with 2 batches
#> 2024-10-29 16:20:31.41336 INFO::Done
#> 2024-10-29 16:20:31.414071 INFO::Acquiring quality metrics after batch effect correction.
#> 2024-10-29 16:20:31.419254 INFO::ASW Batch was 0.492773245691086 prior to batch effect correction and is now -0.0377157224767566 .
#> 2024-10-29 16:20:31.419925 INFO::ASW Label was 0.40854766060101 prior to batch effect correction and is now 0.895560693013661 .
#> 2024-10-29 16:20:31.420797 INFO::Total function execution time is  0.4003746509552  s and adjustment time is  0.365240573883057 s ( 91.22 )

6.5 BERT with references

In rare cases, class distributions across experiments may be severely skewed. In particular, a batch might contain classes that other batches don’t contain. In these cases, samples of common conditions may serve as references (bridges) between the batches (method="ref"). BERT utilizes those samples as references that have a condition specified in the “Reference” column of the input. All other samples are co-adjusted. Please note, that this strategy implicitly uses limma as underlying batch effect correction algorithm.

# import BERT
library(BERT)
# generate data with 4 batches, 6 features, 15 samples per batch, 15% missing values and 2 classes
dataset_raw <- generate_dataset(features=6, batches=4, samplesperbatch=15, mvstmt=0.15, classes=2)
# create reference column with default value 0.  The 0 indicates, that the respective sample should be co-adjusted only.
dataset_raw[, "Reference"] <- 0
# randomly select 2 references per batch and class - in practice, this choice will be determined by external requirements (e.g. class known for only these samples)
batches <- unique(dataset_raw$Batch) # all the batches
for(b in batches){ # iterate over all batches
    # references from class 1
    ref_idx = sample(which((dataset_raw$Batch==b)&(dataset_raw$Label==1)), size=2, replace=FALSE)
    dataset_raw[ref_idx, "Reference"] <- 1
    # references from class 2
    ref_idx = sample(which((dataset_raw$Batch==b)&(dataset_raw$Label==2)), size=2, replace=FALSE)
    dataset_raw[ref_idx, "Reference"] <- 2
}
# BERT
dataset_adjusted <- BERT(dataset_raw, method="ref")
#> 2024-10-29 16:20:31.470564 INFO::Formatting Data.
#> 2024-10-29 16:20:31.471397 INFO::Replacing NaNs with NAs.
#> 2024-10-29 16:20:31.472346 INFO::Removing potential empty rows and columns
#> 2024-10-29 16:20:31.473378 INFO::Found  60  missing values.
#> 2024-10-29 16:20:31.477126 INFO::Introduced 0 missing values due to singular proteins at batch/covariate level.
#> 2024-10-29 16:20:31.477671 INFO::Done
#> 2024-10-29 16:20:31.478181 INFO::Acquiring quality metrics before batch effect correction.
#> 2024-10-29 16:20:31.481143 INFO::Starting hierarchical adjustment
#> 2024-10-29 16:20:31.48187 INFO::Found  4  batches.
#> 2024-10-29 16:20:31.482409 INFO::Cores argument is not defined or BPPARAM has been specified. Argument corereduction will not be used.
#> 2024-10-29 16:20:31.48304 INFO::Using default BPPARAM
#> 2024-10-29 16:20:31.483562 INFO::Processing subtree level 1
#> 2024-10-29 16:20:31.601148 INFO::Adjusting the last 2 batches sequentially
#> 2024-10-29 16:20:31.602822 INFO::Adjusting sequential tree level 1 with 2 batches
#> 2024-10-29 16:20:31.624366 INFO::Done
#> 2024-10-29 16:20:31.625073 INFO::Acquiring quality metrics after batch effect correction.
#> 2024-10-29 16:20:31.628602 INFO::ASW Batch was 0.440355021914032 prior to batch effect correction and is now -0.087480278736629 .
#> 2024-10-29 16:20:31.629246 INFO::ASW Label was 0.373906827748893 prior to batch effect correction and is now 0.919791677398366 .
#> 2024-10-29 16:20:31.630076 INFO::Total function execution time is  0.159577131271362  s and adjustment time is  0.142629861831665 s ( 89.38 )

7 Issues

Issues can be reported in the GitHub forum, the BioConductor forum or directly to the authors.

8 License

This code is published under the GPLv3.0 License and is available for non-commercial academic purposes.

9 Reference

Please cite our manuscript, if you use BERT for your research: Yannis Schumann, Simon Schlumbohm et al., BERT - Batch Effect Reduction Trees with Missing Value Tolerance, 2023

10 Session Info

sessionInfo()
#> R Under development (unstable) (2024-10-21 r87258)
#> Platform: x86_64-pc-linux-gnu
#> Running under: Ubuntu 24.04.1 LTS
#> 
#> Matrix products: default
#> BLAS:   /home/biocbuild/bbs-3.21-bioc/R/lib/libRblas.so 
#> LAPACK: /usr/lib/x86_64-linux-gnu/lapack/liblapack.so.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] stats     graphics  grDevices utils     datasets  methods   base     
#> 
#> other attached packages:
#> [1] BERT_1.3.0       BiocStyle_2.35.0
#> 
#> loaded via a namespace (and not attached):
#>  [1] tidyselect_1.2.1            blob_1.2.4                 
#>  [3] Biostrings_2.75.0           fastmap_1.2.0              
#>  [5] janitor_2.2.0               XML_3.99-0.17              
#>  [7] digest_0.6.37               timechange_0.3.0           
#>  [9] lifecycle_1.0.4             cluster_2.1.6              
#> [11] statmod_1.5.0               survival_3.7-0             
#> [13] KEGGREST_1.47.0             invgamma_1.1               
#> [15] RSQLite_2.3.7               magrittr_2.0.3             
#> [17] genefilter_1.89.0           compiler_4.5.0             
#> [19] rlang_1.1.4                 sass_0.4.9                 
#> [21] tools_4.5.0                 yaml_2.3.10                
#> [23] knitr_1.48                  S4Arrays_1.7.0             
#> [25] bit_4.5.0                   DelayedArray_0.33.0        
#> [27] abind_1.4-8                 BiocParallel_1.41.0        
#> [29] BiocGenerics_0.53.0         grid_4.5.0                 
#> [31] stats4_4.5.0                xtable_1.8-4               
#> [33] edgeR_4.5.0                 iterators_1.0.14           
#> [35] logging_0.10-108            SummarizedExperiment_1.37.0
#> [37] cli_3.6.3                   rmarkdown_2.28             
#> [39] crayon_1.5.3                generics_0.1.3             
#> [41] httr_1.4.7                  DBI_1.2.3                  
#> [43] cachem_1.1.0                stringr_1.5.1              
#> [45] zlibbioc_1.53.0             splines_4.5.0              
#> [47] parallel_4.5.0              AnnotationDbi_1.69.0       
#> [49] BiocManager_1.30.25         XVector_0.47.0             
#> [51] matrixStats_1.4.1           vctrs_0.6.5                
#> [53] Matrix_1.7-1                jsonlite_1.8.9             
#> [55] sva_3.55.0                  bookdown_0.41              
#> [57] comprehenr_0.6.10           IRanges_2.41.0             
#> [59] S4Vectors_0.45.0            bit64_4.5.2                
#> [61] locfit_1.5-9.10             foreach_1.5.2              
#> [63] limma_3.63.0                jquerylib_0.1.4            
#> [65] annotate_1.85.0             glue_1.8.0                 
#> [67] codetools_0.2-20            lubridate_1.9.3            
#> [69] stringi_1.8.4               GenomeInfoDb_1.43.0        
#> [71] GenomicRanges_1.59.0        UCSC.utils_1.3.0           
#> [73] htmltools_0.5.8.1           GenomeInfoDbData_1.2.13    
#> [75] R6_2.5.1                    evaluate_1.0.1             
#> [77] lattice_0.22-6              Biobase_2.67.0             
#> [79] png_0.1-8                   memoise_2.0.1              
#> [81] snakecase_0.11.1            bslib_0.8.0                
#> [83] SparseArray_1.7.0           nlme_3.1-166               
#> [85] mgcv_1.9-1                  xfun_0.48                  
#> [87] MatrixGenerics_1.19.0