Title: | Automated Optimization of XCMS Data Processing parameters |
---|---|
Description: | The outcome of XCMS data processing strongly depends on the parameter settings. IPO (`Isotopologue Parameter Optimization`) is a parameter optimization tool that is applicable for different kinds of samples and liquid chromatography coupled to high resolution mass spectrometry devices, fast and free of labeling steps. IPO uses natural, stable 13C isotopes to calculate a peak picking score. Retention time correction is optimized by minimizing the relative retention time differences within features and grouping parameters are optimized by maximizing the number of features showing exactly one peak from each injection of a pooled sample. The different parameter settings are achieved by design of experiment. The resulting scores are evaluated using response surface models. |
Authors: | Gunnar Libiseller <[email protected]>, Christoph Magnes <[email protected]>, Thomas Lieb <[email protected]> |
Maintainer: | Thomas Lieb <[email protected]> |
License: | GPL (>= 2) + file LICENSE |
Version: | 1.33.0 |
Built: | 2024-10-31 03:48:54 UTC |
Source: | https://github.com/bioc/IPO |
IPO provides a framework for parameter optimization for the software package XCMS. It provides optimisation of peak picking parameters by using natural, stable 13C isotopes. Retention time correction is optimized by minimizing the relative retention time differences within features and grouping parameters are optimized by maximizing the number of features showing exactly one peak from each injection of a pooled sample.
An overview of how to use the package, including the most important functions
Gunnar Libiseller
Maintainer: Thomas Riebenbauer <[email protected]>
Lenth, R. V. (2009). Response-Surface Methods in R , Using rsm. Journal of Statistical Software, 32(7), 1-17. Retrieved from http://www.jstatsoft.org/v32/i07
Smith, C.A. and Want, E.J. and O'Maille, G. and Abagyan,R. and Siuzdak, G.: XCMS: Processing mass spectrometry data for metabolite profiling using nonlinear peak alignment, matching and identification, Analytical Chemistry, 78:779-787 (2006)
Ralf Tautenhahn, Christoph Boettcher, Steffen Neumann: Highly sensitive feature detection for high resolution LC/MS BMC Bioinformatics, 9:504 (2008)
H. Paul Benton, Elizabeth J. Want and Timothy M. D. Ebbels Correction of mass calibration gaps in liquid chromatography-mass spectrometry metabolomics data Bioinformatics, 26:2488 (2010)
Yu, H. (2002). Rmpi: Parallel Statistical Computing in R. R News, 2(2), 10-14. Retrieved from http://cran.r-project.org/doc/Rnews/Rnews_2002-2.pdf
## Not run: mtbls2files <- list.files(file.path(find.package("mtbls2"), "mzML"), full.names=TRUE) paramsPP <- getDefaultXcmsSetStartingParams() paramsPP$mzdiff <- -0.001 #paramsPP$ppm <- 25 paramsPP$min_peakwidth <- c(7,14) paramsPP$max_peakwidth <- c(20,30) paramsPP$noise <- 10000 resultPP <- optimizeXcmsSet(mtbls2files[1:2], paramsPP, subdir="mtbls2") paramsRG <- getDefaultRetGroupStartingParams() paramsRG$gapInit <- 0.2 paramsRG$profStep <- 1 paramsRG$minfrac <- 0.75 resultRG <- optimizeRetGroup(resultPP$best_settings$xset, paramsRG, nSlaves=2) writeRScript(resultPP$best_settings$parameters, resultRG$best_settings, subdir="mtbls2", 4) ## End(Not run)
## Not run: mtbls2files <- list.files(file.path(find.package("mtbls2"), "mzML"), full.names=TRUE) paramsPP <- getDefaultXcmsSetStartingParams() paramsPP$mzdiff <- -0.001 #paramsPP$ppm <- 25 paramsPP$min_peakwidth <- c(7,14) paramsPP$max_peakwidth <- c(20,30) paramsPP$noise <- 10000 resultPP <- optimizeXcmsSet(mtbls2files[1:2], paramsPP, subdir="mtbls2") paramsRG <- getDefaultRetGroupStartingParams() paramsRG$gapInit <- 0.2 paramsRG$profStep <- 1 paramsRG$minfrac <- 0.75 resultRG <- optimizeRetGroup(resultPP$best_settings$xset, paramsRG, nSlaves=2) writeRScript(resultPP$best_settings$parameters, resultRG$best_settings, subdir="mtbls2", 4) ## End(Not run)
This function attaches one list at the end of another list.
attachList(params_1, params_2)
attachList(params_1, params_2)
params_1 |
A List |
params_2 |
A second list which will be attached at the end of the first list. |
This is a convenience funktion, but the implementation is not optimized for speed.
A List composed of the two input lists.
Gunnar Libiseller
a <- list("a"=1, "b"=2) b <- list("c"=4, "d"=4) attachList(a, b)
a <- list("a"=1, "b"=2) b <- list("c"=4, "d"=4) attachList(a, b)
This function calculates PPS by identifying natural, stable 13C isotopes within
an xcmsSet object. Peaks beeing part of an isotopologue are defined
as reliable peaks (RP). Peaks which are not part of an isotopologue
and where the intensity of possible isotopes is below a cutoff are defined as
'low intensity peaks' (LIP). PPS is then calculated by:
PPS = RP^2/(#all_peaks - LIP)
calcPPS(xset, isotopeIdentification, ...)
calcPPS(xset, isotopeIdentification, ...)
xset |
|
isotopeIdentification |
This parameter defines the method for isotope identification. The method IPO was especially implemented for high resolution data. CAMERA is an established isotope and adduct annotation package. |
... |
Additional parameters to CAMERA's findIsotopes function. |
Calculation of a peak picking score (PPS) by using natural, stable 13C isotopes
An array with 5 items:
1. Space for experimentid of the Design of Experiments (0 since not known in calcPPS)
2. Number of peaks
3. Number of peaks without LIP and RP
4. Reliable peaks (RP)
5. Peak picking score (PPS)
Gunnar Libiseller
findIsotopes.IPO
findIsotopes.CAMERA
mzmlfile <- file.path(find.package("msdata"), "microtofq/MM14.mzML") xset <- xcmsSet(mzmlfile, peakwidth=c(5,12), method="centWave") calcPPS(xset)
mzmlfile <- file.path(find.package("msdata"), "microtofq/MM14.mzML") xset <- xcmsSet(mzmlfile, peakwidth=c(5,12), method="centWave") calcPPS(xset)
This function encapsulates xcms::findPeaks-methods for IPO
calculateXcmsSet(files, xcmsSetParameters, scanrange=NULL, task=1, BPPARAM = bpparam(), nSlaves=0)
calculateXcmsSet(files, xcmsSetParameters, scanrange=NULL, task=1, BPPARAM = bpparam(), nSlaves=0)
files |
a vector containing the files for peak picking |
xcmsSetParameters |
a list with all parameters for |
scanrange |
scan range to process. See |
task |
The task-id when using this method in parallel calculations. |
BPPARAM |
a |
nSlaves |
|
Encapsulation of xcms::findPeaks-methods used in IPO.
An xcmsSet-object
Gunnar Libiseller, Thomas Riebenbauer ([email protected])
Smith, C.A. and Want, E.J. and O'Maille, G. and Abagyan,R. and Siuzdak, G.: XCMS: Processing mass spectrometry data for metabolite profiling using nonlinear peak alignment, matching and identification, Analytical Chemistry, 78:779-787 (2006)
Ralf Tautenhahn, Christoph Boettcher, Steffen Neumann: Highly sensitive feature detection for high resolution LC/MS BMC Bioinformatics, 9:504 (2008)
mzmlfile <- file.path(find.package("msdata"), "microtofq/MM14.mzML") params <- list(min_peakwidth=5, max_peakwidth=12, ppm=58, mzdiff=-0.001, snthresh=10, noise=0, prefilter=3, value_of_prefilter=100, mzCenterFun="wMean", integrate=1, fitgauss=FALSE, verbose.columns=FALSE, nSlaves=1) xset <- calculateXcmsSet(mzmlfile, params)
mzmlfile <- file.path(find.package("msdata"), "microtofq/MM14.mzML") params <- list(min_peakwidth=5, max_peakwidth=12, ppm=58, mzdiff=-0.001, snthresh=10, noise=0, prefilter=3, value_of_prefilter=100, mzCenterFun="wMean", integrate=1, fitgauss=FALSE, verbose.columns=FALSE, nSlaves=1) xset <- calculateXcmsSet(mzmlfile, params)
This function combines two lists of parameters. The first is a list of parameters which should be optimized. These parameters have different values set by Design of Experimnt. The second list consists of parameters which should not be optimized, hence only one value is set for each parameter. The parameters of the second list are replicated to have the same length as the number experiments in the DoE. Then the two lists are combined.
combineParams(params_1, params_2)
combineParams(params_1, params_2)
params_1 |
A list holding parameters which should be optimized. Each parameter already has value set for each experiment of an Design of Experiment. |
params_2 |
A list holding parameters which should not be optimized, hence only one value is set. |
Special treatment is needed for the findPeaks.matchedFilter-parameters 'sigma', 'mzdiff'
since these two parameters are defined by default relative to the parameters 'fwhm' or
'step' and 'steps' respectively.sigma=fwhm/2.3548
mzdiff=0.8-step*steps
A list of consting of all parameters needed for an xcms-method (findPeaks.centWave, findPeaks.matchedFilter, retcor.obiwarp or group.density). Each list item has the same length which is equal to the number of experiments within the DoE.
Gunnar Libiseller
params <- getDefaultXcmsSetStartingParams() typ_params <- typeCastParams(params) design <- getBbdParameter(typ_params$to_optimize) xcms_design <- decode.data(design) xcms_design <- combineParams(xcms_design, typ_params$no_optimization) xcms_design
params <- getDefaultXcmsSetStartingParams() typ_params <- typeCastParams(params) design <- getBbdParameter(typ_params$to_optimize) xcms_design <- decode.data(design) xcms_design <- combineParams(xcms_design, typ_params$no_optimization) xcms_design
This function uses a design of experiments, a response for the experiments within the design and the used parameters to create a response surface model
createModel(design, params, resp)
createModel(design, params, resp)
design |
A design of experiments (Box-Behnken-Design or Central-Composite-Design) |
params |
The parameters which were used. |
resp |
The responses achivied for the various experiments. |
This function uses a design of experiments, a response for the experiments within the design and the used parameters to create a response surface model
A response surface model.
getBbdParameter
getCcdParameter
typeCastParams
Gunnar Libiseller
Lenth, R. V. (2009). Response-Surface Methods in R , Using rsm. Journal of Statistical Software, 32(7), 1-17. Retrieved from http://www.jstatsoft.org/v32/i07
params <- getDefaultXcmsSetStartingParams() type_params <- typeCastParams(params) design <- getBbdParameter(type_params$to_optimize) resp <- runif(nrow(design),1,3) model <- createModel(design, type_params$to_optimize, resp) dev.new() par(mfrow=c(3,2)) contour(model, ~ x1*x2*x3*x4, image=TRUE)
params <- getDefaultXcmsSetStartingParams() type_params <- typeCastParams(params) design <- getBbdParameter(type_params$to_optimize) resp <- runif(nrow(design),1,3) model <- createModel(design, type_params$to_optimize, resp) dev.new() par(mfrow=c(3,2)) contour(model, ~ x1*x2*x3*x4, image=TRUE)
Encode and decode values that are in a range of -1 to 1 into a specified range.
encode(value, bounds) decode(value, bounds) decodeAll(values, params)
encode(value, bounds) decode(value, bounds) decodeAll(values, params)
value |
A value |
values |
A vector with values in the range [-1,1] |
bounds |
A vector of two values defining the lower and upper bound of a range. |
params |
A list where evere list-item consist of two values defining a lower and an upper bound. |
Decodes a values from ranges of -1 to 1 to ranges specified.
A function used to decode values that are in a range of -1 to 1 into a specified range. For every value a list item with lower and upper bound has to be supplied.
A function used to encode values that are in a specified range into a range between -1 to 1.
decode: The encoded value. decodeAll: A vector of decoded values.
Gunnar Libiseller
decode(0, c(10, 20)) decode(-0.5, c(10, 20)) decode(1, c(10, 20)) bounds <- c(10, 20) encode(decode(1, bounds), bounds) ## Multiple values: values <- c(-1, -0.25, 0, 0.75) params <- getDefaultXcmsSetStartingParams() type_params <- typeCastParams(params) decodeAll(values, type_params$to_optimize) ## Combination of encode and decode encode(15, c(10, 20)) encode(10, c(10, 20)) encode(5, c(1, 5)) bounds <- c(1,5) decode(encode(5, bounds), bounds)
decode(0, c(10, 20)) decode(-0.5, c(10, 20)) decode(1, c(10, 20)) bounds <- c(10, 20) encode(decode(1, bounds), bounds) ## Multiple values: values <- c(-1, -0.25, 0, 0.75) params <- getDefaultXcmsSetStartingParams() type_params <- typeCastParams(params) decodeAll(values, type_params$to_optimize) ## Combination of encode and decode encode(15, c(10, 20)) encode(10, c(10, 20)) encode(5, c(1, 5)) bounds <- c(1,5) decode(encode(5, bounds), bounds)
This function finds isotopes using CAMERA's find peak function. Isotopes are separately found within each sample.
findIsotopes.CAMERA(xset, ...)
findIsotopes.CAMERA(xset, ...)
xset |
|
... |
Additional parameters to the findIsotopes function of CAMERA |
Identification of 13C isotopes
An matrix with 2 columns. Column one shows the peak id of the 12C, peak column two shows the id of the respective 13C isotope peak.
Gunnar Libiseller
C. Kuhl and R. Tautenhahn and C. Boettcher and T. R. Larson and S. Neumann: CAMERA: an integrated strategy for compound spectra extraction and annotation of liquid chromatography/mass spectrometry data sets Analytical Chemistry 84:283 (2012)
mzmlfile <- file.path(find.package("msdata"), "microtofq/MM14.mzML") xset <- xcmsSet(mzmlfile, peakwidth=c(5,12), method="centWave") isotopes <- findIsotopes.CAMERA(xset, ppm=15, maxcharge=1)
mzmlfile <- file.path(find.package("msdata"), "microtofq/MM14.mzML") xset <- xcmsSet(mzmlfile, peakwidth=c(5,12), method="centWave") isotopes <- findIsotopes.CAMERA(xset, ppm=15, maxcharge=1)
This function identifies natural, stable 13C isotopes within an xcmsSet object of LC-HRMS data. Isotopes have to be within a mass-, retentiontime- and intensitywindow to be recognized as isotopes. If checkBorderIntensity is TRUE the maximum intensity of each peaks has to be at least three times the intensity at rtmin and rtmax.
findIsotopes.IPO(xset, checkPeakShape=c("none", "borderIntensity", "sinusCurve", "normalDistr"))
findIsotopes.IPO(xset, checkPeakShape=c("none", "borderIntensity", "sinusCurve", "normalDistr"))
xset |
|
checkPeakShape |
character to choose if the peakshape should be checked and if so how |
Identification of 13C isotopes
An matrix with 2 columns. Column one shows the peak id of the 12C, peak column two shows the id of the respective 13C isotope peak.
Gunnar Libiseller
mzmlfile <- file.path(find.package("msdata"), "microtofq/MM14.mzML") xset <- xcmsSet(mzmlfile, peakwidth=c(5,12), method="centWave") isotopes <- findIsotopes.IPO(xset, "borderIntensity")
mzmlfile <- file.path(find.package("msdata"), "microtofq/MM14.mzML") xset <- xcmsSet(mzmlfile, peakwidth=c(5,12), method="centWave") isotopes <- findIsotopes.IPO(xset, "borderIntensity")
Creates a Box-Behnken Design of Experiment out of a list of parameters. Each of the list items has to be a pair defining the lower und upper limits of the value-range to test. The method then returns a Center faced Box-Behnken Design of Experiments. The list has to hold a least three pairs.
getBbdParameter(params)
getBbdParameter(params)
params |
A list of value pairs defining lower und upper limits of an optimization range. |
Creates a Box-Behnken Design of Experiment out of a list of parameters. Each of the list items has to be a pair defining the lower und upper limits of the value-range to test. The method then returns a Center faced Box-Behnken Design of Experiments. The list has to hold a least three pairs.
A Box-Behnken Design of Experiments
Gunnar Libiseller
Lenth, R. V. (2009). Response-Surface Methods in R , Using rsm. Journal of Statistical Software, 32(7), 1-17. Retrieved from http://www.jstatsoft.org/v32/i07
params <- getDefaultXcmsSetStartingParams() typ_params <- typeCastParams(params) design <- getBbdParameter(typ_params$to_optimize)
params <- getDefaultXcmsSetStartingParams() typ_params <- typeCastParams(params) design <- getBbdParameter(typ_params$to_optimize)
Creates a Central-Composite Design of Experiment out of a list of parameters. Each of the list items has to be a pair defining the lower und upper limits of the value-range to test. The method then returns a Center faced Central-Composite Design of Experiments.
getCcdParameter(params)
getCcdParameter(params)
params |
A list of value pairs defining lower und upper limits of an optimization range. |
Creates a Central-Composite Design of Experiment out of a list of parameters. Each of the list items has to be a pair defining the lower und upper limits of the value-range to test. The method then returns a Center faced Central-Composite Design of Experiments.
A Central-Composite Design of Experiments
Gunnar Libiseller
Lenth, R. V. (2009). Response-Surface Methods in R , Using rsm. Journal of Statistical Software, 32(7), 1-17. Retrieved from http://www.jstatsoft.org/v32/i07
params <- getDefaultXcmsSetStartingParams() typ_params <- typeCastParams(params) design <- getCcdParameter(typ_params$to_optimize)
params <- getDefaultXcmsSetStartingParams() typ_params <- typeCastParams(params) design <- getCcdParameter(typ_params$to_optimize)
Gets the index of the sample with most peaks in it. This is used if no center sample for retention time correction has been defined by the user.
getDefaultRetCorCenterSample(xset)
getDefaultRetCorCenterSample(xset)
xset |
|
Gets the index of the sample with most peaks in it. This is used if no center sample for retention time correction has been defined by the user.
The file index of the sample with most peaks in it.
Gunnar Libiseller
## The function is currently defined as function (xset) { ret <- NULL for (i in 1:length(filepaths(xset))) { ret <- c(ret, sum(peaks(xset)[, "sample"] == i)) } return(which.max(ret)) }
## The function is currently defined as function (xset) { ret <- NULL for (i in 1:length(filepaths(xset))) { ret <- c(ret, sum(peaks(xset)[, "sample"] == i)) } return(which.max(ret)) }
This function creates a list of parameters used in the xcms-methods
retcor.obiwarp and group.density. Per default the following parameters
have a defined range where optimization should start:
retcor.obiwarp parameters: 'gapInit'; 'gapExtend', 'profStep'
group.density parameters: 'bw', 'minfrac', 'mzwid'
getDefaultRetGroupStartingParams(retcorMethod=c("obiwarp", "loess", "none"), distfunc=c("cor_opt", "cor", "cov", "prd", "euc"), high_resolution=TRUE)
getDefaultRetGroupStartingParams(retcorMethod=c("obiwarp", "loess", "none"), distfunc=c("cor_opt", "cor", "cov", "prd", "euc"), high_resolution=TRUE)
retcorMethod |
The name of the retention time correction method that should be used. The XCMS methods retcor.obiwarp and retcor.loess are supported. If no retention time correction should be done use "none". |
distfunc |
The name of the distance function used by retcor.obiwarp |
high_resolution |
If high_resolution = TRUE starting values for mzwid are set to 0.015 and 0.035; if high_resolution = FALSE to 0.15, 0.35 |
* Do not delete a parameter from the list returned.
* Optimization of qualitative parameters is not supported yet.
* If you want to optimize additional parameter just set an lower and
an upper bound (e.g. params$max <- c(4,8))
* If you dont want to optimize a parameter set a default value
(e.g. params$max <- 10)
A List of parameters used in the xcms-methods retcor.obiwarp or retcor.loess and group.density
Gunnar Libiseller
params <- getDefaultRetGroupStartingParams() params$bw <- 10 params$max <- c(4,8) params
params <- getDefaultRetGroupStartingParams() params$bw <- 10 params$max <- c(4,8) params
This function creates a list of parameters used in the xcmsSet.findPeak-methods
'centWave' and 'matchedFilter'. Per default the following parameters
have a defined range where optimization should start:
'centWave' parameters: 'peakwidth' (split into 'min_peakwidth' and 'max_peakwidth'),
'ppm', 'mzdiff'
'matchedFilter' parameters: 'fwhm', 'snthresh', 'step', 'steps'
getDefaultXcmsSetStartingParams(method = c("centWave", "matchedFilter"))
getDefaultXcmsSetStartingParams(method = c("centWave", "matchedFilter"))
method |
Either parameters for 'centWave' or 'matchedFilter' should be created |
* Do not delete a parameter from the list returned.
* Optimization of qualitative parameters is not supported yet.
* If you want to optimize additional parameter just set an lower and
an upper bound (e.g. params$snthresh <- c(5,20))
* If you dont want to optimize a parameter set a default value
(e.g. params$snthresh <- 10)
A List of parameters for the xcmsSet.findPeak-methods 'centWave' or 'matchedFilter'
Gunnar Libiseller
params <- getDefaultXcmsSetStartingParams() params$ppm <- 10 params$snthresh <- c(5,15) params params <- getDefaultXcmsSetStartingParams("matchedFilter") params
params <- getDefaultXcmsSetStartingParams() params$ppm <- 10 params$snthresh <- c(5,15) params params <- getDefaultXcmsSetStartingParams("matchedFilter") params
This function does unity based normalization on Retention time Correction Scores (RCS) as well as Grouping Scores (GS).
getNormalizedResponse(response)
getNormalizedResponse(response)
response |
A List of all responses calculated by getRGTVValues for all experiments of an Design of Experiment |
Grouping Score (GS) is calculated by:
'good groups'^2/'bad groups
For all RCS and GS values unitiy based normalization is done. For every experiment within the DoE these two values are added together and returned.
A vector with RTGV values
Since RCS and GS can be within completely different ranges, normalization has to be done to prevent an excessive influence of either RCS or GS.
Gunnar Libiseller
mtbls2files <- list.files(file.path(find.package("mtbls2"), "mzML"), full.names=TRUE) params <- list(min_peakwidth=12, max_peakwidth=30, ppm=30, mzdiff=-0.001, snthresh=10, noise=10000, prefilter=3, value_of_prefilter=100, mzCenterFun="wMean", integrate=1, fitgauss=FALSE, verbose.columns=FALSE, nSlaves=2) xset <- calculateXcmsSet(mtbls2files[1:2], params) xset <- retcor(xset, method="obiwarp") xset <- group(xset) result <- getRGTVValues(xset) result
mtbls2files <- list.files(file.path(find.package("mtbls2"), "mzML"), full.names=TRUE) params <- list(min_peakwidth=12, max_peakwidth=30, ppm=30, mzdiff=-0.001, snthresh=10, noise=10000, prefilter=3, value_of_prefilter=100, mzCenterFun="wMean", integrate=1, fitgauss=FALSE, verbose.columns=FALSE, nSlaves=2) xset <- calculateXcmsSet(mtbls2files[1:2], params) xset <- retcor(xset, method="obiwarp") xset <- group(xset) result <- getRGTVValues(xset) result
This function calculates the Retention time Correction Score (RCS) of all features within an xcmsSet-object. Also features having exactly one peak from each sample are defined as 'good groups', all others a 'bad groups'.
getRGTVValues(xset, exp_index = 1, retcor_penalty = 1)
getRGTVValues(xset, exp_index = 1, retcor_penalty = 1)
xset |
|
exp_index |
Experiment-id of the experiment within a Design of Experiments |
retcor_penalty |
Penalty if an error occured with the used retention time correction parameters |
This function calculates the Retention time Correction Score (RCS) of all features within an xcmsSet-object. Also features having exactly one peak from each sample are defined as 'good groups', all others a 'bad groups' which leads to a Grouping Score (GS) by calculating 'good groups'^2/'bad groups'.
a list containing the items exp_index, good_groups, bad_groups, GS and RCS.
Gunnar Libiseller
mtbls2files <- list.files(paste(find.package("mtbls2"), "/mzML", sep=""), full.names=TRUE) xset <- xcmsSet(mtbls2files[1:2], method="centWave", peakwidth=(c(12, 30)), ppm=30, noise=10000) xset <- retcor(xset, method="obiwarp") xset <- group(xset) getRGTVValues(xset)
mtbls2files <- list.files(paste(find.package("mtbls2"), "/mzML", sep=""), full.names=TRUE) xset <- xcmsSet(mtbls2files[1:2], method="centWave", peakwidth=(c(12, 30)), ppm=30, noise=10000) xset <- retcor(xset, method="obiwarp") xset <- group(xset) getRGTVValues(xset)
This function provides optimisation for parameters of the xcms-method retcor.obiwarp and group.density. The retention time correction is optimised by minimizing intra-feature retention time shifts; grouping is optimized by increasing the number of features which have exactly one peak per sample.
optimizeRetGroup(xset, params = getDefaultRetGroupStartingParams(), nSlaves = 4, subdir = "IPO", plot = TRUE)
optimizeRetGroup(xset, params = getDefaultRetGroupStartingParams(), nSlaves = 4, subdir = "IPO", plot = TRUE)
xset |
|
params |
A list of parameters which are needed by xcms-methods retcor.obiwarp and group.density. List-items with two values will be optimized. The first value defines the lower test value, the second one the higher test value. |
nSlaves |
Number of slaves the optimization process should spawn. |
subdir |
The name of the subdirectory which is created and where the figures of the response
surface models will be saved to. |
plot |
Defines if plots should be generated ( |
This function provides optimisation for parameters of the xcms-method retcor.obiwarp and group.density. The retention time correction is optimised by minimizing intra-feature retention time shifts; grouping is optimized by increasing the number of features which have exactly one peak per sample.
A LIST of length n+1 with n beeing the optimization runs needed
comp1-comp(n) |
A LIST containing: |
comp(n+1) |
A LIST containing: |
Gunnar Libiseller
Obiwarp
Prince, J. T., & Marcotte, E. M. (2006). Chromatographic alignment of
ESI-LC-MS proteomics data sets by ordered bijective interpolated warping.
Analytical chemistry, 78(17), 6140-52. doi:10.1021/ac0605344
getDefaultRetGroupStartingParams
mtbls2files <- list.files(file.path(find.package("mtbls2"), "mzML"), full.names=TRUE) params <- list(min_peakwidth=12, max_peakwidth=30, ppm=30, mzdiff=-0.001, snthresh=10, noise=10000, prefilter=3, value_of_prefilter=100, mzCenterFun="wMean", integrate=1, fitgauss=FALSE, verbose.columns=FALSE, nSlaves=2) xset <- calculateXcmsSet(mtbls2files[1:2], params) #optimize the retention time correction and grouping parameters paramsRG <- getDefaultRetGroupStartingParams() paramsRG$profStep <- 1 paramsRG$minfrac <- 0.75 resultRG <- optimizeRetGroup(xset, params=paramsRG, nSlaves=4,subdir="mtbls2") writeRScript(params, resultRG$best_settings, 4)
mtbls2files <- list.files(file.path(find.package("mtbls2"), "mzML"), full.names=TRUE) params <- list(min_peakwidth=12, max_peakwidth=30, ppm=30, mzdiff=-0.001, snthresh=10, noise=10000, prefilter=3, value_of_prefilter=100, mzCenterFun="wMean", integrate=1, fitgauss=FALSE, verbose.columns=FALSE, nSlaves=2) xset <- calculateXcmsSet(mtbls2files[1:2], params) #optimize the retention time correction and grouping parameters paramsRG <- getDefaultRetGroupStartingParams() paramsRG$profStep <- 1 paramsRG$minfrac <- 0.75 resultRG <- optimizeRetGroup(xset, params=paramsRG, nSlaves=4,subdir="mtbls2") writeRScript(params, resultRG$best_settings, 4)
This function provides optimisation of peak picking parameters by using natural, stable 13C isotopes.
optimizeXcmsSet(files, params = getDefaultXcmsSetStartingParams(), isotopeIdentification = c("IPO", "CAMERA"), BPPARAM = bpparam(), nSlaves = 4, subdir = "IPO", plot = TRUE, ...)
optimizeXcmsSet(files, params = getDefaultXcmsSetStartingParams(), isotopeIdentification = c("IPO", "CAMERA"), BPPARAM = bpparam(), nSlaves = 4, subdir = "IPO", plot = TRUE, ...)
files |
A directory or list of files, passed to |
params |
A list of parameters which are needed by XCMS::findPeaks-Methods. List-items with two values will be optimized. The first value defines the lower test value, the second one the higher test value. |
isotopeIdentification |
This parameter defines the method for isotope identification. The method 'IPO' was especially implemented for high resolution data. CAMERA is an established isotope and adduct annotation package. |
BPPARAM |
a |
nSlaves |
Number of slaves the optimization process should spawn. |
subdir |
The name of the subdirectory which is created and where the figures of the response
surface models will be saved to. |
plot |
Defines if plots should be generated ( |
... |
Additional parameters to CAMERA's or IPO's findIsotopes functions |
This function provides optimisation of peak picking parameters by using natural, stable 13C isotopes.
A LIST of length n+1 with n beeing the optimization runs (DoEs) needed
comp1-comp(n) |
A LIST containing: |
comp(n+1) |
A LIST containing: |
Gunnar Libiseller, Thomas Riebenbauer ([email protected])
Smith, C.A. and Want, E.J. and O'Maille, G. and Abagyan,R. and Siuzdak, G.: XCMS: Processing mass spectrometry data for metabolite profiling using nonlinear peak alignment, matching and identification, Analytical Chemistry, 78:779-787 (2006)
Ralf Tautenhahn, Christoph Boettcher, Steffen Neumann: Highly sensitive feature detection for high resolution LC/MS BMC Bioinformatics, 9:504 (2008)
H. Paul Benton, Elizabeth J. Want and Timothy M. D. Ebbels: Correction of mass calibration gaps in liquid chromatography-mass spectrometry metabolomics data Bioinformatics, 26:2488 (2010)
C. Kuhl and R. Tautenhahn and C. Boettcher and T. R. Larson and S. Neumann: CAMERA: an integrated strategy for compound spectra extraction and annotation of liquid chromatography/mass spectrometry data sets Analytical Chemistry 84:283 (2012)
getDefaultXcmsSetStartingParams
,
calcPPS
,
findIsotopes.IPO
,
findIsotopes.CAMERA
#library(IPO) mzmlfile <- file.path(find.package("msdata"), "microtofq/MM14.mzML") paramsPP <- getDefaultXcmsSetStartingParams() paramsPP$mzdiff <- -0.001 paramsPP$min_peakwidth <- c(7,14) paramsPP$max_peakwidth <- c(20,30) #example using IPO isotope identification resultPP <- optimizeXcmsSet(mzmlfile, paramsPP, subdir="mtbls2") #example using CAMERA isotope identification resultPP <- optimizeXcmsSet(mzmlfile, paramsPP, isotopeIdentification="CAMERA", subdir="mtbls2", ppm=15, maxcharge=2)
#library(IPO) mzmlfile <- file.path(find.package("msdata"), "microtofq/MM14.mzML") paramsPP <- getDefaultXcmsSetStartingParams() paramsPP$mzdiff <- -0.001 paramsPP$min_peakwidth <- c(7,14) paramsPP$max_peakwidth <- c(20,30) #example using IPO isotope identification resultPP <- optimizeXcmsSet(mzmlfile, paramsPP, subdir="mtbls2") #example using CAMERA isotope identification resultPP <- optimizeXcmsSet(mzmlfile, paramsPP, isotopeIdentification="CAMERA", subdir="mtbls2", ppm=15, maxcharge=2)
This function converts an array into a matrix. This is useful to counter the implicit casting of matrices into arrays when only one row is selected. If a matrix is passed to the function, the matrix is returned, if an array is passed, a matrix with one row is returned.
toMatrix(data)
toMatrix(data)
data |
An array or a matrix |
A matrix
Gunnar Libsieller
data <- matrix(1:9, nrow=3) colnames(data) <- c("a","b","c") x <- data[1,] is.matrix(x) x <- toMatrix(x) is.matrix(x)
data <- matrix(1:9, nrow=3) colnames(data) <- c("a","b","c") x <- data[1,] is.matrix(x) x <- toMatrix(x) is.matrix(x)
This method takes a list of parameters and returns a list consisting of another two lists; one holding parameters ment for optimization and one holding fixed parameters.
typeCastParams(params)
typeCastParams(params)
params |
A list of parameters for an xcms-method |
This method takes a list of parameters and returns a list consisting of another two lists; one holding parameters ment for optimization and one holding fixed parameters.
A list of:
to_optimize |
A LIST containing all parameters which should be optimized. |
no_optimization |
A LIST containing all parameters which should not be optimized. |
Gunnar Libiseller
optimizeXcmsSet
, optimizeRetGroup
params <- getDefaultXcmsSetStartingParams() typ_params <- typeCastParams(params)
params <- getDefaultXcmsSetStartingParams() typ_params <- typeCastParams(params)
This function writes findPeaks, retcor and grouping parameters to a file using write.table.
writeParamsTable(peakPickingSettings, retCorGroupSettings, file, ...)
writeParamsTable(peakPickingSettings, retCorGroupSettings, file, ...)
peakPickingSettings |
A list of optimized settings for xcms-methods findPeaks.centWave or findPeaks.matchedFilter |
retCorGroupSettings |
A list of optimized settings for xcms-methods for retcor.obiwarp and group.density |
file |
The name of the outputfile for the parameters. |
... |
Additional parameters for write.table. |
This function writes findPeaks, retcor and grouping parameters to a file using write.table.
none
Gunnar Libiseller
#creating list of peak picking parameters paramsPP <- list(min_peakwidth=5, max_peakwidth=12, ppm=58, mzdiff=-0.001, snthresh=10, noise=0, prefilter=3, value_of_prefilter=100, mzCenterFun="wMean", integrate=1, fitgauss=FALSE, verbose.columns=FALSE, nSlaves=1) #creating list of retention time correction and grouping parameters paramsRTCGroup <- list(retcorMethod="obiwarp", distFunc="cor", gapInit=0.2, gapExtend=2.4, profStep=1, plottype="none", response=1, factorDiag=2, factorGap=1, localAlignment=0, initPenalty=0, bw=30, minfrac=0.5, minsamp=1, mzwid=0.25, max=50) #writing parameters to the file "params.tsv" writeParamsTable(paramsPP, paramsRTCGroup, "params.tsv")
#creating list of peak picking parameters paramsPP <- list(min_peakwidth=5, max_peakwidth=12, ppm=58, mzdiff=-0.001, snthresh=10, noise=0, prefilter=3, value_of_prefilter=100, mzCenterFun="wMean", integrate=1, fitgauss=FALSE, verbose.columns=FALSE, nSlaves=1) #creating list of retention time correction and grouping parameters paramsRTCGroup <- list(retcorMethod="obiwarp", distFunc="cor", gapInit=0.2, gapExtend=2.4, profStep=1, plottype="none", response=1, factorDiag=2, factorGap=1, localAlignment=0, initPenalty=0, bw=30, minfrac=0.5, minsamp=1, mzwid=0.25, max=50) #writing parameters to the file "params.tsv" writeParamsTable(paramsPP, paramsRTCGroup, "params.tsv")
This function prints a script of the optimized findPeaks, retcor and grouping parameters to the screen.
writeRScript(peakPickingSettings, retCorGroupSettings, nSlaves = 0)
writeRScript(peakPickingSettings, retCorGroupSettings, nSlaves = 0)
peakPickingSettings |
The optimized settings for xcms-methods findPeaks.centWave or findPeaks.matchedFilter |
retCorGroupSettings |
The optimized settings for xcms-methods for retcor.obiwarp and group.density |
nSlaves |
DEPRECATED |
This function prints a script out of the optimized findPeaks, retcor and grouping parameters to the screen.
The function message
is used to print the script. For capuring the
output capture.output(writeRScript(...), type = "message")
might
be used.
none
Gunnar Libiseller, Thomas Riebenbauer ([email protected])
#creating list of peak picking parameters paramsPP <- list(min_peakwidth=5, max_peakwidth=12, ppm=58, mzdiff=-0.001, snthresh=10, noise=0, prefilter=3, value_of_prefilter=100, mzCenterFun="wMean", integrate=1, fitgauss=FALSE, verbose.columns=FALSE, nSlaves=1) #creating list of retention time correction and grouping parameters paramsRTCGroup <- list(retcorMethod="obiwarp", distFunc="cor", gapInit=0.2, gapExtend=2.4, profStep=1, plottype="none", response=1, factorDiag=2, factorGap=1, localAlignment=0, initPenalty=0, bw=30, minfrac=0.5, minsamp=1, mzwid=0.25, max=50) #outputting an xcms-script to the display writeRScript(paramsPP, paramsRTCGroup, 4)
#creating list of peak picking parameters paramsPP <- list(min_peakwidth=5, max_peakwidth=12, ppm=58, mzdiff=-0.001, snthresh=10, noise=0, prefilter=3, value_of_prefilter=100, mzCenterFun="wMean", integrate=1, fitgauss=FALSE, verbose.columns=FALSE, nSlaves=1) #creating list of retention time correction and grouping parameters paramsRTCGroup <- list(retcorMethod="obiwarp", distFunc="cor", gapInit=0.2, gapExtend=2.4, profStep=1, plottype="none", response=1, factorDiag=2, factorGap=1, localAlignment=0, initPenalty=0, bw=30, minfrac=0.5, minsamp=1, mzwid=0.25, max=50) #outputting an xcms-script to the display writeRScript(paramsPP, paramsRTCGroup, 4)