XCMS Parameter Optimization with IPO

Gunnar Libiseller, Thomas Riebenbauer
JOANNEUM RESEARCH Forschungsgesellschaft m.b.H., Graz, Austria

2019-01-04

Introduction

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:


IPO optimization process

Installation

# try http:// if https:// URLs are not supported
if (!requireNamespace("BiocManager", quietly=TRUE))
    install.packages("BiocManager")
BiocManager::install("IPO")

Installing main suggested packages

# for examples of peak picking parameter optimization:
BiocManager::install("msdata")
# for examples of optimization of retention time correction and grouping
# parameters:
BiocManager::install("faahKO")

Raw data

xcms handles the file processing hence all files can be used that can be processed by xcms.

datapath <- system.file("cdf", package = "faahKO")
datafiles <- list.files(datapath, recursive = TRUE, full.names = TRUE)

Optimize peak picking parameters

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:

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.

library(IPO)
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.

Optimize retention time correction and grouping parameters

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
#> Processing 53233 mz slices ... 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
#> Processing 53233 mz slices ... 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
#> Processing 53233 mz slices ... 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
#> Processing 53233 mz slices ... 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
#> Processing 53233 mz slices ... 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
#> Processing 53233 mz slices ... 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
#> Processing 53233 mz slices ... 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
#> Processing 53233 mz slices ... 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
#> Processing 22815 mz slices ... 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
#> Processing 22815 mz slices ... 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
#> Processing 22815 mz slices ... 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
#> Processing 22815 mz slices ... 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
#> Processing 22815 mz slices ... 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
#> Processing 22815 mz slices ... 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
#> Processing 22815 mz slices ... 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
#> Processing 22815 mz slices ... 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
#> Processing 31940 mz slices ... 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
#> Processing 31940 mz slices ... 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
#> Processing 31940 mz slices ... 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
#> Processing 31940 mz slices ... 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
#> Processing 31940 mz slices ... 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
#> Processing 31940 mz slices ... 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
#> Processing 31940 mz slices ... 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
#> Processing 53233 mz slices ... 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
#> Processing 22815 mz slices ... 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
#> Processing 31940 mz slices ... 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
#> Processing 22815 mz slices ... 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
#> Processing 34718 mz slices ... 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
#> Processing 34718 mz slices ... 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
#> Processing 34718 mz slices ... 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
#> Processing 34718 mz slices ... 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
#> Processing 34718 mz slices ... 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
#> Processing 34718 mz slices ... 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
#> Processing 34718 mz slices ... 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
#> Processing 34718 mz slices ... 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
#> Processing 16990 mz slices ... 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
#> Processing 16990 mz slices ... 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
#> Processing 16990 mz slices ... 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
#> Processing 16990 mz slices ... 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
#> Processing 16990 mz slices ... 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
#> Processing 16990 mz slices ... 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
#> Processing 16990 mz slices ... 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
#> Processing 16990 mz slices ... 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
#> Processing 22815 mz slices ... 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
#> Processing 22815 mz slices ... 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
#> Processing 22815 mz slices ... 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
#> Processing 22815 mz slices ... 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
#> Processing 22815 mz slices ... 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
#> Processing 22815 mz slices ... 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
#> Processing 22815 mz slices ... 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
#> Processing 34718 mz slices ... 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
#> Processing 16990 mz slices ... 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
#> Processing 22815 mz slices ... 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
#> Processing 16990 mz slices ... 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
#> Processing 24494 mz slices ... 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
#> Processing 24494 mz slices ... 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
#> Processing 24494 mz slices ... 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
#> Processing 24494 mz slices ... 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
#> Processing 24494 mz slices ... 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
#> Processing 24494 mz slices ... 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
#> Processing 24494 mz slices ... 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
#> Processing 24494 mz slices ... 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
#> Processing 13006 mz slices ... 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
#> Processing 13006 mz slices ... 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
#> Processing 13006 mz slices ... 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
#> Processing 13006 mz slices ... 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
#> Processing 13006 mz slices ... 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
#> Processing 13006 mz slices ... 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
#> Processing 13006 mz slices ... 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
#> Processing 13006 mz slices ... 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
#> Processing 16990 mz slices ... 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
#> Processing 16990 mz slices ... 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
#> Processing 16990 mz slices ... 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
#> Processing 16990 mz slices ... 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
#> Processing 16990 mz slices ... 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
#> Processing 16990 mz slices ... 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
#> Processing 16990 mz slices ... 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
#> Processing 24494 mz slices ... 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
#> Processing 13006 mz slices ... 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
#> Processing 16990 mz slices ... 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
#> Processing 15135 mz slices ... 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.

Display optimized settings

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)

Running times and session info

Above calculations proceeded with following running times.

time.xcmsSet # time for optimizing peak picking parameters
#>    user  system elapsed 
#> 309.396  41.480 185.205
time.RetGroup # time for optimizing retention time correction and grouping parameters
#>    user  system elapsed 
#> 804.388   2.472 829.740

sessionInfo()
#> R version 3.5.2 (2018-12-20)
#> Platform: x86_64-pc-linux-gnu (64-bit)
#> Running under: Ubuntu 16.04.5 LTS
#> 
#> Matrix products: default
#> BLAS: /home/biocbuild/bbs-3.8-bioc/R/lib/libRblas.so
#> LAPACK: /home/biocbuild/bbs-3.8-bioc/R/lib/libRlapack.so
#> 
#> 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       
#> 
#> attached base packages:
#> [1] stats4    parallel  stats     graphics  grDevices utils     datasets 
#> [8] methods   base     
#> 
#> other attached packages:
#>  [1] IPO_1.8.1           CAMERA_1.38.1       rsm_2.10           
#>  [4] faahKO_1.22.0       xcms_3.4.2          MSnbase_2.8.3      
#>  [7] ProtGenerics_1.14.0 S4Vectors_0.20.1    mzR_2.16.1         
#> [10] Rcpp_1.0.0          BiocParallel_1.16.5 Biobase_2.42.0     
#> [13] BiocGenerics_0.28.0
#> 
#> loaded via a namespace (and not attached):
#>  [1] doParallel_1.0.14      RColorBrewer_1.1-2     backports_1.1.3       
#>  [4] tools_3.5.2            R6_2.3.0               affyio_1.52.0         
#>  [7] rpart_4.1-13           Hmisc_4.1-1            lazyeval_0.2.1        
#> [10] colorspace_1.3-2       nnet_7.3-12            tidyselect_0.2.5      
#> [13] gridExtra_2.3          emmeans_1.3.1          compiler_3.5.2        
#> [16] MassSpecWavelet_1.48.1 preprocessCore_1.44.0  graph_1.60.0          
#> [19] htmlTable_1.13         sandwich_2.5-0         checkmate_1.8.5       
#> [22] scales_1.0.0           DEoptimR_1.0-8         mvtnorm_1.0-8         
#> [25] robustbase_0.93-3      affy_1.60.0            RBGL_1.58.1           
#> [28] stringr_1.3.1          digest_0.6.18          foreign_0.8-71        
#> [31] rmarkdown_1.11         base64enc_0.1-3        pkgconfig_2.0.2       
#> [34] htmltools_0.3.6        limma_3.38.3           htmlwidgets_1.3       
#> [37] rlang_0.3.0.1          rstudioapi_0.8         impute_1.56.0         
#> [40] bindr_0.1.1            zoo_1.8-4              mzID_1.20.1           
#> [43] acepack_1.4.1          dplyr_0.7.8            magrittr_1.5          
#> [46] Formula_1.2-3          MALDIquant_1.18        Matrix_1.2-15         
#> [49] munsell_0.5.0          vsn_3.50.0             stringi_1.2.4         
#> [52] multcomp_1.4-8         yaml_2.2.0             MASS_7.3-51.1         
#> [55] zlibbioc_1.28.0        plyr_1.8.4             grid_3.5.2            
#> [58] crayon_1.3.4           lattice_0.20-38        splines_3.5.2         
#> [61] multtest_2.38.0        knitr_1.21             pillar_1.3.1          
#> [64] igraph_1.2.2           estimability_1.3       codetools_0.2-16      
#> [67] XML_3.98-1.16          glue_1.3.0             evaluate_0.12         
#> [70] latticeExtra_0.6-28    data.table_1.11.8      pcaMethods_1.74.0     
#> [73] BiocManager_1.30.4     foreach_1.4.4          gtable_0.2.0          
#> [76] RANN_2.6               purrr_0.2.5            assertthat_0.2.0      
#> [79] ggplot2_3.1.0          xfun_0.4               xtable_1.8-3          
#> [82] coda_0.19-2            ncdf4_1.16             survival_2.43-3       
#> [85] tibble_2.0.0           iterators_1.0.10       IRanges_2.16.0        
#> [88] bindrcpp_0.2.2         cluster_2.0.7-1        TH.data_1.0-9