BiocStyle::markdown()
knitr::opts_chunk$set(fig.wide = TRUE, fig.retina = 3, error=FALSE, eval=TRUE)
knitr::include_graphics("octopussy.png")
The octopussy `rawDiag` package logo by Lilly van de Venn.

Figure 1: The octopussy rawDiag package logo by Lilly van de Venn

1 Introduction

Over the past two decades, liquid chromatography coupled to mass spectrometry (LC–MS) has evolved into the method of choice in the field of proteomics. (Cox and Mann 2011; Mallick and Kuster 2010) During a typical LC–MS measurement, a complex mixture of analytes is separated by a liquid chromatography system coupled to a mass spectrometer (MS) through an ion source interface. This interface converts the analytes that elute from the chromatography system over time into a beam of ions. The MS records from this ion beam a series of mass spectra containing detailed information on the analyzed sample. (Savaryn, Toby, and Kelleher 2016) The resulting raw data consist of the mass spectra and their metadata, typically recorded in a vendor-specific binary format. During a measurement the mass spectrometer applies internal heuristics, which enables the instrument to adapt to sample properties, for example, sample complexity and amount of ions in near real time. Still, method parameters controlling these heuristics need to be set prior to the measurement. Optimal measurement results require a careful balancing of instrument parameters, but their complex interactions with each other make LC–MS method optimization a challenging task.

Here we present rawDiag, a platform-independent software tool implemented in the R language (Becker, Chambers, and Wilks 1988) that supports LC–MS operators during the process of empirical method optimization. Our work builds on the ideas of the discontinued software rawMeat (VAST Scientific). Our application is currently tailored toward spectral data acquired on Thermo Fisher Scientific instruments (raw format), with a particular focus on Orbitrap (Zubarev and Makarov 2013) mass analyzers (Exactive or Fusion instruments). These instruments are heavily used in the field of bottom-up proteomics (Aebersold and Mann 2003) to analyze complex peptide mixtures derived from enzymatic digests of proteomes.

rawDiag is meant to run after MS acquisition, optimally as an interactive R shiny application, and produces a series of diagnostic plots visualizing the impact of method parameter choices on the acquired data across injections. If static reports are required then pdf files can be generated using rmarkdown. In this vignette, we present the usage of our tool.

rawDiag gains advantages from being part of the Bioconductor ecosystem, such as its ability to utilize the rawrr package and potentially extend its functionality through interaction with the Spectra infrastructure, particularly with the MsBackendRawFileReader.

2 Requirements

rawDiag proides a wrapper function readRaw using the rawrr methods raw::readIndex, rawrr::readTrailer, and rawrr::readChromatogram to read proprietary mass spectrometer generated data by invoking third-party managed methods through a system2 text connection. The rawrr package provides the entire stack below, which rawDiag utilizes.

R>
text connection
system2
Mono Runtime
Managed Assembly (CIL/.NET code)
rawrr.exe
ThermoFisher.CommonCore.*.dll

3 Install

3.1 OS dependencies

3.1.1 Linux (debian/ubuntu)

In case you prefer to compile rawrr.exe from C# source code, please install the mono compiler and xbuild by installing the following Linux packages:

sudo apt-get install mono-mcs mono-xbuild

Otherwise, to execute the precompiled code, the following Linux packages are sufficient:

sudo apt-get install mono-runtime libmono-system-data4.0-cil -y

3.1.2 macOS (Catalina/BigSur)

brew install mono

or install from

https://www.mono-project.com/

3.1.3 Microsoft Windows

Running the rawrr.exe will run out of the box.

If the native C# compiler is not available install mono from:

https://www.mono-project.com/

3.2 R

To install this package, start R (version “>=4.4”) and enter:

if (!require("BiocManager", quietly = TRUE))
    install.packages("BiocManager")

BiocManager::install("rawrr")

assemblies aka Common Intermediate Language bytecode - the download and install can be done on all platforms using the command:

rawDiag::checkRawrr
## function () 
## {
##     if (isFALSE(requireNamespace("BiocManager", quietly = TRUE))) 
##         stop("exec", "install.packages('BiocManager')")
##     if (isFALSE(requireNamespace("rawrr", quietly = TRUE))) 
##         stop("exec", "BiocManager::install('rawrr')")
##     if (isFALSE(rawrr:::.checkRawFileReaderDLLs())) 
##         rawrr::installRawFileReaderDLLs()
##     if (isFALSE(file.exists(rawrr:::.rawrrAssembly()))) 
##         rawrr::installRawrrExe()
##     if (isFALSE(rawrr:::.isAssemblyWorking())) 
##         stop("rawrr assembly is not working. check if mono if available.")
##     TRUE
## }
## <bytecode: 0x55a3831e7ca8>
## <environment: namespace:rawDiag>
rawDiag::checkRawrr()
## [1] TRUE
if (isFALSE(rawrr::.checkDllInMonoPath())){
  rawrr::installRawFileReaderDLLs()
}
## removing files in directory '/home/biocbuild/.cache/R/rawrr/rawrrassembly'
##  ThermoFisher.CommonCore.BackgroundSubtraction.dll 
##                                                  0 
##                   ThermoFisher.CommonCore.Data.dll 
##                                                  0 
## ThermoFisher.CommonCore.MassPrecisionEstimator.dll 
##                                                  0 
##          ThermoFisher.CommonCore.RawFileReader.dll 
##                                                  0
rawrr::installRawrrExe()
## MD5 96e3a4cc1b7caaf92890d85ed4c72f77 /home/biocbuild/.cache/R/rawrr/rawrrassembly/rawrr.exe
## [1] 0

for more information please read the INSTALL file in the rawrr package.

4 Usage

4.1 R command line

4.1.1 Input

fetch example Orbitrap raw files from ExperimentHub’s tartare package.

library(ExperimentHub)
## Loading required package: BiocGenerics
## 
## Attaching package: 'BiocGenerics'
## The following objects are masked from 'package:stats':
## 
##     IQR, mad, sd, var, xtabs
## The following objects are masked from 'package:base':
## 
##     Filter, Find, Map, Position, Reduce, anyDuplicated, aperm, append,
##     as.data.frame, basename, cbind, colnames, dirname, do.call,
##     duplicated, eval, evalq, get, grep, grepl, intersect, is.unsorted,
##     lapply, mapply, match, mget, order, paste, pmax, pmax.int, pmin,
##     pmin.int, rank, rbind, rownames, sapply, setdiff, table, tapply,
##     union, unique, unsplit, which.max, which.min
## Loading required package: AnnotationHub
## Loading required package: BiocFileCache
## Loading required package: dbplyr
ExperimentHub::ExperimentHub() -> eh
normalizePath(eh[["EH3222"]]) -> EH3222
## see ?tartare and browseVignettes('tartare') for documentation
## loading from cache
normalizePath(eh[["EH4547"]]) -> EH4547
## see ?tartare and browseVignettes('tartare') for documentation
## loading from cache
(rawfileEH3222 <- paste0(EH3222, ".raw"))
## [1] "/home/biocbuild/.cache/R/ExperimentHub/21fa175b5a0962_3238.raw"
if (!file.exists(rawfileEH3222)){
  file.copy(EH3222, rawfileEH3222)
}

(rawfileEH4547 <- paste0(EH4547, ".raw"))
## [1] "/home/biocbuild/.cache/R/ExperimentHub/21fa1715d0f23a_4590.raw"
if (!file.exists(rawfileEH4547)){
  file.copy(EH4547, rawfileEH4547)
}

c(rawfileEH3222, rawfileEH4547) -> rawfile

Of note, the proprietary .Net assemblies (Shofstahl 2018) require a file extentention of .raw. Therfore we have to rename the EH files and add the .raw suffix.

list meta data of the raw files.

(rawfile |>
  lapply(FUN = rawrr::readFileHeader) -> rawFileHeader)
## [[1]]
## [[1]]$`RAW file`
## [1] "21fa175b5a0962_3238.raw"
## 
## [[1]]$`RAW file version`
## [1] "66"
## 
## [[1]]$`Creation date`
## [1] "7/16/2019 5:56:24 PM"
## 
## [[1]]$Operator
## [1] "lumos"
## 
## [[1]]$`Number of instruments`
## [1] 2
## 
## [[1]]$Description
## [1] ""
## 
## [[1]]$`Instrument model`
## [1] "Orbitrap Fusion Lumos"
## 
## [[1]]$`Instrument name`
## [1] "Orbitrap Fusion Lumos"
## 
## [[1]]$`Instrument method`
## [1] "C:/Xcalibur/methods/p3181/methMelt_GT20_RT35.meth"
## 
## [[1]]$`Serial number`
## [1] "FSN20583"
## 
## [[1]]$`Software version`
## [1] "3.1.2412.25"
## 
## [[1]]$`Firmware version`
## [1] ""
## 
## [[1]]$Units
## [1] "None"
## 
## [[1]]$`Mass resolution`
## [1] "0.500"
## 
## [[1]]$`Number of scans`
## [1] 8742
## 
## [[1]]$`Number of ms2 scans`
## [1] 2346
## 
## [[1]]$`Scan range`
## [1]    1 8742
## 
## [[1]]$`Time range`
## [1]  0.00 35.01
## 
## [[1]]$`Mass range`
## [1]   91 2000
## 
## [[1]]$`Scan filter (first scan)`
## [1] "FTMS + c ESI Full ms [350.0000-2000.0000]"
## 
## [[1]]$`Scan filter (last scan)`
## [1] "FTMS + c ESI Full ms [350.0000-2000.0000]"
## 
## [[1]]$`Total number of filters`
## [1] "2205"
## 
## [[1]]$`Sample name`
## [1] ""
## 
## [[1]]$`Sample id`
## [1] "1:A,5"
## 
## [[1]]$`Sample type`
## [1] "Unknown"
## 
## [[1]]$`Sample comment`
## [1] ""
## 
## [[1]]$`Sample vial`
## [1] "1:B,1"
## 
## [[1]]$`Sample volume`
## [1] "0"
## 
## [[1]]$`Sample injection volume`
## [1] "2"
## 
## [[1]]$`Sample row number`
## [1] "0"
## 
## [[1]]$`Sample dilution factor`
## [1] "1"
## 
## [[1]]$`Sample barcode`
## [1] ""
## 
## [[1]]$`User text 0`
## [1] ""
## 
## [[1]]$`User text 1`
## [1] ""
## 
## [[1]]$`User text 2`
## [1] ""
## 
## [[1]]$`User text 3`
## [1] ""
## 
## [[1]]$`User text 4`
## [1] ""
## 
## 
## [[2]]
## [[2]]$`RAW file`
## [1] "21fa1715d0f23a_4590.raw"
## 
## [[2]]$`RAW file version`
## [1] "66"
## 
## [[2]]$`Creation date`
## [1] "11/16/2018 9:58:53 AM"
## 
## [[2]]$Operator
## [1] "QexactiveHF"
## 
## [[2]]$`Number of instruments`
## [1] 2
## 
## [[2]]$Description
## [1] ""
## 
## [[2]]$`Instrument model`
## [1] "Q Exactive HF Orbitrap"
## 
## [[2]]$`Instrument name`
## [1] "Q Exactive HF Orbitrap"
## 
## [[2]]$`Instrument method`
## [1] "C:/Xcalibur/methods/__QCloud/current_method/forward trap elute/autoQC01_TRAP_GT20min_RT35min.meth"
## 
## [[2]]$`Serial number`
## [1] "Exactive Series slot #2496"
## 
## [[2]]$`Software version`
## [1] "2.9-290204/2.9.2.2947"
## 
## [[2]]$`Firmware version`
## [1] "rev. 1"
## 
## [[2]]$Units
## [1] "None"
## 
## [[2]]$`Mass resolution`
## [1] "0.500"
## 
## [[2]]$`Number of scans`
## [1] 21881
## 
## [[2]]$`Number of ms2 scans`
## [1] 20885
## 
## [[2]]$`Scan range`
## [1]     1 21881
## 
## [[2]]$`Time range`
## [1]  0 35
## 
## [[2]]$`Mass range`
## [1]  100 1805
## 
## [[2]]$`Scan filter (first scan)`
## [1] "FTMS + c NSI Full ms [350.0000-1800.0000]"
## 
## [[2]]$`Scan filter (last scan)`
## [1] "FTMS + c NSI Full ms2 582.3190@hcd27.00 [100.0000-1205.0000]"
## 
## [[2]]$`Total number of filters`
## [1] "22"
## 
## [[2]]$`Sample name`
## [1] "autoQC01"
## 
## [[2]]$`Sample id`
## [1] "NA"
## 
## [[2]]$`Sample type`
## [1] "Unknown"
## 
## [[2]]$`Sample comment`
## [1] ""
## 
## [[2]]$`Sample vial`
## [1] "1:F,8"
## 
## [[2]]$`Sample volume`
## [1] "0"
## 
## [[2]]$`Sample injection volume`
## [1] "2"
## 
## [[2]]$`Sample row number`
## [1] "0"
## 
## [[2]]$`Sample dilution factor`
## [1] "0"
## 
## [[2]]$`Sample barcode`
## [1] ""
## 
## [[2]]$`User text 0`
## [1] "1000"
## 
## [[2]]$`User text 1`
## [1] ""
## 
## [[2]]$`User text 2`
## [1] "FGCZ"
## 
## [[2]]$`User text 3`
## [1] ""
## 
## [[2]]$`User text 4`
## [1] ""

4.1.2 readRaw - read Orbitrap raw file

read the two instrument raw files by using the rawDiag package.

rawfile |>
  BiocParallel::bplapply(FUN = rawDiag::readRaw) |>
  Reduce(f = rbind) -> x

4.1.3 Output - Visualization

This package provides several plot functions tailored toward MS data. The following list shows all available plot methods.

library(rawDiag)
ls("package:rawDiag") |>
  grep(pattern = '^plot', value = TRUE) -> pm

pm |>
  knitr::kable(col.names = "package:rawDiag plot functions")
package:rawDiag plot functions
plotChargeState
plotCycleLoad
plotCycleTime
plotInjectionTime
plotLockMassCorrection
plotMassDistribution
plotMzDistribution
plotPrecursorHeatmap
plotScanTime
plotTicBasepeak

An inherent problem of visualizing data is the fact that depending on the data at hand, specific visualizations lose their usefulness, e.g., overplotting in a scatter plot if too many data points are present. To address this problem, we implemented most of the plot functions in different versions inspired by the work of Cleveland (1993), Sarkar (2008), and Wickham (2009). The data can be displayed in trellis plot manner using the faceting functionality of ggplot2. Alternatively, overplotting using color coding or violin plots based on descriptive statistics values can be chosen, which allows the user to interactively change the appearance of the plots based on the situation at hand. For instance, a large number of files are best visualized by violin plots, giving the user an idea about the distribution of the data points. On the basis of this, a smaller subset of files can be selected and visualized with another technique.

The code snippet below applies all plot methods on the example data.

pm |> 
  lapply(FUN = function(plotFUN) {
    lapply(c('trellis'), function(method) {
      message("plotting", plotFUN, "using method", method, "...")
      do.call(plotFUN, list(x, method)) 
    })
  }) 
## plottingplotChargeStateusing methodtrellis...
## plottingplotCycleLoadusing methodtrellis...
## plottingplotCycleTimeusing methodtrellis...
## plottingplotInjectionTimeusing methodtrellis...
## plottingplotLockMassCorrectionusing methodtrellis...
## plottingplotMassDistributionusing methodtrellis...
## plottingplotMzDistributionusing methodtrellis...
## plottingplotPrecursorHeatmapusing methodtrellis...
## plottingplotScanTimeusing methodtrellis...
## plottingplotTicBasepeakusing methodtrellis...
## [[1]]
## [[1]][[1]]

## 
## 
## [[2]]
## [[2]][[1]]
## `geom_smooth()` using formula = 'y ~ x'

## 
## 
## [[3]]
## [[3]][[1]]

## 
## 
## [[4]]
## [[4]][[1]]

## 
## 
## [[5]]
## [[5]][[1]]
## Warning: Removed 6164 rows containing non-finite outside the scale range
## (`stat_smooth()`).
## Warning: Removed 6164 rows containing missing values or values outside the scale range
## (`geom_line()`).

## 
## 
## [[6]]
## [[6]][[1]]

## 
## 
## [[7]]
## [[7]][[1]]

## 
## 
## [[8]]
## [[8]][[1]]

## 
## 
## [[9]]
## [[9]][[1]]

## 
## 
## [[10]]
## [[10]][[1]]

The appearance of each plot depends on the instrument, sample, and method used to acquire the data. Therefore, it is hard to say what each ideal plot should look like. In particular, in the example above, we use data generated on an Orbitrap Fusion Lumos, 21fa175b5a0962_3238.raw and Q Exactive HF Orbitrap, 21fa1715d0f23a_4590.raw instrument using data-independent acquisition (DIA) (Bruderer et al. 2017) and data-dependent acquisition (DDA) methods. For more information on the plot methods and its application, please read the package man pages and the application examples in the manuscript (Trachsel et al. 2018).

4.2 Launching the shiny application

The package provides a simple interactive shiny-based graphical user interface for exploring Thermo Fisher Scientific raw data.

If you have a directory containing raw files, you can create a shiny application as follows:

rawfile |>
  dirname() |>
  rawDiag::buildRawDiagShinyApp() -> app

The shiny runApp function launches the app in our browser.

shiny::runApp(app)

By default, the application lets you choose the raw files in the provided directory and provides the visualizations of the raw data as output. The user can interactively change the by the rawDiag the package provided plot functions and arguments.

Additionally, the application provides PDF generation and download buttons. Optionally height and width can be changed in the user interface.

Of note, the rawDiag::rawDiagServer module can be integrated into an existing shinydashboard application, e.g., https://shiny-ms.fgcz.uzh.ch/fgczmsqc-dashboard/.

5 FAQ

5.1 I would like to load multiple files into a single data.frame to do comparisons; what is the preferred method for doing so?

consider all raw files of your working directory, e.g., ~/Downloads and load them.

file.path(Sys.getenv("HOME"), "Downloads") |>
  setwd()

list.files() |>
  grep(pattern = '*.raw$', value = TRUE) |> 
  lapply(FUN = rawDiag::readRaw) |> 
  Reduce(f = rbind) -> x

as alternative to lapply you can utilize the BiocParallel package bplapply function.

References

Aebersold, Ruedi, and Matthias Mann. 2003. “Mass Spectrometry-Based Proteomics.” Nature 422 (6928): 198–207. https://doi.org/10.1038/nature01511.

Becker, Richard A., John M. Chambers, and Allan R. Wilks. 1988. The New S Language. London: Chapman & Hall.

Bruderer, Roland, Oliver M. Bernhardt, Tejas Gandhi, Yue Xuan, Julia Sondermann, Manuela Schmidt, David Gomez-Varela, and Lukas Reiter. 2017. “Optimization of Experimental Parameters in Data-Independent Mass Spectrometry Significantly Increases Depth and Reproducibility of Results.” Molecular &Amp; Cellular Proteomics 16 (12): 2296–2309. https://doi.org/10.1074/mcp.ra117.000314.

Cleveland, William S. 1993. Visualizing Data. 1st ed. Hobart Press, Summit, New Jersey, U.S.A.

Cox, Jürgen, and Matthias Mann. 2011. “Quantitative, High-Resolution Proteomics for Data-Driven Systems Biology.” Annual Review of Biochemistry 80 (1): 273–99. https://doi.org/10.1146/annurev-biochem-061308-093216.

Kockmann, Tobias, and Christian Panse. 2021. “The rawrr R Package: Direct Access to Orbitrap Data and Beyond.” Journal of Proteome Research. https://doi.org/10.1021/acs.jproteome.0c00866.

Mallick, Parag, and Bernhard Kuster. 2010. “Proteomics: A Pragmatic Perspective.” Nature Biotechnology 28 (7): 695–709. https://doi.org/10.1038/nbt.1658.

Sarkar, Deepayan. 2008. Lattice: Multivariate Data Visualization with R. Springer New York. https://doi.org/10.1007/978-0-387-75969-2.

Savaryn, John P., Timothy K. Toby, and Neil L. Kelleher. 2016. “A Researcher’s Guide to Mass Spectrometry‐based Proteomics.” PROTEOMICS 16 (18): 2435–43. https://doi.org/10.1002/pmic.201600113.

Shofstahl, Jim. 2018. “New RawFileReader from Thermo Fisher Scientific, Http://Planetorbitrap.com/Rawfilereader.” 2018. http://planetorbitrap.com/rawfilereader#.WvWESK3QPmE.

Trachsel, Christian, Christian Panse, Tobias Kockmann, Witold E. Wolski, Jonas Grossmann, and Ralph Schlapbach. 2018. “rawDiag: An R Package Supporting Rational LCMS Method Optimization for Bottom-up Proteomics.” Journal of Proteome Research 17 (8): 2908–14. https://doi.org/10.1021/acs.jproteome.8b00173.

Wickham, Hadley. 2009. Ggplot2: Elegant Graphics for Data Analysis. Springer-Verlag New York. https://doi.org/10.1007/978-0-387-98141-3.

Zubarev, Roman A., and Alexander Makarov. 2013. “Orbitrap Mass Spectrometry.” Analytical Chemistry 85 (11): 5288–96. https://doi.org/10.1021/ac4001223.

Appendix

Session information

sessionInfo()
## R version 4.4.0 (2024-04-24)
## Platform: x86_64-pc-linux-gnu
## Running under: Ubuntu 22.04.4 LTS
## 
## Matrix products: default
## BLAS:   /home/biocbuild/bbs-3.19-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] rawDiag_1.0.0        tartare_1.18.0       ExperimentHub_2.12.0
## [4] AnnotationHub_3.12.0 BiocFileCache_2.12.0 dbplyr_2.5.0        
## [7] BiocGenerics_0.50.0  BiocStyle_2.32.0    
## 
## loaded via a namespace (and not attached):
##  [1] tidyselect_1.2.1        dplyr_1.1.4             farver_2.1.2           
##  [4] blob_1.2.4              filelock_1.0.3          Biostrings_2.72.0      
##  [7] fastmap_1.1.1           promises_1.3.0          digest_0.6.35          
## [10] mime_0.12               lifecycle_1.0.4         KEGGREST_1.44.0        
## [13] RSQLite_2.3.6           magrittr_2.0.3          compiler_4.4.0         
## [16] rlang_1.1.3             sass_0.4.9              tools_4.4.0            
## [19] utf8_1.2.4              yaml_2.3.8              knitr_1.46             
## [22] labeling_0.4.3          bit_4.0.5               curl_5.2.1             
## [25] plyr_1.8.9              BiocParallel_1.38.0     withr_3.0.0            
## [28] purrr_1.0.2             grid_4.4.0              stats4_4.4.0           
## [31] fansi_1.0.6             xtable_1.8-4            colorspace_2.1-0       
## [34] ggplot2_3.5.1           scales_1.3.0            tinytex_0.51           
## [37] cli_3.6.2               rmarkdown_2.26          crayon_1.5.2           
## [40] generics_0.1.3          httr_1.4.7              reshape2_1.4.4         
## [43] rawrr_1.12.0            DBI_1.2.2               cachem_1.0.8           
## [46] stringr_1.5.1           zlibbioc_1.50.0         splines_4.4.0          
## [49] parallel_4.4.0          AnnotationDbi_1.66.0    BiocManager_1.30.23    
## [52] XVector_0.44.0          vctrs_0.6.5             Matrix_1.7-0           
## [55] jsonlite_1.8.8          bookdown_0.39           IRanges_2.38.0         
## [58] S4Vectors_0.42.0        bit64_4.0.5             magick_2.8.3           
## [61] fontawesome_0.5.2       jquerylib_0.1.4         hexbin_1.28.3          
## [64] glue_1.7.0              codetools_0.2-20        stringi_1.8.4          
## [67] gtable_0.3.5            BiocVersion_3.19.1      later_1.3.2            
## [70] GenomeInfoDb_1.40.0     UCSC.utils_1.0.0        munsell_0.5.1          
## [73] tibble_3.2.1            pillar_1.9.0            rappdirs_0.3.3         
## [76] htmltools_0.5.8.1       GenomeInfoDbData_1.2.12 R6_2.5.1               
## [79] evaluate_0.23           shiny_1.8.1.1           lattice_0.22-6         
## [82] Biobase_2.64.0          highr_0.10              png_0.1-8              
## [85] memoise_2.0.1           httpuv_1.6.15           bslib_0.7.0            
## [88] Rcpp_1.0.12             nlme_3.1-164            mgcv_1.9-1             
## [91] xfun_0.43               pkgconfig_2.0.3