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 -0.93537581  0.2713915     1
#> 2  1.25412899 -0.8464100     1
#> 3  0.70129767  0.6022083     2
#> 4  0.09470375  1.4189676     2
#> 5 -0.24315261 -1.9582269     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)
#> 2026-04-28 17:26:14.90186 INFO::Formatting Data.
#> 2026-04-28 17:26:14.913654 INFO::Replacing NaNs with NAs.
#> 2026-04-28 17:26:14.925164 INFO::Removing potential empty rows and columns
#> 2026-04-28 17:26:15.210462 INFO::Found  600  missing values.
#> 2026-04-28 17:26:15.222411 INFO::Introduced 0 missing values due to singular proteins at batch/covariate level.
#> 2026-04-28 17:26:15.223101 INFO::Done
#> 2026-04-28 17:26:15.223644 INFO::Acquiring quality metrics before batch effect correction.
#> 2026-04-28 17:26:15.235325 INFO::Starting hierarchical adjustment
#> 2026-04-28 17:26:15.236213 INFO::Found  10  batches.
#> 2026-04-28 17:26:15.23676 INFO::Cores argument is not defined or BPPARAM has been specified. Argument corereduction will not be used.
#> 2026-04-28 17:26:17.271945 INFO::Using default BPPARAM
#> 2026-04-28 17:26:17.272714 INFO::Processing subtree level 1
#> 2026-04-28 17:26:19.189497 INFO::Processing subtree level 2
#> 2026-04-28 17:26:21.201407 INFO::Adjusting the last 1 batches sequentially
#> 2026-04-28 17:26:21.204412 INFO::Done
#> 2026-04-28 17:26:21.205545 INFO::Acquiring quality metrics after batch effect correction.
#> 2026-04-28 17:26:21.21435 INFO::ASW Batch was 0.462382967432971 prior to batch effect correction and is now -0.11787234972041 .
#> 2026-04-28 17:26:21.215655 INFO::ASW Label was 0.352368890394986 prior to batch effect correction and is now 0.889748609538338 .
#> 2026-04-28 17:26:21.217512 INFO::Total function execution time is  6.34682512283325  s and adjustment time is  5.96841764450073 s ( 94.04 )

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")
#> 2026-04-28 17:26:21.302089 INFO::Formatting Data.
#> 2026-04-28 17:26:21.303313 INFO::Replacing NaNs with NAs.
#> 2026-04-28 17:26:21.30499 INFO::Removing potential empty rows and columns
#> 2026-04-28 17:26:21.308243 INFO::Found  2700  missing values.
#> 2026-04-28 17:26:21.33755 INFO::Introduced 0 missing values due to singular proteins at batch/covariate level.
#> 2026-04-28 17:26:21.33833 INFO::Done
#> 2026-04-28 17:26:21.338964 INFO::Acquiring quality metrics before batch effect correction.
#> 2026-04-28 17:26:21.352596 INFO::Starting hierarchical adjustment
#> 2026-04-28 17:26:21.353489 INFO::Found  20  batches.
#> 2026-04-28 17:26:21.354108 INFO::Cores argument is not defined or BPPARAM has been specified. Argument corereduction will not be used.
#> 2026-04-28 17:26:21.354779 INFO::Using default BPPARAM
#> 2026-04-28 17:26:21.355392 INFO::Processing subtree level 1
#> 2026-04-28 17:26:21.860679 INFO::Processing subtree level 2
#> 2026-04-28 17:26:22.272945 INFO::Processing subtree level 3
#> 2026-04-28 17:26:22.725239 INFO::Adjusting the last 1 batches sequentially
#> 2026-04-28 17:26:22.72778 INFO::Done
#> 2026-04-28 17:26:22.728697 INFO::Acquiring quality metrics after batch effect correction.
#> 2026-04-28 17:26:22.743573 INFO::ASW Batch was 0.43650308569708 prior to batch effect correction and is now -0.110179672157025 .
#> 2026-04-28 17:26:22.744638 INFO::ASW Label was 0.347397417777328 prior to batch effect correction and is now 0.805051357146745 .
#> 2026-04-28 17:26:22.745891 INFO::Total function execution time is  1.44388246536255  s and adjustment time is  1.37442493438721 s ( 95.19 )

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)
#> 2026-04-28 17:26:22.807436 INFO::Formatting Data.
#> 2026-04-28 17:26:22.808564 INFO::Replacing NaNs with NAs.
#> 2026-04-28 17:26:22.810041 INFO::Removing potential empty rows and columns
#> 2026-04-28 17:26:22.813011 INFO::Found  2700  missing values.
#> 2026-04-28 17:26:22.843683 INFO::Introduced 0 missing values due to singular proteins at batch/covariate level.
#> 2026-04-28 17:26:22.844547 INFO::Done
#> 2026-04-28 17:26:22.845237 INFO::Acquiring quality metrics before batch effect correction.
#> 2026-04-28 17:26:22.859168 INFO::Starting hierarchical adjustment
#> 2026-04-28 17:26:22.860104 INFO::Found  20  batches.
#> 2026-04-28 17:26:23.578872 INFO::Set up parallel execution backend with 2 workers
#> 2026-04-28 17:26:23.580149 INFO::Processing subtree level 1 with 20 batches using 2 cores.
#> 2026-04-28 17:26:26.283669 INFO::Adjusting the last 2 batches sequentially
#> 2026-04-28 17:26:26.285109 INFO::Adjusting sequential tree level 1 with 2 batches
#> 2026-04-28 17:26:27.979598 INFO::Done
#> 2026-04-28 17:26:27.980464 INFO::Acquiring quality metrics after batch effect correction.
#> 2026-04-28 17:26:27.989588 INFO::ASW Batch was 0.465778128902812 prior to batch effect correction and is now -0.127738330821995 .
#> 2026-04-28 17:26:27.990254 INFO::ASW Label was 0.323364094031014 prior to batch effect correction and is now 0.813590475952723 .
#> 2026-04-28 17:26:27.991066 INFO::Total function execution time is  5.18379402160645  s and adjustment time is  5.11928915977478 s ( 98.76 )

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")
#> 2026-04-28 17:26:28.047589 INFO::Formatting Data.
#> 2026-04-28 17:26:28.048374 INFO::Recognized SummarizedExperiment
#> 2026-04-28 17:26:28.048929 INFO::Typecasting input to dataframe.
#> 2026-04-28 17:26:28.081824 INFO::Replacing NaNs with NAs.
#> 2026-04-28 17:26:28.083008 INFO::Removing potential empty rows and columns
#> 2026-04-28 17:26:28.086268 INFO::Found  0  missing values.
#> 2026-04-28 17:26:28.09236 INFO::Introduced 0 missing values due to singular proteins at batch/covariate level.
#> 2026-04-28 17:26:28.09291 INFO::Done
#> 2026-04-28 17:26:28.093437 INFO::Acquiring quality metrics before batch effect correction.
#> 2026-04-28 17:26:28.09717 INFO::Starting hierarchical adjustment
#> 2026-04-28 17:26:28.097853 INFO::Found  2  batches.
#> 2026-04-28 17:26:28.0984 INFO::Cores argument is not defined or BPPARAM has been specified. Argument corereduction will not be used.
#> 2026-04-28 17:26:28.099 INFO::Using default BPPARAM
#> 2026-04-28 17:26:28.099524 INFO::Adjusting the last 2 batches sequentially
#> 2026-04-28 17:26:28.100457 INFO::Adjusting sequential tree level 1 with 2 batches
#> 2026-04-28 17:26:28.144145 INFO::Done
#> 2026-04-28 17:26:28.144835 INFO::Acquiring quality metrics after batch effect correction.
#> 2026-04-28 17:26:28.148568 INFO::ASW Batch was -0.00705317567571381 prior to batch effect correction and is now -0.093615315447412 .
#> 2026-04-28 17:26:28.149147 INFO::ASW Label was 0.0106206821097985 prior to batch effect correction and is now 0.0278813778229404 .
#> 2026-04-28 17:26:28.149896 INFO::Total function execution time is  0.102298974990845  s and adjustment time is  0.0464060306549072 s ( 45.36 )

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)
#> 2026-04-28 17:26:28.196873 INFO::Formatting Data.
#> 2026-04-28 17:26:28.197686 INFO::Replacing NaNs with NAs.
#> 2026-04-28 17:26:28.198644 INFO::Removing potential empty rows and columns
#> 2026-04-28 17:26:28.200443 INFO::Found  1350  missing values.
#> 2026-04-28 17:26:28.20132 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.
#> 2026-04-28 17:26:28.21806 INFO::Introduced 0 missing values due to singular proteins at batch/covariate level.
#> 2026-04-28 17:26:28.218749 INFO::Done
#> 2026-04-28 17:26:28.219337 INFO::Acquiring quality metrics before batch effect correction.
#> 2026-04-28 17:26:28.224098 INFO::Starting hierarchical adjustment
#> 2026-04-28 17:26:28.224803 INFO::Found  5  batches.
#> 2026-04-28 17:26:28.225332 INFO::Cores argument is not defined or BPPARAM has been specified. Argument corereduction will not be used.
#> 2026-04-28 17:26:28.225922 INFO::Using default BPPARAM
#> 2026-04-28 17:26:28.226481 INFO::Processing subtree level 1
#> 2026-04-28 17:26:28.456777 INFO::Adjusting the last 2 batches sequentially
#> 2026-04-28 17:26:28.458472 INFO::Adjusting sequential tree level 1 with 2 batches
#> 2026-04-28 17:26:28.512192 INFO::Done
#> 2026-04-28 17:26:28.512999 INFO::Acquiring quality metrics after batch effect correction.
#> 2026-04-28 17:26:28.518444 INFO::ASW Batch was 0.492773245691086 prior to batch effect correction and is now -0.0377157224767566 .
#> 2026-04-28 17:26:28.519132 INFO::ASW Label was 0.40854766060101 prior to batch effect correction and is now 0.895560693013661 .
#> 2026-04-28 17:26:28.520018 INFO::Total function execution time is  0.323177814483643  s and adjustment time is  0.287485122680664 s ( 88.96 )

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")
#> 2026-04-28 17:26:28.570582 INFO::Formatting Data.
#> 2026-04-28 17:26:28.571504 INFO::Replacing NaNs with NAs.
#> 2026-04-28 17:26:28.572506 INFO::Removing potential empty rows and columns
#> 2026-04-28 17:26:28.573622 INFO::Found  60  missing values.
#> 2026-04-28 17:26:28.577767 INFO::Introduced 0 missing values due to singular proteins at batch/covariate level.
#> 2026-04-28 17:26:28.578376 INFO::Done
#> 2026-04-28 17:26:28.578985 INFO::Acquiring quality metrics before batch effect correction.
#> 2026-04-28 17:26:28.582221 INFO::Starting hierarchical adjustment
#> 2026-04-28 17:26:28.583066 INFO::Found  4  batches.
#> 2026-04-28 17:26:28.583724 INFO::Cores argument is not defined or BPPARAM has been specified. Argument corereduction will not be used.
#> 2026-04-28 17:26:28.58444 INFO::Using default BPPARAM
#> 2026-04-28 17:26:28.58515 INFO::Processing subtree level 1
#> 2026-04-28 17:26:28.697091 INFO::Adjusting the last 2 batches sequentially
#> 2026-04-28 17:26:28.699015 INFO::Adjusting sequential tree level 1 with 2 batches
#> 2026-04-28 17:26:28.72175 INFO::Done
#> 2026-04-28 17:26:28.722556 INFO::Acquiring quality metrics after batch effect correction.
#> 2026-04-28 17:26:28.726007 INFO::ASW Batch was 0.440355021914032 prior to batch effect correction and is now -0.087480278736629 .
#> 2026-04-28 17:26:28.72665 INFO::ASW Label was 0.373906827748893 prior to batch effect correction and is now 0.919791677398366 .
#> 2026-04-28 17:26:28.727577 INFO::Total function execution time is  0.157063961029053  s and adjustment time is  0.138847827911377 s ( 88.4 )

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: Schumann Y, Gocke A, Neumann J (2024). Computational Methods for Data Integration and Imputation of Missing Values in Omics Datasets. PROTEOMICS. ISSN 1615-9861, doi:10.1002/pmic.202400100

10 Session Info

sessionInfo()
#> R version 4.6.0 RC (2026-04-17 r89917)
#> Platform: x86_64-pc-linux-gnu
#> Running under: Ubuntu 24.04.4 LTS
#> 
#> Matrix products: default
#> BLAS:   /home/biocbuild/bbs-3.24-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] stats     graphics  grDevices utils     datasets  methods   base     
#> 
#> other attached packages:
#> [1] BERT_1.9.0       BiocStyle_2.41.0
#> 
#> loaded via a namespace (and not attached):
#>  [1] KEGGREST_1.53.0             SummarizedExperiment_1.43.0
#>  [3] xfun_0.57                   bslib_0.10.0               
#>  [5] Biobase_2.73.0              lattice_0.22-9             
#>  [7] vctrs_0.7.3                 tools_4.6.0                
#>  [9] generics_0.1.4              stats4_4.6.0               
#> [11] parallel_4.6.0              AnnotationDbi_1.75.0       
#> [13] RSQLite_2.4.6               cluster_2.1.8.2            
#> [15] blob_1.3.0                  logging_0.10-108           
#> [17] Matrix_1.7-5                S4Vectors_0.51.0           
#> [19] lifecycle_1.0.5             compiler_4.6.0             
#> [21] stringr_1.6.0               Biostrings_2.81.0          
#> [23] statmod_1.5.1               janitor_2.2.1              
#> [25] Seqinfo_1.3.0               codetools_0.2-20           
#> [27] snakecase_0.11.1            htmltools_0.5.9            
#> [29] sass_0.4.10                 yaml_2.3.12                
#> [31] crayon_1.5.3                jquerylib_0.1.4            
#> [33] comprehenr_0.6.10           BiocParallel_1.47.0        
#> [35] limma_3.69.0                DelayedArray_0.39.0        
#> [37] cachem_1.1.0                iterators_1.0.14           
#> [39] abind_1.4-8                 foreach_1.5.2              
#> [41] nlme_3.1-169                sva_3.61.0                 
#> [43] genefilter_1.95.0           locfit_1.5-9.12            
#> [45] tidyselect_1.2.1            digest_0.6.39              
#> [47] stringi_1.8.7               bookdown_0.46              
#> [49] splines_4.6.0               fastmap_1.2.0              
#> [51] grid_4.6.0                  cli_3.6.6                  
#> [53] invgamma_1.2                SparseArray_1.13.0         
#> [55] magrittr_2.0.5              S4Arrays_1.13.0            
#> [57] survival_3.8-6              XML_3.99-0.23              
#> [59] edgeR_4.11.0                bit64_4.8.0                
#> [61] lubridate_1.9.5             timechange_0.4.0           
#> [63] rmarkdown_2.31              XVector_0.53.0             
#> [65] httr_1.4.8                  matrixStats_1.5.0          
#> [67] bit_4.6.0                   otel_0.2.0                 
#> [69] png_0.1-9                   memoise_2.0.1              
#> [71] evaluate_1.0.5              knitr_1.51                 
#> [73] GenomicRanges_1.65.0        IRanges_2.47.0             
#> [75] mgcv_1.9-4                  rlang_1.2.0                
#> [77] xtable_1.8-8                glue_1.8.1                 
#> [79] DBI_1.3.0                   BiocManager_1.30.27        
#> [81] BiocGenerics_0.59.0         annotate_1.91.0            
#> [83] jsonlite_2.0.0              R6_2.6.1                   
#> [85] MatrixGenerics_1.25.0