adjustRtime-peakGroups: Retention time correction based on alignment of house keeping...

Description Usage Arguments Value Slots Note Author(s) References See Also Examples

Description

This method performs retention time adjustment based on the alignment of chromatographic peak groups present in all/most samples (hence corresponding to house keeping compounds). First the retention time deviation of these peak groups is described by fitting either a polynomial (smooth = "loess") or a linear ( smooth = "linear") model to the data points. These models are subsequently used to adjust the retention time of each spectrum in each sample.

It is also possible to exclude certain samples within an experiment from the estimation of the alignment models. The parameter subset allows to define the indices of samples within object that should be aligned. Samples not part of this subset are left out in the estimation of the alignment models, but their retention times are subsequently adjusted based on the alignment results of the closest sample in subset (close in terms of position within the object). Alignment could thus be performed on only real samples leaving out e.g. blanks, which are then in turn adjusted based on the closest real sample. Here it is up to the user to ensure that the samples within object are ordered correctly (e.g. by injection index).

How the non-subset samples are adjusted bases also on the parameter subsetAdjust: with subsetAdjust = "previous", each non-subset sample is adjusted based on the closest previous subset sample which results in most cases with adjusted retention times of the non-subset sample being identical to the subset sample on which the adjustment bases. The second, default, option is to use subsetAdjust = "average" in which case each non subset sample is adjusted based on the average retention time adjustment from the previous and following subset sample. For the average a weighted mean is used with weights being the inverse of the distance of the non-subset sample to the subset samples used for alignment.

See also section Alignment of experiments including blanks in the xcms vignette for an example.

The PeakGroupsParam class allows to specify all settings for the retention time adjustment based on house keeping peak groups present in most samples. Instances should be created with the PeakGroupsParam constructor.

adjustRtimePeakGroups returns the features (peak groups) which would, depending on the provided PeakGroupsParam, be selected for alignment/retention time correction.

minFraction,minFraction<-: getter and setter for the minFraction slot of the object.

extraPeaks,extraPeaks<-: getter and setter for the extraPeaks slot of the object.

smooth,smooth<-: getter and setter for the smooth slot of the object.

span,span<-: getter and setter for the span slot of the object.

family,family<-: getter and setter for the family slot of the object.

peakGroupsMatrix,peakGroupsMatrix<-: getter and setter for the peakGroupsMatrix slot of the object.

subset,subset<-: getter and setter for the subset slot of the object.

subsetAdjust,subsetAdjust<-: getter and setter for the subsetAdjust slot of the object.

adjustRtime,XCMSnExp,PeakGroupsParam: performs retention time correction based on the alignment of peak groups (features) found in all/most samples. The correction function identified on these peak groups is applied to the retention time of all spectra in the object, i.e. retention times of all spectra, also MS level > 1 are adjusted.

Usage

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
PeakGroupsParam(minFraction = 0.9, extraPeaks = 1, smooth = "loess",
  span = 0.2, family = "gaussian", peakGroupsMatrix = matrix(nrow =
  0, ncol = 0), subset = integer(), subsetAdjust = c("average",
  "previous"))

adjustRtimePeakGroups(object, param = PeakGroupsParam())

## S4 method for signature 'PeakGroupsParam'
show(object)

## S4 method for signature 'PeakGroupsParam'
minFraction(object)

## S4 replacement method for signature 'PeakGroupsParam'
minFraction(object) <- value

## S4 method for signature 'PeakGroupsParam'
extraPeaks(object)

## S4 replacement method for signature 'PeakGroupsParam'
extraPeaks(object) <- value

## S4 method for signature 'PeakGroupsParam'
smooth(x)

## S4 replacement method for signature 'PeakGroupsParam'
smooth(object) <- value

## S4 method for signature 'PeakGroupsParam'
span(object)

## S4 replacement method for signature 'PeakGroupsParam'
span(object) <- value

## S4 method for signature 'PeakGroupsParam'
family(object)

## S4 replacement method for signature 'PeakGroupsParam'
family(object) <- value

## S4 method for signature 'PeakGroupsParam'
peakGroupsMatrix(object)

## S4 replacement method for signature 'PeakGroupsParam'
peakGroupsMatrix(object) <- value

## S4 method for signature 'PeakGroupsParam'
subset(x)

## S4 replacement method for signature 'PeakGroupsParam'
subset(object) <- value

## S4 method for signature 'PeakGroupsParam'
subsetAdjust(object)

## S4 replacement method for signature 'PeakGroupsParam'
subsetAdjust(object) <- value

## S4 method for signature 'XCMSnExp,PeakGroupsParam'
adjustRtime(object, param)

Arguments

minFraction

numeric(1) between 0 and 1 defining the minimum required fraction of samples in which peaks for the peak group were identified. Peak groups passing this criteria will aligned across samples and retention times of individual spectra will be adjusted based on this alignment. For minFraction = 1 the peak group has to contain peaks in all samples of the experiment. Note that if subset is provided, the specified fraction is relative to the defined subset of samples and not to the total number of samples within the experiment (i.e. a peak has to be present in the specified proportion of subset samples).

extraPeaks

numeric(1) defining the maximal number of additional peaks for all samples to be assigned to a peak group (i.e. feature) for retention time correction. For a data set with 6 samples, extraPeaks = 1 uses all peak groups with a total peak count <= 6 + 1. The total peak count is the total number of peaks being assigned to a peak group and considers also multiple peaks within a sample being assigned to the group.

smooth

character defining the function to be used, to interpolate corrected retention times for all peak groups. Either "loess" or "linear".

span

numeric(1) defining the degree of smoothing (if smooth = "loess"). This parameter is passed to the internal call to loess.

family

character defining the method to be used for loess smoothing. Allowed values are "gaussian" and "symmetric".See loess for more information.

peakGroupsMatrix

optional matrix of (raw) retention times for the peak groups on which the alignment should be performed. Each column represents a sample, each row a feature/peak group. Such a matrix is for example returned by the adjustRtimePeakGroups method.

subset

integer with the indices of samples within the experiment on which the alignment models should be estimated. Samples not part of the subset are adjusted based on the closest subset sample. See description above for more details.

subsetAdjust

character specifying the method with which non-subset samples should be adjusted. Supported options are "previous" and "average" (default). See description above for more information.

object

For adjustRtime: an XCMSnExp object containing the results from a previous chromatographic peak detection (see findChromPeaks) and alignment analysis (see groupChromPeaks).

For all other methods: a PeakGroupsParam object.

param

A PeakGroupsParam object containing all settings for the retention time correction method..

value

The value for the slot.

x

a PeakGroupsParam object.

Value

The PeakGroupsParam function returns a PeakGroupsParam class instance with all of the settings specified for retention time adjustment based on house keeping features/peak groups.

For adjustRtimePeakGroups: a matrix, rows being features, columns samples, of retention times. The features are ordered by the median retention time across columns.

For adjustRtime: a XCMSnExp object with the results of the retention time adjustment step. These can be accessed with the adjustedRtime method. Retention time correction does also adjust the retention time of the identified chromatographic peaks (accessed via chromPeaks. Note that retention time correction drops all previous alignment results from the result object.

Slots

.__classVersion__,minFraction,extraPeaks,smooth,span,family,peakGroupsMatrix,subset,subsetAdjust

See corresponding parameter above. .__classVersion__ stores the version from the class. Slots values should exclusively be accessed via the corresponding getter and setter methods listed above.

Note

These methods and classes are part of the updated and modernized xcms user interface which will eventually replace the group methods. All of the settings to the alignment algorithm can be passed with a PeakGroupsParam object.

The matrix with the (raw) retention times of the peak groups used in the alignment is added to the peakGroupsMatrix slot of the PeakGroupsParam object that is stored into the corresponding process history step (see processHistory for how to access the process history).

adjustRtimePeakGroups is supposed to be called before the sample alignment, but after a correspondence (peak grouping).

This method requires that a correspondence analysis has been performed on the data, i.e. that grouped chromatographic peaks/features are present (see groupChromPeaks for details).

Calling adjustRtime on an XCMSnExp object will cause all peak grouping (correspondence) results and any previous retention time adjustments to be dropped. In some instances, the adjustRtime,XCMSnExp,PeakGroupsParam re-adjusts adjusted retention times to ensure them being in the same order than the raw (original) retention times.

Author(s)

Colin Smith, Johannes Rainer

References

Colin A. Smith, Elizabeth J. Want, Grace O'Maille, Ruben Abagyan and Gary Siuzdak. "XCMS: Processing Mass Spectrometry Data for Metabolite Profiling Using Nonlinear Peak Alignment, Matching, and Identification" Anal. Chem. 2006, 78:779-787.

See Also

The do_adjustRtime_peakGroups core API function and retcor.peakgroups for the old user interface. plotAdjustedRtime for visualization of alignment results.

XCMSnExp for the object containing the results of the alignment.

Other retention time correction methods: adjustRtime-obiwarp, adjustRtime

Examples

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
##############################
## Chromatographic peak detection and grouping.
##
## Below we perform first a peak detection (using the matchedFilter
## method) on some of the test files from the faahKO package followed by
## a peak grouping.
library(faahKO)
library(xcms)
fls <- dir(system.file("cdf/KO", package = "faahKO"), recursive = TRUE,
           full.names = TRUE)

## Reading 2 of the KO samples
raw_data <- readMSData(fls[1:2], mode = "onDisk")

## Perform the peak detection using the matchedFilter method.
mfp <- MatchedFilterParam(snthresh = 20, binSize = 1)
res <- findChromPeaks(raw_data, param = mfp)

head(chromPeaks(res))
## The number of peaks identified per sample:
table(chromPeaks(res)[, "sample"])

## Performing the peak grouping using the "peak density" method.
p <- PeakDensityParam(sampleGroups = c(1, 1))
res <- groupChromPeaks(res, param = p)

## Perform the retention time adjustment using peak groups found in both
## files.
fgp <- PeakGroupsParam(minFraction = 1)

## Before running the alignment we can evaluate which features (peak groups)
## would be used based on the specified parameters.
pkGrps <- adjustRtimePeakGroups(res, param = fgp)

## We can also plot these to evaluate if the peak groups span a large portion
## of the retention time range.
plot(x = pkGrps[, 1], y = rep(1, nrow(pkGrps)), xlim = range(rtime(res)),
    ylim = c(1, 2), xlab = "rt", ylab = "", yaxt = "n")
points(x = pkGrps[, 2], y = rep(2, nrow(pkGrps)))
segments(x0 = pkGrps[, 1], x1 = pkGrps[, 2],
    y0 = rep(1, nrow(pkGrps)), y1 = rep(2, nrow(pkGrps)))
grid()
axis(side = 2, at = c(1, 2), labels = colnames(pkGrps))

## Next we perform the alignment.
res <- adjustRtime(res, param = fgp)

## Any grouping information was dropped
hasFeatures(res)

## Plot the raw against the adjusted retention times.
plot(rtime(raw_data), rtime(res), pch = 16, cex = 0.25, col = fromFile(res))

## Adjusterd retention times can be accessed using
## rtime(object, adjusted = TRUE) and adjustedRtime
all.equal(rtime(res), adjustedRtime(res))

## To get the raw, unadjusted retention times:
all.equal(rtime(res, adjusted = FALSE), rtime(raw_data))

## To extract the retention times grouped by sample/file:
rts <- rtime(res, bySample = TRUE)

anupbharade09/xcms_test documentation built on May 14, 2019, 4:07 a.m.