This document describes how to use the R-package IPO
to
optimize xcms
parameters. Code examples on how to use
IPO
are provided. Additional to IPO
the
R-packages xcms
and rsm
are required. The
R-package msdata
andmtbls2
are recommended.
The optimization process looks as following:
# try http:// if https:// URLs are not supported
if (!requireNamespace("BiocManager", quietly=TRUE))
install.packages("BiocManager")
BiocManager::install("IPO")
Installing main suggested packages
xcms
handles the file processing hence all files can be
used that can be processed by xcms
.
To optimize parameters different values (levels) have to tested for
these parameters. To efficiently test many different levels design of
experiment (DoE) is used. Box-Behnken and central composite designs set
three evenly spaced levels for each parameter. The method
getDefaultXcmsSetStartingParams
provides default values for
the lower and upper levels defining a range. Since the levels are evenly
spaced the middle level or center point is calculated automatically. To
edit the starting levels of a parameter set the lower and upper level as
desired. If a parameter should not be optimized, set a single default
value for xcms
processing, do not set this parameter to
NULL.
The method getDefaultXcmsSetStartingParams
creates a
list with default values for the optimization of the peak picking
methods centWave
or matchedFilter
. To choose
between these two method set the parameter accordingly.
The method optimizeXcmsSet
has the following
parameters:
xcms
peak picking methods parameters. A default list is
created by getDefaultXcmsSetStartingParams()
.BiocParallelParam
-object (see
?BiocParallel::BiocParallelParam
) to controll the use of
parallelisation of xcms
. Defaults to
bpparam()
.NULL
if no rsm’s should be saved.The optimization process starts at the specified levels. After the calculation of the DoE is finished the result is evaluated and the levels automatically set accordingly. Then a new DoE is generated and processed. This continues until an optimum is found.
The result of peak picking optimization is a list consisting of all
calculated DoEs including the used levels, design, response, rsm and
best setting. Additionally the last list item is a list
(\$best_settings
) providing the optimized parameters
(\$parameters
), an xcmsSet object (\$xset
)
calculated with these parameters and the response this
xcms
-object gives.
peakpickingParameters <- getDefaultXcmsSetStartingParams('matchedFilter')
#setting levels for step to 0.2 and 0.3 (hence 0.25 is the center point)
peakpickingParameters$step <- c(0.2, 0.3)
peakpickingParameters$fwhm <- c(40, 50)
#setting only one value for steps therefore this parameter is not optimized
peakpickingParameters$steps <- 2
time.xcmsSet <- system.time({ # measuring time
resultPeakpicking <-
optimizeXcmsSet(files = datafiles[1:2],
params = peakpickingParameters,
nSlaves = 1,
subdir = NULL,
plot = TRUE)
})
#>
#> starting new DoE with:
#> fwhm: c(40, 50)
#> snthresh: c(3, 17)
#> step: c(0.2, 0.3)
#> steps: 2
#> sigma: 0
#> max: 5
#> mzdiff: 0
#> index: FALSE
#> 1
#> 2
#> 3
#> 4
#> 5
#> 6
#> 7
#> 8
#> 9
#> 10
#> 11
#> 12
#> 13
#> 14
#> 15
#> 16
#>
#> starting new DoE with:
#> fwhm: c(45, 55)
#> snthresh: c(1, 15)
#> step: c(0.22, 0.3)
#> steps: 2
#> sigma: 0
#> max: 5
#> mzdiff: 0
#> index: FALSE
#> 1
#> 2
#> 3
#> 4
#> 5
#> 6
#> 7
#> 8
#> 9
#> 10
#> 11
#> 12
#> 13
#> 14
#> 15
#> 16
#>
#> starting new DoE with:
#> fwhm: c(50, 60)
#> snthresh: c(1, 15)
#> step: c(0.26, 0.34)
#> steps: 2
#> sigma: 0
#> max: 5
#> mzdiff: 0
#> index: FALSE
#> 1
#> 2
#> 3
#> 4
#> 5
#> 6
#> 7
#> 8
#> 9
#> 10
#> 11
#> 12
#> 13
#> 14
#> 15
#> 16
#> no increase, stopping
#> best parameter settings:
#> fwhm: 50
#> snthresh: 3
#> step: 0.26
#> steps: 2
#> sigma: 21.2332257516562
#> max: 5
#> mzdiff: 0.28
#> index: FALSE
resultPeakpicking$best_settings$result
#> ExpId #peaks #NonRP #RP PPS
#> 0.000 3228.000 2264.000 569.000 143.004
optimizedXcmsSetObject <- resultPeakpicking$best_settings$xset
The response surface models of all optimization steps for the parameter optimization of peak picking are shown above.
Currently the xcms
peak picking methods
centWave
and matchedFilter
are supported. The
parameter peakwidth
of the peak picking method
centWave
needs two values defining a minimum and maximum
peakwidth. These two values need separate optimization and are therefore
split into min_peakwidth
and max_peakwidth
in
getDefaultXcmsSetStartingParams
. Also for the
centWave
parameter prefilter two values have to be set. To
optimize these use set prefilter
to optimize the first
value and prefilter_value
to optimize the second value
respectively.
Optimization of retention time correction and grouping parameters is
done simultaneously. The method
getDefaultRetGroupStartingParams
provides default
optimization levels for the xcms
retention time correction
method obiwarp
and the grouping method
density
. Modifying these levels should be done the same way
done for the peak picking parameter optimization.
The method getDefaultRetGroupStartingParams
only
supports one retention time correction method (obiwarp
) and
one grouping method (density
) at the moment.
The method optimizeRetGroup
provides the following
parameter: - xset: an xcmsSet
-object used as basis for
retention time correction and grouping. - params: a list consisting of
items named according to xcms
retention time correction and
grouping methods parameters. A default list is created by
getDefaultRetGroupStartingParams
. - nSlaves: the number of
experiments of an DoE processed in parallel - subdir: a directory where
the response surface models are stored. Can also be NULL if no rsm’s
should be saved.
A list is returned similar to the one returned from peak picking
optimization. The last list item consists of the optimized retention
time correction and grouping parameters
(\$best_settings
).
retcorGroupParameters <- getDefaultRetGroupStartingParams()
retcorGroupParameters$profStep <- 1
retcorGroupParameters$gapExtend <- 2.7
time.RetGroup <- system.time({ # measuring time
resultRetcorGroup <-
optimizeRetGroup(xset = optimizedXcmsSetObject,
params = retcorGroupParameters,
nSlaves = 1,
subdir = NULL,
plot = TRUE)
})
#>
#> starting new DoE with:
#> distFunc: cor_opt
#> gapInit: c(0, 0.4)
#> gapExtend: 2.7
#> profStep: 1
#> plottype: none
#> response: 1
#> factorDiag: 2
#> factorGap: 1
#> localAlignment: 0
#> retcorMethod: obiwarp
#> bw: c(22, 38)
#> minfrac: c(0.3, 0.7)
#> mzwid: c(0.015, 0.035)
#> minsamp: 1
#> max: 50
#> center: 2
#> center sample: ko16
#> Processing: ko15
#> Create profile matrix with method 'bin' and step 1 ...
#> OK
#> Create profile matrix with method 'bin' and step 1 ... OK
#>
#> center sample: ko16
#> Processing: ko15
#> Create profile matrix with method 'bin' and step 1 ... OK
#> Create profile matrix with method 'bin' and step 1 ... OK
#>
#> center sample: ko16
#> Processing: ko15
#> Create profile matrix with method 'bin' and step 1 ... OK
#> Create profile matrix with method 'bin' and step 1 ... OK
#>
#> center sample: ko16
#> Processing: ko15
#> Create profile matrix with method 'bin' and step 1 ... OK
#> Create profile matrix with method 'bin' and step 1 ... OK
#>
#> center sample: ko16
#> Processing: ko15
#> Create profile matrix with method 'bin' and step 1 ... OK
#> Create profile matrix with method 'bin' and step 1 ... OK
#>
#> center sample: ko16
#> Processing: ko15
#> Create profile matrix with method 'bin' and step 1 ... OK
#> Create profile matrix with method 'bin' and step 1 ... OK
#>
#> center sample: ko16
#> Processing: ko15
#> Create profile matrix with method 'bin' and step 1 ... OK
#> Create profile matrix with method 'bin' and step 1 ... OK
#>
#> center sample: ko16
#> Processing: ko15
#> Create profile matrix with method 'bin' and step 1 ... OK
#> Create profile matrix with method 'bin' and step 1 ... OK
#>
#> center sample: ko16
#> Processing: ko15
#> Create profile matrix with method 'bin' and step 1 ... OK
#> Create profile matrix with method 'bin' and step 1 ... OK
#>
#> center sample: ko16
#> Processing: ko15
#> Create profile matrix with method 'bin' and step 1 ... OK
#> Create profile matrix with method 'bin' and step 1 ... OK
#>
#> center sample: ko16
#> Processing: ko15
#> Create profile matrix with method 'bin' and step 1 ... OK
#> Create profile matrix with method 'bin' and step 1 ... OK
#>
#> center sample: ko16
#> Processing: ko15
#> Create profile matrix with method 'bin' and step 1 ... OK
#> Create profile matrix with method 'bin' and step 1 ... OK
#>
#> center sample: ko16
#> Processing: ko15
#> Create profile matrix with method 'bin' and step 1 ... OK
#> Create profile matrix with method 'bin' and step 1 ... OK
#>
#> center sample: ko16
#> Processing: ko15
#> Create profile matrix with method 'bin' and step 1 ... OK
#> Create profile matrix with method 'bin' and step 1 ... OK
#>
#> center sample: ko16
#> Processing: ko15
#> Create profile matrix with method 'bin' and step 1 ... OK
#> Create profile matrix with method 'bin' and step 1 ... OK
#>
#> center sample: ko16
#> Processing: ko15
#> Create profile matrix with method 'bin' and step 1 ... OK
#> Create profile matrix with method 'bin' and step 1 ... OK
#>
#> center sample: ko16
#> Processing: ko15
#> Create profile matrix with method 'bin' and step 1 ... OK
#> Create profile matrix with method 'bin' and step 1 ... OK
#>
#> center sample: ko16
#> Processing: ko15
#> Create profile matrix with method 'bin' and step 1 ... OK
#> Create profile matrix with method 'bin' and step 1 ... OK
#>
#> center sample: ko16
#> Processing: ko15
#> Create profile matrix with method 'bin' and step 1 ... OK
#> Create profile matrix with method 'bin' and step 1 ... OK
#>
#> center sample: ko16
#> Processing: ko15
#> Create profile matrix with method 'bin' and step 1 ... OK
#> Create profile matrix with method 'bin' and step 1 ... OK
#>
#> center sample: ko16
#> Processing: ko15
#> Create profile matrix with method 'bin' and step 1 ... OK
#> Create profile matrix with method 'bin' and step 1 ... OK
#>
#> center sample: ko16
#> Processing: ko15
#> Create profile matrix with method 'bin' and step 1 ... OK
#> Create profile matrix with method 'bin' and step 1 ... OK
#>
#> center sample: ko16
#> Processing: ko15
#> Create profile matrix with method 'bin' and step 1 ... OK
#> Create profile matrix with method 'bin' and step 1 ... OK
#>
#> center sample: ko16
#> Processing: ko15
#> Create profile matrix with method 'bin' and step 1 ... OK
#> Create profile matrix with method 'bin' and step 1 ... OK
#>
#> center sample: ko16
#> Processing: ko15
#> Create profile matrix with method 'bin' and step 1 ... OK
#> Create profile matrix with method 'bin' and step 1 ... OK
#>
#> center sample: ko16
#> Processing: ko15
#> Create profile matrix with method 'bin' and step 1 ... OK
#> Create profile matrix with method 'bin' and step 1 ... OK
#> center sample: ko16
#> Processing: ko15
#> Create profile matrix with method 'bin' and step 1 ... OK
#> Create profile matrix with method 'bin' and step 1 ... OK
#>
#>
#>
#> starting new DoE with:
#>
#> gapInit: c(0.16, 0.64)
#> bw: c(12.4, 31.6)
#> minfrac: c(0.46, 0.94)
#> mzwid: c(0.023, 0.047)
#> distFunc: cor_opt
#> gapExtend: 2.7
#> profStep: 1
#> plottype: none
#> response: 1
#> factorDiag: 2
#> factorGap: 1
#> localAlignment: 0
#> retcorMethod: obiwarp
#> minsamp: 1
#> max: 50
#> center: 2
#> center sample: ko16
#> Processing: ko15
#> Create profile matrix with method 'bin' and step 1 ... OK
#> Create profile matrix with method 'bin' and step 1 ... OK
#>
#> center sample: ko16
#> Processing: ko15
#> Create profile matrix with method 'bin' and step 1 ... OK
#> Create profile matrix with method 'bin' and step 1 ... OK
#>
#> center sample: ko16
#> Processing: ko15
#> Create profile matrix with method 'bin' and step 1 ... OK
#> Create profile matrix with method 'bin' and step 1 ... OK
#>
#> center sample: ko16
#> Processing: ko15
#> Create profile matrix with method 'bin' and step 1 ... OK
#> Create profile matrix with method 'bin' and step 1 ... OK
#>
#> center sample: ko16
#> Processing: ko15
#> Create profile matrix with method 'bin' and step 1 ... OK
#> Create profile matrix with method 'bin' and step 1 ... OK
#>
#> center sample: ko16
#> Processing: ko15
#> Create profile matrix with method 'bin' and step 1 ... OK
#> Create profile matrix with method 'bin' and step 1 ... OK
#>
#> center sample: ko16
#> Processing: ko15
#> Create profile matrix with method 'bin' and step 1 ... OK
#> Create profile matrix with method 'bin' and step 1 ... OK
#>
#> center sample: ko16
#> Processing: ko15
#> Create profile matrix with method 'bin' and step 1 ... OK
#> Create profile matrix with method 'bin' and step 1 ... OK
#>
#> center sample: ko16
#> Processing: ko15
#> Create profile matrix with method 'bin' and step 1 ... OK
#> Create profile matrix with method 'bin' and step 1 ... OK
#>
#> center sample: ko16
#> Processing: ko15
#> Create profile matrix with method 'bin' and step 1 ... OK
#> Create profile matrix with method 'bin' and step 1 ... OK
#>
#> center sample: ko16
#> Processing: ko15
#> Create profile matrix with method 'bin' and step 1 ... OK
#> Create profile matrix with method 'bin' and step 1 ... OK
#>
#> center sample: ko16
#> Processing: ko15
#> Create profile matrix with method 'bin' and step 1 ... OK
#> Create profile matrix with method 'bin' and step 1 ... OK
#>
#> center sample: ko16
#> Processing: ko15
#> Create profile matrix with method 'bin' and step 1 ... OK
#> Create profile matrix with method 'bin' and step 1 ... OK
#>
#> center sample: ko16
#> Processing: ko15
#> Create profile matrix with method 'bin' and step 1 ... OK
#> Create profile matrix with method 'bin' and step 1 ... OK
#>
#> center sample: ko16
#> Processing: ko15
#> Create profile matrix with method 'bin' and step 1 ... OK
#> Create profile matrix with method 'bin' and step 1 ... OK
#>
#> center sample: ko16
#> Processing: ko15
#> Create profile matrix with method 'bin' and step 1 ... OK
#> Create profile matrix with method 'bin' and step 1 ... OK
#>
#> center sample: ko16
#> Processing: ko15
#> Create profile matrix with method 'bin' and step 1 ... OK
#> Create profile matrix with method 'bin' and step 1 ... OK
#>
#> center sample: ko16
#> Processing: ko15
#> Create profile matrix with method 'bin' and step 1 ... OK
#> Create profile matrix with method 'bin' and step 1 ... OK
#>
#> center sample: ko16
#> Processing: ko15
#> Create profile matrix with method 'bin' and step 1 ... OK
#> Create profile matrix with method 'bin' and step 1 ... OK
#>
#> center sample: ko16
#> Processing: ko15
#> Create profile matrix with method 'bin' and step 1 ... OK
#> Create profile matrix with method 'bin' and step 1 ... OK
#>
#> center sample: ko16
#> Processing: ko15
#> Create profile matrix with method 'bin' and step 1 ... OK
#> Create profile matrix with method 'bin' and step 1 ... OK
#>
#> center sample: ko16
#> Processing: ko15
#> Create profile matrix with method 'bin' and step 1 ... OK
#> Create profile matrix with method 'bin' and step 1 ... OK
#>
#> center sample: ko16
#> Processing: ko15
#> Create profile matrix with method 'bin' and step 1 ... OK
#> Create profile matrix with method 'bin' and step 1 ... OK
#>
#> center sample: ko16
#> Processing: ko15
#> Create profile matrix with method 'bin' and step 1 ... OK
#> Create profile matrix with method 'bin' and step 1 ... OK
#>
#> center sample: ko16
#> Processing: ko15
#> Create profile matrix with method 'bin' and step 1 ... OK
#> Create profile matrix with method 'bin' and step 1 ... OK
#>
#> center sample: ko16
#> Processing: ko15
#> Create profile matrix with method 'bin' and step 1 ... OK
#> Create profile matrix with method 'bin' and step 1 ... OK
#> center sample: ko16
#> Processing: ko15
#> Create profile matrix with method 'bin' and step 1 ... OK
#> Create profile matrix with method 'bin' and step 1 ... OK
#> profStep or minfrac greater 1, decreasing to 0.54 and 1
#>
#>
#>
#> starting new DoE with:
#>
#> gapInit: c(0.352, 0.928)
#> bw: c(0.879999999999999, 23.92)
#> minfrac: c(0.54, 1)
#> mzwid: c(0.0326, 0.0614)
#> distFunc: cor_opt
#> gapExtend: 2.7
#> profStep: 1
#> plottype: none
#> response: 1
#> factorDiag: 2
#> factorGap: 1
#> localAlignment: 0
#> retcorMethod: obiwarp
#> minsamp: 1
#> max: 50
#> center: 2
#> center sample: ko16
#> Processing: ko15
#> Create profile matrix with method 'bin' and step 1 ... OK
#> Create profile matrix with method 'bin' and step 1 ... OK
#>
#> center sample: ko16
#> Processing: ko15
#> Create profile matrix with method 'bin' and step 1 ... OK
#> Create profile matrix with method 'bin' and step 1 ... OK
#>
#> center sample: ko16
#> Processing: ko15
#> Create profile matrix with method 'bin' and step 1 ... OK
#> Create profile matrix with method 'bin' and step 1 ... OK
#>
#> center sample: ko16
#> Processing: ko15
#> Create profile matrix with method 'bin' and step 1 ... OK
#> Create profile matrix with method 'bin' and step 1 ... OK
#>
#> center sample: ko16
#> Processing: ko15
#> Create profile matrix with method 'bin' and step 1 ... OK
#> Create profile matrix with method 'bin' and step 1 ... OK
#>
#> center sample: ko16
#> Processing: ko15
#> Create profile matrix with method 'bin' and step 1 ... OK
#> Create profile matrix with method 'bin' and step 1 ... OK
#>
#> center sample: ko16
#> Processing: ko15
#> Create profile matrix with method 'bin' and step 1 ... OK
#> Create profile matrix with method 'bin' and step 1 ... OK
#>
#> center sample: ko16
#> Processing: ko15
#> Create profile matrix with method 'bin' and step 1 ... OK
#> Create profile matrix with method 'bin' and step 1 ... OK
#>
#> center sample: ko16
#> Processing: ko15
#> Create profile matrix with method 'bin' and step 1 ... OK
#> Create profile matrix with method 'bin' and step 1 ... OK
#>
#> center sample: ko16
#> Processing: ko15
#> Create profile matrix with method 'bin' and step 1 ... OK
#> Create profile matrix with method 'bin' and step 1 ... OK
#>
#> center sample: ko16
#> Processing: ko15
#> Create profile matrix with method 'bin' and step 1 ... OK
#> Create profile matrix with method 'bin' and step 1 ... OK
#>
#> center sample: ko16
#> Processing: ko15
#> Create profile matrix with method 'bin' and step 1 ... OK
#> Create profile matrix with method 'bin' and step 1 ... OK
#>
#> center sample: ko16
#> Processing: ko15
#> Create profile matrix with method 'bin' and step 1 ... OK
#> Create profile matrix with method 'bin' and step 1 ... OK
#>
#> center sample: ko16
#> Processing: ko15
#> Create profile matrix with method 'bin' and step 1 ... OK
#> Create profile matrix with method 'bin' and step 1 ... OK
#>
#> center sample: ko16
#> Processing: ko15
#> Create profile matrix with method 'bin' and step 1 ... OK
#> Create profile matrix with method 'bin' and step 1 ... OK
#>
#> center sample: ko16
#> Processing: ko15
#> Create profile matrix with method 'bin' and step 1 ... OK
#> Create profile matrix with method 'bin' and step 1 ... OK
#>
#> center sample: ko16
#> Processing: ko15
#> Create profile matrix with method 'bin' and step 1 ... OK
#> Create profile matrix with method 'bin' and step 1 ... OK
#>
#> center sample: ko16
#> Processing: ko15
#> Create profile matrix with method 'bin' and step 1 ... OK
#> Create profile matrix with method 'bin' and step 1 ... OK
#>
#> center sample: ko16
#> Processing: ko15
#> Create profile matrix with method 'bin' and step 1 ... OK
#> Create profile matrix with method 'bin' and step 1 ... OK
#>
#> center sample: ko16
#> Processing: ko15
#> Create profile matrix with method 'bin' and step 1 ... OK
#> Create profile matrix with method 'bin' and step 1 ... OK
#>
#> center sample: ko16
#> Processing: ko15
#> Create profile matrix with method 'bin' and step 1 ... OK
#> Create profile matrix with method 'bin' and step 1 ... OK
#>
#> center sample: ko16
#> Processing: ko15
#> Create profile matrix with method 'bin' and step 1 ... OK
#> Create profile matrix with method 'bin' and step 1 ... OK
#>
#> center sample: ko16
#> Processing: ko15
#> Create profile matrix with method 'bin' and step 1 ... OK
#> Create profile matrix with method 'bin' and step 1 ... OK
#>
#> center sample: ko16
#> Processing: ko15
#> Create profile matrix with method 'bin' and step 1 ... OK
#> Create profile matrix with method 'bin' and step 1 ... OK
#>
#> center sample: ko16
#> Processing: ko15
#> Create profile matrix with method 'bin' and step 1 ... OK
#> Create profile matrix with method 'bin' and step 1 ... OK
#>
#> center sample: ko16
#> Processing: ko15
#> Create profile matrix with method 'bin' and step 1 ... OK
#> Create profile matrix with method 'bin' and step 1 ... OK
#> center sample: ko16
#> Processing: ko15
#> Create profile matrix with method 'bin' and step 1 ... OK
#> Create profile matrix with method 'bin' and step 1 ... OK
#> no increase stopping
The response surface models of all optimization steps for the retention time correction and grouping parameters are shown above.
Currently the xcms
retention time correction method
obiwarp
and grouping method density
are
supported.
A script which you can use to process your raw data can be generated
by using the function writeRScript
.
writeRScript(resultPeakpicking$best_settings$parameters,
resultRetcorGroup$best_settings)
#> library(xcms)
#> library(Rmpi)
#> xset <- xcmsSet(
#> method = "matchedFilter",
#> fwhm = 50,
#> snthresh = 3,
#> step = 0.26,
#> steps = 2,
#> sigma = 21.2332257516562,
#> max = 5,
#> mzdiff = 0.28,
#> index = FALSE)
#> xset <- retcor(
#> xset,
#> method = "obiwarp",
#> plottype = "none",
#> distFunc = "cor_opt",
#> profStep = 1,
#> center = 2,
#> response = 1,
#> gapInit = 0.64,
#> gapExtend = 2.7,
#> factorDiag = 2,
#> factorGap = 1,
#> localAlignment = 0)
#> xset <- group(
#> xset,
#> method = "density",
#> bw = 12.4,
#> mzwid = 0.047,
#> minfrac = 0.94,
#> minsamp = 1,
#> max = 50)
#> xset <- fillPeaks(xset)
Above calculations proceeded with following running times.
time.xcmsSet # time for optimizing peak picking parameters
#> user system elapsed
#> 280.672 64.830 173.831
time.RetGroup # time for optimizing retention time correction and grouping parameters
#> user system elapsed
#> 494.455 2.143 495.591
sessionInfo()
#> R version 4.4.2 (2024-10-31)
#> Platform: x86_64-pc-linux-gnu
#> Running under: Ubuntu 24.04.1 LTS
#>
#> Matrix products: default
#> BLAS: /usr/lib/x86_64-linux-gnu/openblas-pthread/libblas.so.3
#> LAPACK: /usr/lib/x86_64-linux-gnu/openblas-pthread/libopenblasp-r0.3.26.so; LAPACK version 3.12.0
#>
#> locale:
#> [1] LC_CTYPE=en_US.UTF-8 LC_NUMERIC=C
#> [3] LC_TIME=en_US.UTF-8 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: Etc/UTC
#> tzcode source: system (glibc)
#>
#> attached base packages:
#> [1] stats graphics grDevices utils datasets methods base
#>
#> other attached packages:
#> [1] IPO_1.33.0 CAMERA_1.63.0 Biobase_2.67.0
#> [4] BiocGenerics_0.53.3 generics_0.1.3 rsm_2.10.5
#> [7] faahKO_1.46.0 xcms_4.5.1 BiocParallel_1.41.0
#> [10] rmarkdown_2.29
#>
#> loaded via a namespace (and not attached):
#> [1] DBI_1.2.3 RBGL_1.83.0
#> [3] gridExtra_2.3 rlang_1.1.4
#> [5] magrittr_2.0.3 clue_0.3-66
#> [7] MassSpecWavelet_1.73.0 matrixStats_1.4.1
#> [9] compiler_4.4.2 vctrs_0.6.5
#> [11] reshape2_1.4.4 stringr_1.5.1
#> [13] ProtGenerics_1.39.0 MetaboCoreUtils_1.15.0
#> [15] pkgconfig_2.0.3 crayon_1.5.3
#> [17] fastmap_1.2.0 backports_1.5.0
#> [19] XVector_0.47.0 utf8_1.2.4
#> [21] graph_1.85.0 UCSC.utils_1.3.0
#> [23] preprocessCore_1.69.0 purrr_1.0.2
#> [25] xfun_0.49 MultiAssayExperiment_1.33.1
#> [27] zlibbioc_1.52.0 cachem_1.1.0
#> [29] GenomeInfoDb_1.43.2 jsonlite_1.8.9
#> [31] progress_1.2.3 DelayedArray_0.33.2
#> [33] prettyunits_1.2.0 parallel_4.4.2
#> [35] cluster_2.1.6 R6_2.5.1
#> [37] bslib_0.8.0 stringi_1.8.4
#> [39] RColorBrewer_1.1-3 limma_3.63.2
#> [41] rpart_4.1.23 GenomicRanges_1.59.1
#> [43] jquerylib_0.1.4 Rcpp_1.0.13-1
#> [45] SummarizedExperiment_1.37.0 iterators_1.0.14
#> [47] knitr_1.49 base64enc_0.1-3
#> [49] IRanges_2.41.1 nnet_7.3-19
#> [51] Matrix_1.7-1 igraph_2.1.1
#> [53] tidyselect_1.2.1 rstudioapi_0.17.1
#> [55] abind_1.4-8 yaml_2.3.10
#> [57] doParallel_1.0.17 codetools_0.2-20
#> [59] affy_1.85.0 lattice_0.22-6
#> [61] tibble_3.2.1 plyr_1.8.9
#> [63] evaluate_1.0.1 foreign_0.8-87
#> [65] Spectra_1.17.1 pillar_1.9.0
#> [67] affyio_1.77.0 BiocManager_1.30.25
#> [69] MatrixGenerics_1.19.0 checkmate_2.3.2
#> [71] foreach_1.5.2 stats4_4.4.2
#> [73] MSnbase_2.33.2 MALDIquant_1.22.3
#> [75] ncdf4_1.23 hms_1.1.3
#> [77] S4Vectors_0.45.2 ggplot2_3.5.1
#> [79] munsell_0.5.1 scales_1.3.0
#> [81] MsExperiment_1.9.0 glue_1.8.0
#> [83] Hmisc_5.2-0 MsFeatures_1.15.0
#> [85] lazyeval_0.2.2 maketools_1.3.1
#> [87] tools_4.4.2 data.table_1.16.2
#> [89] mzID_1.45.0 sys_3.4.3
#> [91] QFeatures_1.17.0 vsn_3.75.0
#> [93] mzR_2.41.1 fs_1.6.5
#> [95] buildtools_1.0.0 XML_3.99-0.17
#> [97] grid_4.4.2 impute_1.81.0
#> [99] tidyr_1.3.1 MsCoreUtils_1.19.0
#> [101] colorspace_2.1-1 GenomeInfoDbData_1.2.13
#> [103] PSMatch_1.11.0 htmlTable_2.4.3
#> [105] Formula_1.2-5 cli_3.6.3
#> [107] fansi_1.0.6 S4Arrays_1.7.1
#> [109] dplyr_1.1.4 AnnotationFilter_1.31.0
#> [111] pcaMethods_1.99.0 gtable_0.3.6
#> [113] sass_0.4.9 digest_0.6.37
#> [115] SparseArray_1.7.2 htmlwidgets_1.6.4
#> [117] htmltools_0.5.8.1 lifecycle_1.0.4
#> [119] httr_1.4.7 statmod_1.5.0
#> [121] MASS_7.3-61