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 -2.1132695 -0.184651533     1
#> 2 -0.5738389 -0.008044459     1
#> 3 -0.4600935  1.412635906     2
#> 4 -2.3411000 -1.776417156     2
#> 5  1.0293645  0.603843088     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-11-08 16:07:09.587428 INFO::Formatting Data.
#> 2024-11-08 16:07:09.590909 INFO::Replacing NaNs with NAs.
#> 2024-11-08 16:07:09.593834 INFO::Removing potential empty rows and columns
#> 2024-11-08 16:07:09.67868 INFO::Found  600  missing values.
#> 2024-11-08 16:07:09.682325 INFO::Introduced 0 missing values due to singular proteins at batch/covariate level.
#> 2024-11-08 16:07:09.682491 INFO::Done
#> 2024-11-08 16:07:09.68263 INFO::Acquiring quality metrics before batch effect correction.
#> 2024-11-08 16:07:09.686975 INFO::Starting hierarchical adjustment
#> 2024-11-08 16:07:09.687241 INFO::Found  10  batches.
#> 2024-11-08 16:07:09.687388 INFO::Cores argument is not defined or BPPARAM has been specified. Argument corereduction will not be used.
#> 2024-11-08 16:07:10.462109 INFO::Using default BPPARAM
#> 2024-11-08 16:07:10.462348 INFO::Processing subtree level 1
#> 2024-11-08 16:07:11.176469 INFO::Processing subtree level 2
#> 2024-11-08 16:07:11.872709 INFO::Adjusting the last 1 batches sequentially
#> 2024-11-08 16:07:11.874028 INFO::Done
#> 2024-11-08 16:07:11.874247 INFO::Acquiring quality metrics after batch effect correction.
#> 2024-11-08 16:07:11.880722 INFO::ASW Batch was 0.511027924139007 prior to batch effect correction and is now -0.131209130966805 .
#> 2024-11-08 16:07:11.881145 INFO::ASW Label was 0.299355003730648 prior to batch effect correction and is now 0.799824739927783 .
#> 2024-11-08 16:07:11.882175 INFO::Total function execution time is  2.30353498458862  s and adjustment time is  2.18687915802002 s ( 94.94 )

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-11-08 16:07:11.931217 INFO::Formatting Data.
#> 2024-11-08 16:07:11.931499 INFO::Replacing NaNs with NAs.
#> 2024-11-08 16:07:11.931974 INFO::Removing potential empty rows and columns
#> 2024-11-08 16:07:11.933249 INFO::Found  2700  missing values.
#> 2024-11-08 16:07:11.943369 INFO::Introduced 0 missing values due to singular proteins at batch/covariate level.
#> 2024-11-08 16:07:11.943543 INFO::Done
#> 2024-11-08 16:07:11.943667 INFO::Acquiring quality metrics before batch effect correction.
#> 2024-11-08 16:07:11.948305 INFO::Starting hierarchical adjustment
#> 2024-11-08 16:07:11.94851 INFO::Found  20  batches.
#> 2024-11-08 16:07:11.94864 INFO::Cores argument is not defined or BPPARAM has been specified. Argument corereduction will not be used.
#> 2024-11-08 16:07:11.948812 INFO::Using default BPPARAM
#> 2024-11-08 16:07:11.948937 INFO::Processing subtree level 1
#> 2024-11-08 16:07:12.104617 INFO::Processing subtree level 2
#> 2024-11-08 16:07:12.250178 INFO::Processing subtree level 3
#> 2024-11-08 16:07:12.415528 INFO::Adjusting the last 1 batches sequentially
#> 2024-11-08 16:07:12.41765 INFO::Done
#> 2024-11-08 16:07:12.418104 INFO::Acquiring quality metrics after batch effect correction.
#> 2024-11-08 16:07:12.422808 INFO::ASW Batch was 0.432699232338476 prior to batch effect correction and is now -0.139354730988 .
#> 2024-11-08 16:07:12.423016 INFO::ASW Label was 0.343823202049922 prior to batch effect correction and is now 0.819916170768385 .
#> 2024-11-08 16:07:12.423471 INFO::Total function execution time is  0.492261171340942  s and adjustment time is  0.469177961349487 s ( 95.31 )

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-11-08 16:07:12.451357 INFO::Formatting Data.
#> 2024-11-08 16:07:12.451599 INFO::Replacing NaNs with NAs.
#> 2024-11-08 16:07:12.451969 INFO::Removing potential empty rows and columns
#> 2024-11-08 16:07:12.452665 INFO::Found  2700  missing values.
#> 2024-11-08 16:07:12.460479 INFO::Introduced 0 missing values due to singular proteins at batch/covariate level.
#> 2024-11-08 16:07:12.460652 INFO::Done
#> 2024-11-08 16:07:12.460801 INFO::Acquiring quality metrics before batch effect correction.
#> 2024-11-08 16:07:12.465386 INFO::Starting hierarchical adjustment
#> 2024-11-08 16:07:12.46561 INFO::Found  20  batches.
#> 2024-11-08 16:07:12.713541 INFO::Set up parallel execution backend with 2 workers
#> 2024-11-08 16:07:12.713791 INFO::Processing subtree level 1 with 20 batches using 2 cores.
#> 2024-11-08 16:07:13.791776 INFO::Adjusting the last 2 batches sequentially
#> 2024-11-08 16:07:13.792269 INFO::Adjusting sequential tree level 1 with 2 batches
#> 2024-11-08 16:07:14.297923 INFO::Done
#> 2024-11-08 16:07:14.298127 INFO::Acquiring quality metrics after batch effect correction.
#> 2024-11-08 16:07:14.302012 INFO::ASW Batch was 0.505255168870644 prior to batch effect correction and is now -0.103500128143621 .
#> 2024-11-08 16:07:14.302174 INFO::ASW Label was 0.28794630407003 prior to batch effect correction and is now 0.814480536962776 .
#> 2024-11-08 16:07:14.302366 INFO::Total function execution time is  1.8510468006134  s and adjustment time is  1.83226299285889 s ( 98.99 )

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-11-08 16:07:14.331847 INFO::Formatting Data.
#> 2024-11-08 16:07:14.332113 INFO::Recognized SummarizedExperiment
#> 2024-11-08 16:07:14.332266 INFO::Typecasting input to dataframe.
#> 2024-11-08 16:07:14.342947 INFO::Replacing NaNs with NAs.
#> 2024-11-08 16:07:14.343296 INFO::Removing potential empty rows and columns
#> 2024-11-08 16:07:14.344346 INFO::Found  0  missing values.
#> 2024-11-08 16:07:14.346175 INFO::Introduced 0 missing values due to singular proteins at batch/covariate level.
#> 2024-11-08 16:07:14.346324 INFO::Done
#> 2024-11-08 16:07:14.346461 INFO::Acquiring quality metrics before batch effect correction.
#> 2024-11-08 16:07:14.347683 INFO::Starting hierarchical adjustment
#> 2024-11-08 16:07:14.347875 INFO::Found  2  batches.
#> 2024-11-08 16:07:14.348018 INFO::Cores argument is not defined or BPPARAM has been specified. Argument corereduction will not be used.
#> 2024-11-08 16:07:14.348181 INFO::Using default BPPARAM
#> 2024-11-08 16:07:14.348313 INFO::Adjusting the last 2 batches sequentially
#> 2024-11-08 16:07:14.348562 INFO::Adjusting sequential tree level 1 with 2 batches
#> 2024-11-08 16:07:14.366843 INFO::Done
#> 2024-11-08 16:07:14.367038 INFO::Acquiring quality metrics after batch effect correction.
#> 2024-11-08 16:07:14.368248 INFO::ASW Batch was 0.0182173584368604 prior to batch effect correction and is now -0.0858028241380684 .
#> 2024-11-08 16:07:14.368407 INFO::ASW Label was 0.00770208313591983 prior to batch effect correction and is now 0.0253785873145086 .
#> 2024-11-08 16:07:14.368615 INFO::Total function execution time is  0.0367729663848877  s and adjustment time is  0.0189998149871826 s ( 51.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)
#> 2024-11-08 16:07:14.393989 INFO::Formatting Data.
#> 2024-11-08 16:07:14.394233 INFO::Replacing NaNs with NAs.
#> 2024-11-08 16:07:14.39452 INFO::Removing potential empty rows and columns
#> 2024-11-08 16:07:14.395067 INFO::Found  1350  missing values.
#> 2024-11-08 16:07:14.395309 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-11-08 16:07:14.400332 INFO::Introduced 0 missing values due to singular proteins at batch/covariate level.
#> 2024-11-08 16:07:14.400493 INFO::Done
#> 2024-11-08 16:07:14.400631 INFO::Acquiring quality metrics before batch effect correction.
#> 2024-11-08 16:07:14.402375 INFO::Starting hierarchical adjustment
#> 2024-11-08 16:07:14.402565 INFO::Found  5  batches.
#> 2024-11-08 16:07:14.402703 INFO::Cores argument is not defined or BPPARAM has been specified. Argument corereduction will not be used.
#> 2024-11-08 16:07:14.402857 INFO::Using default BPPARAM
#> 2024-11-08 16:07:14.402989 INFO::Processing subtree level 1
#> 2024-11-08 16:07:14.482791 INFO::Adjusting the last 2 batches sequentially
#> 2024-11-08 16:07:14.484371 INFO::Adjusting sequential tree level 1 with 2 batches
#> 2024-11-08 16:07:14.503773 INFO::Done
#> 2024-11-08 16:07:14.503993 INFO::Acquiring quality metrics after batch effect correction.
#> 2024-11-08 16:07:14.506227 INFO::ASW Batch was 0.492773245691086 prior to batch effect correction and is now -0.0377157224767566 .
#> 2024-11-08 16:07:14.506429 INFO::ASW Label was 0.40854766060101 prior to batch effect correction and is now 0.895560693013661 .
#> 2024-11-08 16:07:14.506797 INFO::Total function execution time is  0.112818956375122  s and adjustment time is  0.101233959197998 s ( 89.73 )

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-11-08 16:07:14.536284 INFO::Formatting Data.
#> 2024-11-08 16:07:14.536629 INFO::Replacing NaNs with NAs.
#> 2024-11-08 16:07:14.537014 INFO::Removing potential empty rows and columns
#> 2024-11-08 16:07:14.537354 INFO::Found  60  missing values.
#> 2024-11-08 16:07:14.53855 INFO::Introduced 0 missing values due to singular proteins at batch/covariate level.
#> 2024-11-08 16:07:14.53871 INFO::Done
#> 2024-11-08 16:07:14.538863 INFO::Acquiring quality metrics before batch effect correction.
#> 2024-11-08 16:07:14.539808 INFO::Starting hierarchical adjustment
#> 2024-11-08 16:07:14.540013 INFO::Found  4  batches.
#> 2024-11-08 16:07:14.540165 INFO::Cores argument is not defined or BPPARAM has been specified. Argument corereduction will not be used.
#> 2024-11-08 16:07:14.540352 INFO::Using default BPPARAM
#> 2024-11-08 16:07:14.540509 INFO::Processing subtree level 1
#> 2024-11-08 16:07:14.580359 INFO::Adjusting the last 2 batches sequentially
#> 2024-11-08 16:07:14.581739 INFO::Adjusting sequential tree level 1 with 2 batches
#> 2024-11-08 16:07:14.590039 INFO::Done
#> 2024-11-08 16:07:14.590245 INFO::Acquiring quality metrics after batch effect correction.
#> 2024-11-08 16:07:14.591358 INFO::ASW Batch was 0.440355021914032 prior to batch effect correction and is now -0.0874802787366291 .
#> 2024-11-08 16:07:14.591545 INFO::ASW Label was 0.373906827748893 prior to batch effect correction and is now 0.919791677398366 .
#> 2024-11-08 16:07:14.591808 INFO::Total function execution time is  0.0555551052093506  s and adjustment time is  0.0500631332397461 s ( 90.11 )

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: aarch64-apple-darwin20
#> Running under: macOS Ventura 13.6.7
#> 
#> Matrix products: default
#> BLAS:   /Library/Frameworks/R.framework/Versions/4.4-arm64/Resources/lib/libRblas.0.dylib 
#> LAPACK: /Library/Frameworks/R.framework/Versions/4.4-arm64/Resources/lib/libRlapack.dylib;  LAPACK version 3.12.0
#> 
#> locale:
#> [1] C/en_US.UTF-8/en_US.UTF-8/C/en_US.UTF-8/en_US.UTF-8
#> 
#> time zone: America/New_York
#> tzcode source: internal
#> 
#> attached base packages:
#> [1] stats     graphics  grDevices utils     datasets  methods   base     
#> 
#> other attached packages:
#> [1] BERT_1.2.0       BiocStyle_2.34.0
#> 
#> loaded via a namespace (and not attached):
#>  [1] tidyselect_1.2.1            blob_1.2.4                 
#>  [3] Biostrings_2.74.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.46.0             invgamma_1.1               
#> [15] RSQLite_2.3.7               magrittr_2.0.3             
#> [17] genefilter_1.88.0           compiler_4.4.1             
#> [19] rlang_1.1.4                 sass_0.4.9                 
#> [21] tools_4.4.1                 yaml_2.3.10                
#> [23] knitr_1.49                  S4Arrays_1.6.0             
#> [25] bit_4.5.0                   DelayedArray_0.32.0        
#> [27] abind_1.4-8                 BiocParallel_1.40.0        
#> [29] BiocGenerics_0.52.0         grid_4.4.1                 
#> [31] stats4_4.4.1                xtable_1.8-4               
#> [33] edgeR_4.4.0                 iterators_1.0.14           
#> [35] logging_0.10-108            SummarizedExperiment_1.36.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] zlibbioc_1.52.0             splines_4.4.1              
#> [47] parallel_4.4.1              AnnotationDbi_1.68.0       
#> [49] BiocManager_1.30.25         XVector_0.46.0             
#> [51] matrixStats_1.4.1           vctrs_0.6.5                
#> [53] Matrix_1.7-1                jsonlite_1.8.9             
#> [55] sva_3.54.0                  bookdown_0.41              
#> [57] comprehenr_0.6.10           IRanges_2.40.0             
#> [59] S4Vectors_0.44.0            bit64_4.5.2                
#> [61] locfit_1.5-9.10             foreach_1.5.2              
#> [63] limma_3.62.1                jquerylib_0.1.4            
#> [65] annotate_1.84.0             glue_1.8.0                 
#> [67] codetools_0.2-20            lubridate_1.9.3            
#> [69] stringi_1.8.4               GenomeInfoDb_1.42.0        
#> [71] GenomicRanges_1.58.0        UCSC.utils_1.2.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.66.0             
#> [79] png_0.1-8                   memoise_2.0.1              
#> [81] snakecase_0.11.1            bslib_0.8.0                
#> [83] SparseArray_1.6.0           nlme_3.1-166               
#> [85] mgcv_1.9-1                  xfun_0.49                  
#> [87] MatrixGenerics_1.18.0