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.6698590  0.5651831     1
#> 2  0.2616950 -0.5436112     1
#> 3 -0.9850240 -0.1895470     2
#> 4 -0.5652119  0.2611247     2
#> 5  0.5807294  1.4047859     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-07-22 16:18:07.594312 INFO::Formatting Data.
#> 2024-07-22 16:18:07.604389 INFO::Replacing NaNs with NAs.
#> 2024-07-22 16:18:07.614344 INFO::Removing potential empty rows and columns
#> 2024-07-22 16:18:07.917954 INFO::Found  600  missing values.
#> 2024-07-22 16:18:07.932011 INFO::Introduced 0 missing values due to singular proteins at batch/covariate level.
#> 2024-07-22 16:18:07.932763 INFO::Done
#> 2024-07-22 16:18:07.933287 INFO::Acquiring quality metrics before batch effect correction.
#> 2024-07-22 16:18:07.948 INFO::Starting hierarchical adjustment
#> 2024-07-22 16:18:07.949078 INFO::Found  10  batches.
#> 2024-07-22 16:18:07.949643 INFO::Cores argument is not defined or BPPARAM has been specified. Argument corereduction will not be used.
#> 2024-07-22 16:18:10.829286 INFO::Using default BPPARAM
#> 2024-07-22 16:18:10.830021 INFO::Processing subtree level 1
#> 2024-07-22 16:18:13.046508 INFO::Processing subtree level 2
#> 2024-07-22 16:18:15.334677 INFO::Adjusting the last 1 batches sequentially
#> 2024-07-22 16:18:15.338541 INFO::Done
#> 2024-07-22 16:18:15.339927 INFO::Acquiring quality metrics after batch effect correction.
#> 2024-07-22 16:18:15.351106 INFO::ASW Batch was 0.505515563453848 prior to batch effect correction and is now -0.121018956695757 .
#> 2024-07-22 16:18:15.352603 INFO::ASW Label was 0.328316405680818 prior to batch effect correction and is now 0.747326631179575 .
#> 2024-07-22 16:18:15.35499 INFO::Total function execution time is  7.78552770614624  s and adjustment time is  7.38971924781799 s ( 94.92 )

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-07-22 16:18:15.443256 INFO::Formatting Data.
#> 2024-07-22 16:18:15.444118 INFO::Replacing NaNs with NAs.
#> 2024-07-22 16:18:15.445331 INFO::Removing potential empty rows and columns
#> 2024-07-22 16:18:15.447924 INFO::Found  2700  missing values.
#> 2024-07-22 16:18:15.477226 INFO::Introduced 0 missing values due to singular proteins at batch/covariate level.
#> 2024-07-22 16:18:15.478055 INFO::Done
#> 2024-07-22 16:18:15.478604 INFO::Acquiring quality metrics before batch effect correction.
#> 2024-07-22 16:18:15.493276 INFO::Starting hierarchical adjustment
#> 2024-07-22 16:18:15.494131 INFO::Found  20  batches.
#> 2024-07-22 16:18:15.494668 INFO::Cores argument is not defined or BPPARAM has been specified. Argument corereduction will not be used.
#> 2024-07-22 16:18:15.495315 INFO::Using default BPPARAM
#> 2024-07-22 16:18:15.495829 INFO::Processing subtree level 1
#> 2024-07-22 16:18:16.022643 INFO::Processing subtree level 2
#> 2024-07-22 16:18:16.442018 INFO::Processing subtree level 3
#> 2024-07-22 16:18:16.882033 INFO::Adjusting the last 1 batches sequentially
#> 2024-07-22 16:18:16.883852 INFO::Done
#> 2024-07-22 16:18:16.884475 INFO::Acquiring quality metrics after batch effect correction.
#> 2024-07-22 16:18:16.895134 INFO::ASW Batch was 0.463488509343107 prior to batch effect correction and is now -0.115387495045219 .
#> 2024-07-22 16:18:16.895808 INFO::ASW Label was 0.306896821573641 prior to batch effect correction and is now 0.829604428333249 .
#> 2024-07-22 16:18:16.896701 INFO::Total function execution time is  1.45353531837463  s and adjustment time is  1.38989043235779 s ( 95.62 )

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-07-22 16:18:16.935653 INFO::Formatting Data.
#> 2024-07-22 16:18:16.936425 INFO::Replacing NaNs with NAs.
#> 2024-07-22 16:18:16.937499 INFO::Removing potential empty rows and columns
#> 2024-07-22 16:18:16.939821 INFO::Found  2700  missing values.
#> 2024-07-22 16:18:16.964715 INFO::Introduced 0 missing values due to singular proteins at batch/covariate level.
#> 2024-07-22 16:18:16.965413 INFO::Done
#> 2024-07-22 16:18:16.965905 INFO::Acquiring quality metrics before batch effect correction.
#> 2024-07-22 16:18:16.979561 INFO::Starting hierarchical adjustment
#> 2024-07-22 16:18:16.980374 INFO::Found  20  batches.
#> 2024-07-22 16:18:17.781549 INFO::Set up parallel execution backend with 2 workers
#> 2024-07-22 16:18:17.782573 INFO::Processing subtree level 1 with 20 batches using 2 cores.
#> 2024-07-22 16:18:21.108607 INFO::Adjusting the last 2 batches sequentially
#> 2024-07-22 16:18:21.110009 INFO::Adjusting sequential tree level 1 with 2 batches
#> 2024-07-22 16:18:22.916483 INFO::Done
#> 2024-07-22 16:18:22.917198 INFO::Acquiring quality metrics after batch effect correction.
#> 2024-07-22 16:18:22.926368 INFO::ASW Batch was 0.444446758665952 prior to batch effect correction and is now -0.140031368878606 .
#> 2024-07-22 16:18:22.926919 INFO::ASW Label was 0.331329997813052 prior to batch effect correction and is now 0.831136387252085 .
#> 2024-07-22 16:18:22.927619 INFO::Total function execution time is  5.99206519126892  s and adjustment time is  5.93558883666992 s ( 99.06 )

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-07-22 16:18:22.998842 INFO::Formatting Data.
#> 2024-07-22 16:18:22.999464 INFO::Recognized SummarizedExperiment
#> 2024-07-22 16:18:22.999914 INFO::Typecasting input to dataframe.
#> 2024-07-22 16:18:23.031972 INFO::Replacing NaNs with NAs.
#> 2024-07-22 16:18:23.033029 INFO::Removing potential empty rows and columns
#> 2024-07-22 16:18:23.035845 INFO::Found  0  missing values.
#> 2024-07-22 16:18:23.041154 INFO::Introduced 0 missing values due to singular proteins at batch/covariate level.
#> 2024-07-22 16:18:23.041622 INFO::Done
#> 2024-07-22 16:18:23.042049 INFO::Acquiring quality metrics before batch effect correction.
#> 2024-07-22 16:18:23.045447 INFO::Starting hierarchical adjustment
#> 2024-07-22 16:18:23.046039 INFO::Found  2  batches.
#> 2024-07-22 16:18:23.046476 INFO::Cores argument is not defined or BPPARAM has been specified. Argument corereduction will not be used.
#> 2024-07-22 16:18:23.046967 INFO::Using default BPPARAM
#> 2024-07-22 16:18:23.047397 INFO::Adjusting the last 2 batches sequentially
#> 2024-07-22 16:18:23.048243 INFO::Adjusting sequential tree level 1 with 2 batches
#> 2024-07-22 16:18:23.088573 INFO::Done
#> 2024-07-22 16:18:23.089148 INFO::Acquiring quality metrics after batch effect correction.
#> 2024-07-22 16:18:23.092539 INFO::ASW Batch was 0.0079492752092062 prior to batch effect correction and is now -0.105217453217764 .
#> 2024-07-22 16:18:23.093026 INFO::ASW Label was -0.0224052552563604 prior to batch effect correction and is now -0.00863404231166552 .
#> 2024-07-22 16:18:23.093673 INFO::Total function execution time is  0.0948295593261719  s and adjustment time is  0.0426311492919922 s ( 44.96 )

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-07-22 16:18:23.125782 INFO::Formatting Data.
#> 2024-07-22 16:18:23.126481 INFO::Replacing NaNs with NAs.
#> 2024-07-22 16:18:23.127375 INFO::Removing potential empty rows and columns
#> 2024-07-22 16:18:23.12922 INFO::Found  1350  missing values.
#> 2024-07-22 16:18:23.13005 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-07-22 16:18:23.145467 INFO::Introduced 0 missing values due to singular proteins at batch/covariate level.
#> 2024-07-22 16:18:23.146034 INFO::Done
#> 2024-07-22 16:18:23.146477 INFO::Acquiring quality metrics before batch effect correction.
#> 2024-07-22 16:18:23.150973 INFO::Starting hierarchical adjustment
#> 2024-07-22 16:18:23.15159 INFO::Found  5  batches.
#> 2024-07-22 16:18:23.152065 INFO::Cores argument is not defined or BPPARAM has been specified. Argument corereduction will not be used.
#> 2024-07-22 16:18:23.152566 INFO::Using default BPPARAM
#> 2024-07-22 16:18:23.153001 INFO::Processing subtree level 1
#> 2024-07-22 16:18:23.375064 INFO::Adjusting the last 2 batches sequentially
#> 2024-07-22 16:18:23.376853 INFO::Adjusting sequential tree level 1 with 2 batches
#> 2024-07-22 16:18:23.435636 INFO::Done
#> 2024-07-22 16:18:23.436344 INFO::Acquiring quality metrics after batch effect correction.
#> 2024-07-22 16:18:23.441837 INFO::ASW Batch was 0.492773245691086 prior to batch effect correction and is now -0.0377157224767566 .
#> 2024-07-22 16:18:23.442483 INFO::ASW Label was 0.40854766060101 prior to batch effect correction and is now 0.895560693013661 .
#> 2024-07-22 16:18:23.443344 INFO::Total function execution time is  0.317615747451782  s and adjustment time is  0.284123182296753 s ( 89.46 )

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-07-22 16:18:23.496067 INFO::Formatting Data.
#> 2024-07-22 16:18:23.496827 INFO::Replacing NaNs with NAs.
#> 2024-07-22 16:18:23.497759 INFO::Removing potential empty rows and columns
#> 2024-07-22 16:18:23.498809 INFO::Found  60  missing values.
#> 2024-07-22 16:18:23.50303 INFO::Introduced 0 missing values due to singular proteins at batch/covariate level.
#> 2024-07-22 16:18:23.503603 INFO::Done
#> 2024-07-22 16:18:23.504174 INFO::Acquiring quality metrics before batch effect correction.
#> 2024-07-22 16:18:23.507486 INFO::Starting hierarchical adjustment
#> 2024-07-22 16:18:23.508271 INFO::Found  4  batches.
#> 2024-07-22 16:18:23.508832 INFO::Cores argument is not defined or BPPARAM has been specified. Argument corereduction will not be used.
#> 2024-07-22 16:18:23.509463 INFO::Using default BPPARAM
#> 2024-07-22 16:18:23.510031 INFO::Processing subtree level 1
#> 2024-07-22 16:18:23.622165 INFO::Adjusting the last 2 batches sequentially
#> 2024-07-22 16:18:23.624007 INFO::Adjusting sequential tree level 1 with 2 batches
#> 2024-07-22 16:18:23.646866 INFO::Done
#> 2024-07-22 16:18:23.647539 INFO::Acquiring quality metrics after batch effect correction.
#> 2024-07-22 16:18:23.65074 INFO::ASW Batch was 0.440355021914032 prior to batch effect correction and is now -0.087480278736629 .
#> 2024-07-22 16:18:23.651344 INFO::ASW Label was 0.373906827748893 prior to batch effect correction and is now 0.919791677398366 .
#> 2024-07-22 16:18:23.652126 INFO::Total function execution time is  0.156152248382568  s and adjustment time is  0.138733148574829 s ( 88.84 )

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 version 4.4.1 (2024-06-14)
#> Platform: x86_64-pc-linux-gnu
#> Running under: Ubuntu 22.04.4 LTS
#> 
#> Matrix products: default
#> BLAS:   /home/biocbuild/bbs-3.20-bioc/R/lib/libRblas.so 
#> LAPACK: /usr/lib/x86_64-linux-gnu/lapack/liblapack.so.3.10.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.1.2       BiocStyle_2.33.1
#> 
#> loaded via a namespace (and not attached):
#>  [1] tidyselect_1.2.1            blob_1.2.4                 
#>  [3] Biostrings_2.73.1           fastmap_1.2.0              
#>  [5] janitor_2.2.0               XML_3.99-0.17              
#>  [7] digest_0.6.36               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.45.1             invgamma_1.1               
#> [15] RSQLite_2.3.7               magrittr_2.0.3             
#> [17] genefilter_1.87.0           compiler_4.4.1             
#> [19] rlang_1.1.4                 sass_0.4.9                 
#> [21] tools_4.4.1                 yaml_2.3.9                 
#> [23] knitr_1.48                  S4Arrays_1.5.5             
#> [25] bit_4.0.5                   DelayedArray_0.31.9        
#> [27] abind_1.4-5                 BiocParallel_1.39.0        
#> [29] BiocGenerics_0.51.0         grid_4.4.1                 
#> [31] stats4_4.4.1                xtable_1.8-4               
#> [33] edgeR_4.3.5                 iterators_1.0.14           
#> [35] logging_0.10-108            SummarizedExperiment_1.35.1
#> [37] cli_3.6.3                   rmarkdown_2.27             
#> [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.51.1             splines_4.4.1              
#> [47] parallel_4.4.1              AnnotationDbi_1.67.0       
#> [49] BiocManager_1.30.23         XVector_0.45.0             
#> [51] matrixStats_1.3.0           vctrs_0.6.5                
#> [53] Matrix_1.7-0                jsonlite_1.8.8             
#> [55] sva_3.53.0                  bookdown_0.40              
#> [57] comprehenr_0.6.10           IRanges_2.39.2             
#> [59] S4Vectors_0.43.2            bit64_4.0.5                
#> [61] locfit_1.5-9.10             foreach_1.5.2              
#> [63] limma_3.61.5                jquerylib_0.1.4            
#> [65] annotate_1.83.0             glue_1.7.0                 
#> [67] codetools_0.2-20            lubridate_1.9.3            
#> [69] stringi_1.8.4               GenomeInfoDb_1.41.1        
#> [71] GenomicRanges_1.57.1        UCSC.utils_1.1.0           
#> [73] htmltools_0.5.8.1           GenomeInfoDbData_1.2.12    
#> [75] R6_2.5.1                    evaluate_0.24.0            
#> [77] lattice_0.22-6              Biobase_2.65.0             
#> [79] png_0.1-8                   memoise_2.0.1              
#> [81] snakecase_0.11.1            bslib_0.7.0                
#> [83] Rcpp_1.0.13                 SparseArray_1.5.25         
#> [85] nlme_3.1-165                mgcv_1.9-1                 
#> [87] xfun_0.46                   MatrixGenerics_1.17.0