| XCMSnExp-class {xcms} | R Documentation |
The XCMSnExp object is designed to contain all results
from metabolomics data preprocessing (chromatographic peak detection,
peak grouping (correspondence) and retention time correction). The
corresponding elements in the msFeatureData slot are
"chromPeaks" (a matrix), "featureDefinitions"
(a DataFrame) and "adjustedRtime" (a list of
numeric vectors). Note that these should not be accessed directly but
rather via their accessor methods.
Along with the results, the object contains the processing history that
allows to track each processing step along with the used settings. This
can be extracted with the processHistory method.
The XCMSnExp object directly extends the
OnDiskMSnExp object and provides thus an easy
access to the full raw data at any stage of an analysis.
Objects from this class should not be created directly, they are
returned as result from the findChromPeaks method.
XCMSnExp objects can be coerced into xcmsSet
objects using the as method.
processHistoryTypes returns the available types of
process histories. These can be passed with argument type to the
processHistory method to extract specific process step(s).
profMat: creates a profile matrix, which
is a n x m matrix, n (rows) representing equally spaced m/z values (bins)
and m (columns) the retention time of the corresponding scans. Each cell
contains the maximum intensity measured for the specific scan and m/z
values. See profMat for more details and description of
the various binning methods.
hasAdjustedRtime: whether the object provides adjusted
retention times.
hasFeatures: whether the object contains correspondence
results (i.e. features).
hasChromPeaks: whether the object contains peak
detection results.
adjustedRtime,adjustedRtime<-:
extract/set adjusted retention times. adjustedRtime<- should not
be called manually, it is called internally by the
adjustRtime methods. For XCMSnExp objects,
adjustedRtime<- does also apply retention time adjustments to
eventually present chromatographic peaks. The bySample parameter
allows to specify whether the adjusted retention time should be grouped
by sample (file).
featureDefinitions, featureDefinitions<-: extract
or set the correspondence results, i.e. the mz-rt features (peak groups).
Similar to the chromPeaks it is possible to extract features for
specified m/z and/or rt ranges (see chromPeaks for more details).
chromPeaks, chromPeaks<-: extract or set
the matrix containing the information on identified chromatographic
peaks. Parameter bySample allows to specify whether peaks should
be returned ungrouped (default bySample = FALSE) or grouped by
sample (bySample = TRUE). The chromPeaks<- method for
XCMSnExp objects removes also all correspondence (peak grouping)
and retention time correction (alignment) results. The optional
arguments rt, mz and ppm allow to extract only
chromatographic peaks overlapping (if type = "any") or completely
within (if type = "within") the defined retention time and mz
ranges.
See description of the return value for details on the returned matrix.
Users usually don't have to use the chromPeaks<- method directly
as detected chromatographic peaks are added to the object by the
findChromPeaks method.
rtime: extracts the retention time for each
scan. The bySample parameter allows to return the values grouped
by sample/file and adjusted whether adjusted or raw retention
times should be returned. By default the method returns adjusted
retention times, if they are available (i.e. if retention times were
adjusted using the adjustRtime method).
mz: extracts the mz values from each scan of
all files within an XCMSnExp object. These values are extracted
from the original data files and eventual processing steps are applied
on the fly. Using the bySample parameter it is possible to
switch from the default grouping of mz values by spectrum/scan to a
grouping by sample/file.
intensity: extracts the intensity values from
each scan of all files within an XCMSnExp object. These values are
extracted from the original data files and eventual processing steps are
applied on the fly. Using the bySample parameter it is
possible to switch from the default grouping of intensity values by
spectrum/scan to a grouping by sample/file.
spectra: extracts the
Spectrum objects containing all data from
object. The values are extracted from the original data files and
eventual processing steps are applied on the fly. By setting
bySample = TRUE, the spectra are returned grouped by sample/file.
If the XCMSnExp object contains adjusted retention times, these
are returned by default in the Spectrum objects (can be
overwritten by setting adjusted = FALSE).
processHistory: returns a list with
ProcessHistory objects (or objects inheriting from this
base class) representing the individual processing steps that have been
performed, eventually along with their settings (Param parameter
class). Optional arguments fileIndex, type and
msLevel allow to restrict to process steps of a certain type or
performed on a certain file or MS level.
dropChromPeaks: drops any identified chromatographic
peaks and returns the object without that information. Note that for
XCMSnExp objects the method drops by default also results from a
correspondence (peak grouping) analysis. Adjusted retention times are
removed if the alignment has been performed after peak detection.
This can be overruled with keepAdjustedRtime = TRUE.
dropFeatureDefinitions: drops the results from a
correspondence (peak grouping) analysis, i.e. the definition of the mz-rt
features and returns the object without that information. Note that for
XCMSnExp objects the method will also by default drop retention
time adjustment results, if these were performed after the last peak
grouping (i.e. which base on the results from the peak grouping that are
going to be removed). All related process history steps are
removed too as well as eventually filled in peaks
(by fillChromPeaks). The parameter keepAdjustedRtime
can be used to avoid removal of adjusted retention times.
dropAdjustedRtime: drops any retention time
adjustment information and returns the object without adjusted retention
time. For XCMSnExp objects, this also reverts the retention times
reported for the chromatographic peaks in the peak matrix to the
original, raw, ones (after chromatographic peak detection). Note that
for XCMSnExp objects the method drops also all peak grouping
results if these were performed after the retention time
adjustment. All related process history steps are removed too.
findChromPeaks performs chromatographic peak detection
on the provided XCMSnExp objects. For more details see the method
for XCMSnExp. Note that the findChromPeaks
method for XCMSnExp objects removes previously identified
chromatographic peaks and aligned features. Previous alignment (retention
time adjustment) results are kept, i.e. chromatographic peak detection
is performed using adjusted retention times if the data was first
aligned using e.g. obiwarp (adjustRtime-obiwarp).
dropFilledChromPeaks: drops any filled-in chromatographic
peaks (filled in by the fillChromPeaks method) and all
related process history steps.
description spectrapply applies the provided function to each
Spectrum in the object and returns its
results. If no function is specified the function simply returns the
list of Spectrum objects.
processHistoryTypes() ## S4 method for signature 'OnDiskMSnExp' profMat(object, method = "bin", step = 0.1, baselevel = NULL, basespace = NULL, mzrange. = NULL, fileIndex, ...) ## S4 method for signature 'XCMSnExp' show(object) ## S4 method for signature 'XCMSnExp' hasAdjustedRtime(object) ## S4 method for signature 'XCMSnExp' hasFeatures(object) ## S4 method for signature 'XCMSnExp' hasChromPeaks(object) ## S4 method for signature 'XCMSnExp' adjustedRtime(object, bySample = FALSE) ## S4 replacement method for signature 'XCMSnExp' adjustedRtime(object) <- value ## S4 method for signature 'XCMSnExp' featureDefinitions(object, mz = numeric(), rt = numeric(), ppm = 0, type = "any") ## S4 replacement method for signature 'XCMSnExp' featureDefinitions(object) <- value ## S4 method for signature 'XCMSnExp' chromPeaks(object, bySample = FALSE, rt = numeric(), mz = numeric(), ppm = 0, type = "any") ## S4 replacement method for signature 'XCMSnExp' chromPeaks(object) <- value ## S4 method for signature 'XCMSnExp' rtime(object, bySample = FALSE, adjusted = hasAdjustedRtime(object)) ## S4 method for signature 'XCMSnExp' mz(object, bySample = FALSE, BPPARAM = bpparam()) ## S4 method for signature 'XCMSnExp' intensity(object, bySample = FALSE, BPPARAM = bpparam()) ## S4 method for signature 'XCMSnExp' spectra(object, bySample = FALSE, adjusted = hasAdjustedRtime(object), BPPARAM = bpparam()) ## S4 method for signature 'XCMSnExp' processHistory(object, fileIndex, type, msLevel) ## S4 method for signature 'XCMSnExp' dropChromPeaks(object, keepAdjustedRtime = FALSE) ## S4 method for signature 'XCMSnExp' dropFeatureDefinitions(object, keepAdjustedRtime = FALSE, dropLastN = -1) ## S4 method for signature 'XCMSnExp' dropAdjustedRtime(object) ## S4 method for signature 'XCMSnExp' profMat(object, method = "bin", step = 0.1, baselevel = NULL, basespace = NULL, mzrange. = NULL, fileIndex, ...) ## S4 method for signature 'XCMSnExp,Param' findChromPeaks(object, param, BPPARAM = bpparam(), return.type = "XCMSnExp", msLevel = 1L) ## S4 method for signature 'XCMSnExp' dropFilledChromPeaks(object) ## S4 method for signature 'XCMSnExp' spectrapply(object, FUN = NULL, BPPARAM = bpparam(), ...)
object |
For |
method |
The profile matrix generation method. Allowed are |
step |
numeric(1) representing the m/z bin size. |
baselevel |
numeric(1) representing the base value to which
empty elements (i.e. m/z bins without a measured intensity) should be
set. Only considered if |
basespace |
numeric(1) representing the m/z length after
which the signal will drop to the base level. Linear interpolation will
be used between consecutive data points falling within
|
mzrange. |
Optional numeric(2) manually specifying the mz value range to be used for binnind. If not provided, the whole mz value range is used. |
fileIndex |
For |
... |
Additional parameters. |
bySample |
logical(1) specifying whether results should be grouped by sample. |
value |
For For For |
mz |
optional |
rt |
optional |
ppm |
optional |
type |
For |
adjusted |
logical(1) whether adjusted or raw (i.e. the original retention times reported in the files) should be returned. |
BPPARAM |
Parameter class for parallel processing. See
|
msLevel |
|
keepAdjustedRtime |
For |
dropLastN |
For |
param |
A |
return.type |
Character specifying what type of object the method should
return. Can be either |
FUN |
For |
For profMat: a list with a the profile matrix
matrix (or matrices if fileIndex was not specified or if
length(fileIndex) > 1). See profile-matrix for
general help and information about the profile matrix.
For adjustedRtime: if bySample = FALSE a numeric
vector with the adjusted retention for each spectrum of all files/samples
within the object. If bySample = TRUE a list (length equal
to the number of samples) with adjusted retention times grouped by
sample. Returns NULL if no adjusted retention times are present.
For featureDefinitions: a DataFrame with peak grouping
information, each row corresponding to one mz-rt feature (grouped peaks
within and across samples) and columns "mzmed" (median mz value),
"mzmin" (minimal mz value), "mzmax" (maximum mz value),
"rtmed" (median retention time), "rtmin" (minimal retention
time), "rtmax" (maximal retention time) and "peakidx".
Column "peakidx" contains a list with indices of
chromatographic peaks (rows) in the matrix returned by the
chromPeaks method that belong to that feature group. The method
returns NULL if no feature definitions are present.
For chromPeaks: if bySample = FALSE a matrix
with at least the following columns:
"mz" (intensity-weighted mean of mz values of the peak across
scans/retention times),
"mzmin" (minimal mz value),
"mzmax" (maximal mz value),
"rt" (retention time for the peak apex),
"rtmin" (minimal retention time),
"rtmax" (maximal retention time),
"into" (integrated, original, intensity of the peak),
"maxo" (maximum intentity of the peak),
"sample" (sample index in which the peak was identified) and
"is_filled" defining whether the chromatographic peak was
identified by the peak picking algorithm (0) or was added by the
fillChromPeaks method (1).
Depending on the employed peak detection algorithm and the
verboseColumns parameter of it additional columns might be
returned. For bySample = TRUE the chronatographic peaks are
returned as a list of matrices, each containing the
chromatographic peaks of a specific sample. For samples in which no
peaks were detected a matrix with 0 rows is returned.
For rtime: if bySample = FALSE a numeric vector with
the retention times of each scan, if bySample = TRUE a
list of numeric vectors with the retention times per sample.
For mz: if bySample = FALSE a list with the mz
values (numeric vectors) of each scan. If bySample = TRUE a
list with the mz values per sample.
For intensity: if bySample = FALSE a list with
the intensity values (numeric vectors) of each scan. If
bySample = TRUE a list with the intensity values per
sample.
For spectra: if bySample = FALSE a list with
Spectrum objects. If bySample = TRUE the
result is grouped by sample, i.e. as a list of lists, each
element in the outer list being the list of spectra
of the specific file.
For processHistory: a list of
ProcessHistory objects providing the details of the
individual data processing steps that have been performed.
.processHistorylist with XProcessHistory objects
tracking all individual analysis steps that have been performed.
msFeatureDataMsFeatureData class extending environment
and containing the results from a chromatographic peak detection (element
"chromPeaks"), peak grouping (element "featureDefinitions")
and retention time correction (element "adjustedRtime") steps.
This object should not be manipulated directly.
The "chromPeaks" element in the msFeatureData slot is
equivalent to the @peaks slot of the xcmsSet object, the
"featureDefinitions" contains information from the @groups
and @groupidx slots from an xcmsSet object.
Johannes Rainer
xcmsSet for the old implementation.
OnDiskMSnExp, MSnExp
and pSet for a complete list of inherited methods.
findChromPeaks for available peak detection methods
returning a XCMSnExp object as a result.
groupChromPeaks for available peak grouping
methods and featureDefinitions for the method to extract
the feature definitions representing the peak grouping results.
adjustRtime for retention time adjustment methods.
chromatogram to extract MS data as
Chromatogram objects.
extractMsData for the method to extract MS data as
data.frames.
fillChromPeaks for the method to fill-in eventually
missing chromatographic peaks for a feature in some samples.
## Loading the data from 2 files of the faahKO package.
library(faahKO)
od <- readMSData(c(system.file("cdf/KO/ko15.CDF", package = "faahKO"),
system.file("cdf/KO/ko16.CDF", package = "faahKO")),
mode = "onDisk")
## Now we perform a chromatographic peak detection on this data set using the
## matched filter method. We are tuning the settings such that it performs
## faster.
mfp <- MatchedFilterParam(binSize = 6)
xod <- findChromPeaks(od, param = mfp)
## The results from the peak detection are now stored in the XCMSnExp
## object
xod
## The detected peaks can be accessed with the chromPeaks method.
head(chromPeaks(xod))
## The settings of the chromatographic peak detection can be accessed with
## the processHistory method
processHistory(xod)
## Also the parameter class for the peak detection can be accessed
processParam(processHistory(xod)[[1]])
## The XCMSnExp inherits all methods from the pSet and OnDiskMSnExp classes
## defined in Bioconductor's MSnbase package. To access the (raw) retention
## time for each spectrum we can use the rtime method. Setting bySample = TRUE
## would cause the retention times to be grouped by sample
head(rtime(xod))
## Similarly it is possible to extract the mz values or the intensity values
## using the mz and intensity method, respectively, also with the option to
## return the results grouped by sample instead of the default, which is
## grouped by spectrum. Finally, to extract all of the data we can use the
## spectra method which returns Spectrum objects containing all raw data.
## Note that all these methods read the information from the original input
## files and subsequently apply eventual data processing steps to them.
mzs <- mz(xod, bySample = TRUE)
length(mzs)
lengths(mzs)
## The full data could also be read using the spectra data, which returns
## a list of Spectrum object containing the mz, intensity and rt values.
## spctr <- spectra(xod)
## To get all spectra of the first file we can split them by file
## head(split(spctr, fromFile(xod))[[1]])
############
## Filtering
##
## XCMSnExp objects can be filtered by file, retention time, mz values or
## MS level. For some of these filter preprocessing results (mostly
## retention time correction and peak grouping results) will be dropped.
## Below we filter the XCMSnExp object by file to extract the results for
## only the second file.
xod_2 <- filterFile(xod, file = 2)
xod_2
## Now the objects contains only the idenfified peaks for the second file
head(chromPeaks(xod_2))
head(chromPeaks(xod)[chromPeaks(xod)[, "sample"] == 2, ])
##########
## Coercing to an xcmsSet object
##
## We can also coerce the XCMSnExp object into an xcmsSet object:
xs <- as(xod, "xcmsSet")
head(peaks(xs))