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.0913860 -1.9432641     1
#> 2 -0.6408452 -1.7986956     1
#> 3  2.2912134  1.6272827     2
#> 4  2.0210501  0.1256088     2
#> 5 -0.8585913  0.7856018     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)
#> 2025-01-29 16:16:49.126329 INFO::Formatting Data.
#> 2025-01-29 16:16:49.138131 INFO::Replacing NaNs with NAs.
#> 2025-01-29 16:16:49.149923 INFO::Removing potential empty rows and columns
#> 2025-01-29 16:16:49.478523 INFO::Found  600  missing values.
#> 2025-01-29 16:16:49.491273 INFO::Introduced 0 missing values due to singular proteins at batch/covariate level.
#> 2025-01-29 16:16:49.492005 INFO::Done
#> 2025-01-29 16:16:49.492553 INFO::Acquiring quality metrics before batch effect correction.
#> 2025-01-29 16:16:49.505832 INFO::Starting hierarchical adjustment
#> 2025-01-29 16:16:49.506808 INFO::Found  10  batches.
#> 2025-01-29 16:16:49.507365 INFO::Cores argument is not defined or BPPARAM has been specified. Argument corereduction will not be used.
#> 2025-01-29 16:16:52.667855 INFO::Using default BPPARAM
#> 2025-01-29 16:16:52.668701 INFO::Processing subtree level 1
#> 2025-01-29 16:16:54.642114 INFO::Processing subtree level 2
#> 2025-01-29 16:16:56.621754 INFO::Adjusting the last 1 batches sequentially
#> 2025-01-29 16:16:56.624384 INFO::Done
#> 2025-01-29 16:16:56.62526 INFO::Acquiring quality metrics after batch effect correction.
#> 2025-01-29 16:16:56.632087 INFO::ASW Batch was 0.451095164317474 prior to batch effect correction and is now -0.0752636523392637 .
#> 2025-01-29 16:16:56.632958 INFO::ASW Label was 0.367513743669761 prior to batch effect correction and is now 0.827875044172967 .
#> 2025-01-29 16:16:56.6344 INFO::Total function execution time is  7.53445792198181  s and adjustment time is  7.11786866188049 s ( 94.47 )

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")
#> 2025-01-29 16:16:56.709496 INFO::Formatting Data.
#> 2025-01-29 16:16:56.710316 INFO::Replacing NaNs with NAs.
#> 2025-01-29 16:16:56.711489 INFO::Removing potential empty rows and columns
#> 2025-01-29 16:16:56.713957 INFO::Found  2700  missing values.
#> 2025-01-29 16:16:56.739703 INFO::Introduced 0 missing values due to singular proteins at batch/covariate level.
#> 2025-01-29 16:16:56.740381 INFO::Done
#> 2025-01-29 16:16:56.740904 INFO::Acquiring quality metrics before batch effect correction.
#> 2025-01-29 16:16:56.755536 INFO::Starting hierarchical adjustment
#> 2025-01-29 16:16:56.756343 INFO::Found  20  batches.
#> 2025-01-29 16:16:56.756892 INFO::Cores argument is not defined or BPPARAM has been specified. Argument corereduction will not be used.
#> 2025-01-29 16:16:56.757524 INFO::Using default BPPARAM
#> 2025-01-29 16:16:56.758037 INFO::Processing subtree level 1
#> 2025-01-29 16:16:57.319492 INFO::Processing subtree level 2
#> 2025-01-29 16:16:57.783073 INFO::Processing subtree level 3
#> 2025-01-29 16:16:58.258726 INFO::Adjusting the last 1 batches sequentially
#> 2025-01-29 16:16:58.26117 INFO::Done
#> 2025-01-29 16:16:58.261979 INFO::Acquiring quality metrics after batch effect correction.
#> 2025-01-29 16:16:58.275673 INFO::ASW Batch was 0.458534670516476 prior to batch effect correction and is now -0.144051469540801 .
#> 2025-01-29 16:16:58.276549 INFO::ASW Label was 0.31095315944725 prior to batch effect correction and is now 0.825183693360462 .
#> 2025-01-29 16:16:58.277653 INFO::Total function execution time is  1.5681631565094  s and adjustment time is  1.50496411323547 s ( 95.97 )

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)
#> 2025-01-29 16:16:58.325869 INFO::Formatting Data.
#> 2025-01-29 16:16:58.326882 INFO::Replacing NaNs with NAs.
#> 2025-01-29 16:16:58.328355 INFO::Removing potential empty rows and columns
#> 2025-01-29 16:16:58.331312 INFO::Found  2700  missing values.
#> 2025-01-29 16:16:58.359482 INFO::Introduced 0 missing values due to singular proteins at batch/covariate level.
#> 2025-01-29 16:16:58.360178 INFO::Done
#> 2025-01-29 16:16:58.360732 INFO::Acquiring quality metrics before batch effect correction.
#> 2025-01-29 16:16:58.37422 INFO::Starting hierarchical adjustment
#> 2025-01-29 16:16:58.374996 INFO::Found  20  batches.
#> 2025-01-29 16:16:59.029207 INFO::Set up parallel execution backend with 2 workers
#> 2025-01-29 16:16:59.030599 INFO::Processing subtree level 1 with 20 batches using 2 cores.
#> 2025-01-29 16:17:02.183244 INFO::Adjusting the last 2 batches sequentially
#> 2025-01-29 16:17:02.184897 INFO::Adjusting sequential tree level 1 with 2 batches
#> 2025-01-29 16:17:03.898302 INFO::Done
#> 2025-01-29 16:17:03.899143 INFO::Acquiring quality metrics after batch effect correction.
#> 2025-01-29 16:17:03.909549 INFO::ASW Batch was 0.507962914856041 prior to batch effect correction and is now -0.107048275856229 .
#> 2025-01-29 16:17:03.910168 INFO::ASW Label was 0.260395642117206 prior to batch effect correction and is now 0.829866237323456 .
#> 2025-01-29 16:17:03.910998 INFO::Total function execution time is  5.58526468276978  s and adjustment time is  5.52294206619263 s ( 98.88 )

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")
#> 2025-01-29 16:17:03.992851 INFO::Formatting Data.
#> 2025-01-29 16:17:03.993601 INFO::Recognized SummarizedExperiment
#> 2025-01-29 16:17:03.994117 INFO::Typecasting input to dataframe.
#> 2025-01-29 16:17:04.030907 INFO::Replacing NaNs with NAs.
#> 2025-01-29 16:17:04.0321 INFO::Removing potential empty rows and columns
#> 2025-01-29 16:17:04.035321 INFO::Found  0  missing values.
#> 2025-01-29 16:17:04.041341 INFO::Introduced 0 missing values due to singular proteins at batch/covariate level.
#> 2025-01-29 16:17:04.041869 INFO::Done
#> 2025-01-29 16:17:04.042351 INFO::Acquiring quality metrics before batch effect correction.
#> 2025-01-29 16:17:04.046172 INFO::Starting hierarchical adjustment
#> 2025-01-29 16:17:04.046854 INFO::Found  2  batches.
#> 2025-01-29 16:17:04.047357 INFO::Cores argument is not defined or BPPARAM has been specified. Argument corereduction will not be used.
#> 2025-01-29 16:17:04.047904 INFO::Using default BPPARAM
#> 2025-01-29 16:17:04.048387 INFO::Adjusting the last 2 batches sequentially
#> 2025-01-29 16:17:04.049299 INFO::Adjusting sequential tree level 1 with 2 batches
#> 2025-01-29 16:17:04.099772 INFO::Done
#> 2025-01-29 16:17:04.100562 INFO::Acquiring quality metrics after batch effect correction.
#> 2025-01-29 16:17:04.104978 INFO::ASW Batch was -0.0131035864611407 prior to batch effect correction and is now -0.087808557020308 .
#> 2025-01-29 16:17:04.105637 INFO::ASW Label was -0.00372094715413252 prior to batch effect correction and is now 0.00922609598208422 .
#> 2025-01-29 16:17:04.106467 INFO::Total function execution time is  0.113599538803101  s and adjustment time is  0.0530133247375488 s ( 46.67 )

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)
#> 2025-01-29 16:17:04.1409 INFO::Formatting Data.
#> 2025-01-29 16:17:04.141735 INFO::Replacing NaNs with NAs.
#> 2025-01-29 16:17:04.142824 INFO::Removing potential empty rows and columns
#> 2025-01-29 16:17:04.145054 INFO::Found  1350  missing values.
#> 2025-01-29 16:17:04.146106 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.
#> 2025-01-29 16:17:04.203836 INFO::Introduced 0 missing values due to singular proteins at batch/covariate level.
#> 2025-01-29 16:17:04.204552 INFO::Done
#> 2025-01-29 16:17:04.205058 INFO::Acquiring quality metrics before batch effect correction.
#> 2025-01-29 16:17:04.209658 INFO::Starting hierarchical adjustment
#> 2025-01-29 16:17:04.210366 INFO::Found  5  batches.
#> 2025-01-29 16:17:04.210858 INFO::Cores argument is not defined or BPPARAM has been specified. Argument corereduction will not be used.
#> 2025-01-29 16:17:04.211418 INFO::Using default BPPARAM
#> 2025-01-29 16:17:04.211896 INFO::Processing subtree level 1
#> 2025-01-29 16:17:04.456633 INFO::Adjusting the last 2 batches sequentially
#> 2025-01-29 16:17:04.458681 INFO::Adjusting sequential tree level 1 with 2 batches
#> 2025-01-29 16:17:04.522429 INFO::Done
#> 2025-01-29 16:17:04.52314 INFO::Acquiring quality metrics after batch effect correction.
#> 2025-01-29 16:17:04.528571 INFO::ASW Batch was 0.492773245691086 prior to batch effect correction and is now -0.0377157224767566 .
#> 2025-01-29 16:17:04.529219 INFO::ASW Label was 0.40854766060101 prior to batch effect correction and is now 0.895560693013661 .
#> 2025-01-29 16:17:04.530086 INFO::Total function execution time is  0.389231204986572  s and adjustment time is  0.31217360496521 s ( 80.2 )

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")
#> 2025-01-29 16:17:04.585952 INFO::Formatting Data.
#> 2025-01-29 16:17:04.586899 INFO::Replacing NaNs with NAs.
#> 2025-01-29 16:17:04.588005 INFO::Removing potential empty rows and columns
#> 2025-01-29 16:17:04.589225 INFO::Found  60  missing values.
#> 2025-01-29 16:17:04.595011 INFO::Introduced 0 missing values due to singular proteins at batch/covariate level.
#> 2025-01-29 16:17:04.595677 INFO::Done
#> 2025-01-29 16:17:04.596326 INFO::Acquiring quality metrics before batch effect correction.
#> 2025-01-29 16:17:04.600054 INFO::Starting hierarchical adjustment
#> 2025-01-29 16:17:04.600974 INFO::Found  4  batches.
#> 2025-01-29 16:17:04.601642 INFO::Cores argument is not defined or BPPARAM has been specified. Argument corereduction will not be used.
#> 2025-01-29 16:17:04.602387 INFO::Using default BPPARAM
#> 2025-01-29 16:17:04.603059 INFO::Processing subtree level 1
#> 2025-01-29 16:17:04.721526 INFO::Adjusting the last 2 batches sequentially
#> 2025-01-29 16:17:04.723184 INFO::Adjusting sequential tree level 1 with 2 batches
#> 2025-01-29 16:17:04.743386 INFO::Done
#> 2025-01-29 16:17:04.744053 INFO::Acquiring quality metrics after batch effect correction.
#> 2025-01-29 16:17:04.747287 INFO::ASW Batch was 0.440355021914032 prior to batch effect correction and is now -0.087480278736629 .
#> 2025-01-29 16:17:04.747879 INFO::ASW Label was 0.373906827748893 prior to batch effect correction and is now 0.919791677398366 .
#> 2025-01-29 16:17:04.748674 INFO::Total function execution time is  0.16278862953186  s and adjustment time is  0.142592430114746 s ( 87.59 )

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 Under development (unstable) (2025-01-20 r87609)
#> 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  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.3.6       BiocStyle_2.35.0
#> 
#> loaded via a namespace (and not attached):
#>  [1] tidyselect_1.2.1            blob_1.2.4                 
#>  [3] Biostrings_2.75.3           fastmap_1.2.0              
#>  [5] janitor_2.2.1               XML_3.99-0.18              
#>  [7] digest_0.6.37               timechange_0.3.0           
#>  [9] lifecycle_1.0.4             cluster_2.1.8              
#> [11] survival_3.8-3              statmod_1.5.0              
#> [13] KEGGREST_1.47.0             invgamma_1.1               
#> [15] RSQLite_2.3.9               magrittr_2.0.3             
#> [17] genefilter_1.89.0           compiler_4.5.0             
#> [19] rlang_1.1.5                 sass_0.4.9                 
#> [21] tools_4.5.0                 yaml_2.3.10                
#> [23] knitr_1.49                  S4Arrays_1.7.1             
#> [25] bit_4.5.0.1                 DelayedArray_0.33.4        
#> [27] abind_1.4-8                 BiocParallel_1.41.0        
#> [29] BiocGenerics_0.53.6         grid_4.5.0                 
#> [31] stats4_4.5.0                xtable_1.8-4               
#> [33] edgeR_4.5.2                 iterators_1.0.14           
#> [35] logging_0.10-108            SummarizedExperiment_1.37.0
#> [37] cli_3.6.3                   rmarkdown_2.29             
#> [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] splines_4.5.0               parallel_4.5.0             
#> [47] AnnotationDbi_1.69.0        BiocManager_1.30.25        
#> [49] XVector_0.47.2              matrixStats_1.5.0          
#> [51] vctrs_0.6.5                 Matrix_1.7-2               
#> [53] jsonlite_1.8.9              sva_3.55.0                 
#> [55] bookdown_0.42               comprehenr_0.6.10          
#> [57] IRanges_2.41.2              S4Vectors_0.45.2           
#> [59] bit64_4.6.0-1               locfit_1.5-9.10            
#> [61] foreach_1.5.2               limma_3.63.3               
#> [63] jquerylib_0.1.4             annotate_1.85.0            
#> [65] glue_1.8.0                  codetools_0.2-20           
#> [67] lubridate_1.9.4             stringi_1.8.4              
#> [69] GenomeInfoDb_1.43.4         GenomicRanges_1.59.1       
#> [71] UCSC.utils_1.3.1            htmltools_0.5.8.1          
#> [73] GenomeInfoDbData_1.2.13     R6_2.5.1                   
#> [75] evaluate_1.0.3              lattice_0.22-6             
#> [77] Biobase_2.67.0              png_0.1-8                  
#> [79] memoise_2.0.1               snakecase_0.11.1           
#> [81] bslib_0.8.0                 SparseArray_1.7.4          
#> [83] nlme_3.1-167                mgcv_1.9-1                 
#> [85] xfun_0.50                   MatrixGenerics_1.19.1