Skip to contents

The XcmsExperiment is a data container for xcms preprocessing results (i.e. results from chromatographic peak detection, alignment and correspondence analysis).

It provides the same functionality than the XCMSnExp object, but uses the more advanced and modern MS infrastructure provided by the MsExperiment and Spectra Bioconductor packages. With this comes a higher flexibility on how and where to store the data.

Documentation of the various functions for XcmsExperiment objects are grouped by topic and provided in the sections below.

The default xcms workflow is to perform

  • chromatographic peak detection using findChromPeaks()

  • optionally refine identified chromatographic peaks using refineChromPeaks()

  • perform an alignment (retention time adjustment) using adjustRtime(). Depending on the method used this requires to run a correspondence analysis first

  • perform a correspondence analysis using the groupChromPeaks() function to group chromatographic peaks across samples to define the LC-MS features.

  • optionally perform a gap-filling to rescue signal in samples in which no chromatographic peak was identified and hence a missing value would be reported. This can be performed using the fillChromPeaks() function.

Usage

filterFeatureDefinitions(object, ...)

# S4 method for class 'MsExperiment'
filterRt(object, rt = numeric(), ...)

# S4 method for class 'MsExperiment'
filterMzRange(object, mz = numeric(), msLevel. = uniqueMsLevels(object))

# S4 method for class 'MsExperiment'
filterMz(object, mz = numeric(), msLevel. = uniqueMsLevels(object))

# S4 method for class 'MsExperiment'
filterMsLevel(object, msLevel. = uniqueMsLevels(object))

# S4 method for class 'MsExperiment'
uniqueMsLevels(object)

# S4 method for class 'MsExperiment'
filterFile(object, file = integer(), ...)

# S4 method for class 'MsExperiment'
rtime(object)

# S4 method for class 'MsExperiment'
fromFile(object)

# S4 method for class 'MsExperiment'
fileNames(object)

# S4 method for class 'MsExperiment'
polarity(object)

# S4 method for class 'MsExperiment'
filterIsolationWindow(object, mz = numeric())

# S4 method for class 'MsExperiment'
chromatogram(
  object,
  rt = matrix(nrow = 0, ncol = 2),
  mz = matrix(nrow = 0, ncol = 2),
  aggregationFun = "sum",
  msLevel = 1L,
  isolationWindowTargetMz = NULL,
  chunkSize = 2L,
  return.type = "MChromatograms",
  BPPARAM = bpparam()
)

# S4 method for class 'XcmsExperiment,ANY,ANY,ANY'
[(x, i, j, ..., drop = TRUE)

# S4 method for class 'XcmsExperiment'
filterIsolationWindow(object, mz = numeric())

# S4 method for class 'XcmsExperiment'
filterRt(object, rt, msLevel.)

# S4 method for class 'XcmsExperiment'
filterMzRange(object, mz = numeric(), msLevel. = uniqueMsLevels(object))

# S4 method for class 'XcmsExperiment'
filterMsLevel(object, msLevel. = uniqueMsLevels(object))

# S4 method for class 'XcmsExperiment'
hasChromPeaks(object, msLevel = integer())

# S4 method for class 'XcmsExperiment'
dropChromPeaks(object, keepAdjustedRtime = FALSE)

# S4 method for class 'XcmsExperiment'
chromPeaks(object) <- value

# S4 method for class 'XcmsExperiment'
chromPeaks(
  object,
  rt = numeric(),
  mz = numeric(),
  ppm = 0,
  msLevel = integer(),
  type = c("any", "within", "apex_within"),
  isFilledColumn = FALSE
)

# S4 method for class 'XcmsExperiment'
chromPeakData(object) <- value

# S4 method for class 'XcmsExperiment'
chromPeakData(
  object,
  msLevel = integer(),
  return.type = c("DataFrame", "data.frame")
)

# S4 method for class 'XcmsExperiment'
filterChromPeaks(
  object,
  keep = rep(TRUE, nrow(.chromPeaks(object))),
  method = "keep",
  ...
)

# S4 method for class 'XcmsExperiment'
dropAdjustedRtime(object)

# S4 method for class 'MsExperiment'
hasAdjustedRtime(object)

# S4 method for class 'XcmsExperiment'
rtime(object, adjusted = hasAdjustedRtime(object))

# S4 method for class 'XcmsExperiment'
adjustedRtime(object)

# S4 method for class 'XcmsExperiment'
hasFeatures(object, msLevel = integer())

# S4 method for class 'XcmsExperiment'
featureDefinitions(object) <- value

# S4 method for class 'XcmsExperiment'
featureDefinitions(
  object,
  mz = numeric(),
  rt = numeric(),
  ppm = 0,
  type = c("any", "within", "apex_within"),
  msLevel = integer()
)

# S4 method for class 'XcmsExperiment'
dropFeatureDefinitions(object, keepAdjustedRtime = FALSE)

# S4 method for class 'XcmsExperiment'
filterFeatureDefinitions(object, features = integer())

# S4 method for class 'XcmsExperiment'
hasFilledChromPeaks(object)

# S4 method for class 'XcmsExperiment'
dropFilledChromPeaks(object)

# S4 method for class 'XcmsExperiment'
quantify(object, ...)

# S4 method for class 'XcmsExperiment'
featureValues(
  object,
  method = c("medret", "maxint", "sum"),
  value = "into",
  intensity = "into",
  filled = TRUE,
  missing = NA_real_,
  msLevel = integer()
)

# S4 method for class 'XcmsExperiment'
chromatogram(
  object,
  rt = matrix(nrow = 0, ncol = 2),
  mz = matrix(nrow = 0, ncol = 2),
  aggregationFun = "sum",
  msLevel = 1L,
  chunkSize = 2L,
  isolationWindowTargetMz = NULL,
  return.type = c("XChromatograms", "MChromatograms"),
  include = character(),
  chromPeaks = c("apex_within", "any", "none"),
  BPPARAM = bpparam()
)

# S4 method for class 'XcmsExperiment'
processHistory(object, type)

# S4 method for class 'XcmsExperiment'
filterFile(
  object,
  file,
  keepAdjustedRtime = hasAdjustedRtime(object),
  keepFeatures = FALSE,
  ...
)

featureArea(
  object,
  mzmin = min,
  mzmax = max,
  rtmin = min,
  rtmax = max,
  features = character()
)

# S4 method for class 'MsExperiment,missing'
plot(x, y, msLevel = 1L, peakCol = "#ff000060", ...)

Arguments

object

An XcmsExperiment object.

...

Additional optional parameters. For quantify: any parameter for the featureValues call used to extract the feature value matrix.

rt

For chromPeaks and featureDefinitions: numeric(2) defining the retention time range for which chromatographic peaks or features should be returned. The full range is used by default. For chromatogram: two column numerical matrix with each row representing the lower and upper retention time window(s) for the chromatograms. If not provided the full retention time range is used.

mz

For chromPeaks and featureDefinitions: numeric(2) optionally defining the m/z range for which chromatographic peaks or feature definitions should be returned. The full m/z range is used by default. For chromatogram: two-column numerical matrix with each row representing m/z range that should be aggregated into a chromatogram. If not provided the full m/z range of the data will be used (and hence a total ion chromatogram will be returned if aggregationFun = "sum" is used). For filterIsolationWindow: numeric(1) defining the m/z that should be contained within the spectra's isolation window.

msLevel.

For filterRt: ignored. filterRt will always filter by retention times on all MS levels regardless of this parameter. For chromatogram: integer with the MS level from which the chromatogram(s) should be extracted. Has to be either of length 1 or length equal to the numer of rows of the parameters mz and rt defining the m/z and rt regions from which the chromatograms should be created. Defaults to msLevel = 1L. for filterMsLevel: integer defining the MS level(s) to which the data should be subset.

file

For filterFile: integer with the indices of the samples (files) to which the data should be subsetted.

aggregationFun

For chromatogram: character(1) defining the function that should be used to aggregate intensities for retention time (i.e. each spectrum) along the specified m/z range (parameter mz). Defaults to aggregationFun = "sum" and hence all intensities will be summed up. Alternatively, use aggregationFun = "max" to use the maximal intensity per m/z range to create a base peak chromatogram (BPC).

msLevel

integer defining the MS level (or multiple MS level if the function supports it).

isolationWindowTargetMz

For chromatogram: numeric (of length equal to the number of rows of rt and mz) with the isolation window target m/z of the MS2 spectra from which the chromatgrom should be generated. For MS1 data (msLevel = 1L, the default), this parameter is ignored. See examples on chromatogram below for further information.

chunkSize

For chromatogram: integer(1) defining the number of files from which the data should be loaded at a time into memory. Defaults to chunkSize = 2L.

return.type

For chromPeakData: character(1) defining the class of the returned object. Can be either "DataFrame" (the default) or "data.frame". For chromatogram: character(1) defining the type of the returned object. Currently only return.type = "MChromatograms" is supported.

BPPARAM

For chromatogram: parallel processing setup. Defaults to BPPARAM = bpparam(). See bpparam() for more information.

x

An XcmsExperiment object.

i

For [: integer or logical defining the samples/files to subset.

j

For [: not supported.

drop

For [: ignored.

keepAdjustedRtime

logical(1): whether adjusted retention times (if present) should be retained.

value

For featureValues: character(1) defining which value should be reported for each feature in each sample. Can be any column of the chromPeaks matrix or "index" if simply the index of the assigned peak should be returned. Defaults to value = "into" thus the integrated peak area is reported.

ppm

For chromPeaks and featureDefinitions: optional numeric(1) specifying the ppm by which the m/z range (defined by mz should be extended. For a value of ppm = 10, all peaks within mz[1] - ppm / 1e6 and mz[2] + ppm / 1e6 are returned.

type

For chromPeaks and featureDefinitions and only if either mz and rt are defined too: character(1): defining which peaks (or features) should be returned. For type = "any": returns all chromatographic peaks or features also only partially overlapping any of the provided ranges. For type = "within": returns only peaks or features completely within the region defined by mz and/or rt. For type = "apex_within": returns peaks or features for which the m/z and retention time of the peak's apex is within the region defined by mz and/or rt. For processHistory: restrict returned processing steps to specific types. Use processHistoryTypes() to list all supported values.

isFilledColumn

For chromPeaks: logical(1) whether a column "is_filled" should be included in the returned matrix with the information whether a peak was detected or only filled-in. Note that this information is also provided in the chromPeakData data frame.

keep

For filterChromPeaks: logical, integer or character specifying which chromatographic peaks to keep. If logical the length of keep needs to match the number of rows of chromPeaks. Alternatively, keep allows to specify the index (row) of peaks to keep or their ID (i.e. row name in chromPeaks).

method

For featureValues: character(1) specifying the method to resolve multi-peak mappings within the same sample (correspondence analysis can assign more than one chromatographic peak within a sample to the same feature, e.g. if they are close in retention time). Options: method = "medret": report the value for the chromatographic peak closest to the feature's median retention time. method = "maxint": report the value for the chromatographic peak with the largest signal (parameter intensity allows to select the column in chromPeaks that should be used for signal). method = "sum": sum the value for all chromatographic peaks in a sample assigned to the same feature. The default is method = "medret". For filterChromPeaks: currently only method = "keep" is supported.

adjusted

For rtime,XcmsExperiment: whether adjusted or raw retention times should be returned. The default is to return adjusted retention times, if available.

features

For filterFeatureDefinitions and featureArea: logical, integer or character defining the features to keep or from which to extract the feature area, respectively. See function description for more information.

intensity

For featureValues: character(1) specifying the name of the column in the chromPeaks(objects) matrix containing the intensity value of the peak that should be used for the conflict resolution if method = "maxint".

filled

For featureValues: logical(1) specifying whether values for filled-in peaks should be reported. For filled = TRUE (the default) filled peak values are returned, otherwise NA is reported for the respective features in the samples in which no peak was detected.

missing

For featureValues: default value for missing values. Allows to define the value that should be reported for a missing peak intensity. Defaults to missing = NA_real_.

include

For chromatogram: deprecated; use parameter chromPeaks instead.

chromPeaks

For chromatogram: character(1) defining which chromatographic peaks should be returned. Can be either chromPeaks = "apex_within" (default) to return all chromatographic peaks with the m/z and RT of their apex within the m/z and retention time window, chromPeaks = "any" for all chromatographic peaks that are overlapping with the m/z - retention time window or chromPeaks = "none" to not include any chromatographic peaks. See also parameter type below for additional information.

keepFeatures

for most subsetting functions ([, filterFile): logical(1): wheter eventually present feature definitions should be retained in the returned (filtered) object.

mzmin

For featureArea: function to calculate the "mzmin" of a feature based on the "mzmin" values of the individual chromatographic peaks assigned to that feature. Defaults to mzmin = min.

mzmax

For featureArea: function to calculate the "mzmax" of a feature based on the "mzmax" values of the individual chromatographic peaks assigned to that feature. Defaults to mzmax = max.

rtmin

For featureArea: function to calculate the "rtmin" of a feature based on the "rtmin" values of the individual chromatographic peaks assigned to that feature. Defaults to rtmin = min.

rtmax

For featureArea: function to calculate the "rtmax" of a feature based on the "rtmax" values of the individual chromatographic peaks assigned to that feature. Defaults to rtmax = max.

y

For plot: should not be defined as it is not supported.

peakCol

For plot: defines the border color of the rectangles indicating the identified chromatographic peaks. Only a single color is supported. Defaults to `peakCol = "#ff000060".

Subsetting and filtering

  • [: subset an XcmsExperiment by sample (parameter i). Subsetting will by default drop correspondence results (as subsetting by samples will obviously affect the feature definition) and alignment results (adjusted retention times) while identified chromatographic peaks (for the selected samples) will be retained. Which preprocessing results should be kept or dropped can also be configured with optional parameters keepChromPeaks (by default TRUE), keepAdjustedRtime (by default FALSE) and keepFeatures (by default FALSE).

  • filterChromPeaks: filter chromatographic peaks of an XcmsExperiment keeping only those specified with parameter keep. Returns the XcmsExperiment with the filtered data. Chromatographic peaks to retain can be specified either by providing their index in the chromPeaks matrix, their ID (rowname in chromPeaks) or with a logical vector with the same length than number of rows of chromPeaks. Assignment of chromatographic peaks are updated to eventually present feature definitions after filtering.

  • filterFeatureDefinitions: filter feature definitions of an XcmsExperiment keeping only those defined with parameter features, which can be a logical of length equal to the number of features, an integer with the index of the features in featureDefinitions(object) to keep or a character with the feature IDs (i.e. row names in featureDefinitions(object)).

  • filterFile: filter an XcmsExperiment (or MsExperiment) by file (sample). The index of the samples to which the data should be subsetted can be specified with parameter file. The sole purpose of this function is to provide backward compatibility with the MSnbase package. Wherever possible, the [ function should be used instead for any sample-based subsetting. Parameters keepChromPeaks, keepAdjustedRtime and keepChromPeaks can be passed using .... Note also that in contrast to [, filterFile does not support subsetting in arbitrary order.

  • filterIsolationWindow: filter the spectra within an MsExperiment or XcmsExperiment object keeping only those with an isolation window containing the specified m/z (i.e., keeping spectra with an "isolationWindowLowerMz" smaller than the user-provided mz and an "isolationWindowUpperMz" larger than mz). For an XcmsExperiment also all chromatographic peaks (and subsequently also features) are removed for which the range of their "isolationWindowLowerMz" and "isolationWindowUpperMz" (columns in chromPeakData) do not contain the user provided mz.

  • filterMsLevel: filter the data of the XcmsExperiment or MsExperiment to keep only data of the MS level(s) specified with parameter msLevel..

  • filterMz, filterMzRange: filter the spectra within an XcmsExperiment or MsExperiment to the specified m/z range (parameter mz). For XcmsExperiment also identified chromatographic peaks and features are filtered keeping only those that are within the specified m/z range (i.e. for which the m/z of the peak apex is within the m/z range). Parameter msLevels. allows to restrict the filtering to only specified MS levels. By default data from all MS levels are filtered.

  • filterRt: filter an XcmsExperiment keeping only data within the specified retention time range (parameter rt). This function will keep all preprocessing results present within the retention time range: all identified chromatographic peaks with the retention time of the apex position within the retention time range rt are retained along, if present, with the associated features. Parameter msLevel. is currently ignored, i.e. filtering will always performed on all MS levels of the object.

  • chromatogram: extract chromatographic data from a data set. Parameters mz and rt allow to define specific m/z - retention time regions to extract the data from (to e.g. for extracted ion chromatograms EICs). Both parameters are expected to be numerical two-column matrices with the first column defining the lower and the second the upper margin. Each row can define a separate m/z - retention time region. Currently the function returns a MChromatograms() object for object being a MsExperiment or, for object being an XcmsExperiment, either a MChromatograms or XChromatograms() depending on parameter return.type (can be either "MChromatograms" or "XChromatograms"). For the latter also chromatographic peaks detected within the provided m/z and retention times are returned. Parameter chromPeaks allows to specify which chromatographic peaks should be reported. See documentation on the chromPeaks parameter for more information. If the XcmsExperiment contains correspondence results, also the associated feature definitions will be included in the returned XChromatograms. By default the function returns chromatograms from MS1 data, but by setting parameter msLevel = 2L it is possible to e.g. extract also MS2 chromatograms. For msLevel other than 1 it is in addition important to also specify the isolationWindowTargetMz for which MS2 data should be extracted (e.g. for SWATH data MS2 spectra are created for different m/z isolation windows and the isolationWindowTargetMz parameter allows to define from which of these the MS2 chromatogram should be extracted. Note that in future more efficient data structures for chromatographic data will be available as well.

  • chromPeaks: returns a numeric matrix with the identified chromatographic peaks. Each row represents a chromatographic peak identified in one sample (file). The number of columns depends on the peak detection algorithm (see findChromPeaks()) but most methods return the following columns: "mz" (intensity-weighted mean of the m/z values of all mass peaks included in the chromatographic peak), "mzmin" ( smallest m/z value of any mass peak in the chromatographic peak), "mzmax" (largest m/z value of any mass peak in the chromatographic peak), "rt" (retention time of the peak apex), "rtmin" (retention time of the first scan/mass peak of the chromatographic peak), "rtmax" (retention time of the last scan/mass peak of the chromatographic peak), "into" (integrated intensity of the chromatographic peak), "maxo" (maximal intensity of any mass peak of the chromatographic peak), "sample" (index of the sample in object in which the peak was identified). Parameters rt, mz, ppm, msLevel and type allow to extract subsets of identified chromatographic peaks from the object. See parameter description below for details.

  • chromPeakData: returns a DataFrame with potential additional annotations for the identified chromatographic peaks. Each row in this DataFrame corresponds to a row (same index and row name) in the chromPeaks matrix. The default annotations are "ms_level" (the MS level in which the peak was identified) and "is_filled" (whether the chromatographic peak was detected (by findChromPeaks) or filled-in (by fillChromPeaks).

  • chromPeakSpectra: extract MS spectra for identified chromatographic peaks. This can be either all (full scan) MS1 spectra with retention times between the retention time range of a chromatographic peak, all MS2 spectra (if present) with a retention time within the retention time range of a (MS1) chromatographic peak and a precursor m/z within the m/z range of the chromatographic peak or single, selected spectra depending on their total signal or highest signal. Parameter msLevel allows to define from which MS level spectra should be extracted, parameter method allows to define if all or selected spectra should be returned. See chromPeakSpectra() for details.

  • dropChromPeaks: removes (all) chromatographic peak detection results from object. This will also remove any correspondence results (i.e. features) and eventually present adjusted retention times from the object if the alignment was performed after the peak detection. Alignment results (adjusted retention times) can be retained if parameter keepAdjustedRtime is set to TRUE.

  • dropFilledChromPeaks: removes chromatographic peaks added by gap filling with fillChromPeaks.

  • fillChromPeaks: perform gap filling to integrate signal missing values in samples in which no chromatographic peak was found. This depends on correspondence results, hence groupChromPeaks needs to be called first. For details and options see fillChromPeaks().

  • findChromPeaks: perform chromatographic peak detection. See findChromPeaks() for details.

  • hasChromPeaks: whether the object contains peak detection results. Parameter msLevel allows to check whether peak detection results are available for the specified MS level(s).

  • hasFilledChromPeaks: whether gap-filling results (i.e., filled-in chromatographic peaks) are present.

  • manualChromPeaks: manually add chromatographic peaks by defining their m/z and retention time ranges. See manualChromPeaks() for details and examples.

  • plotChromPeakImage: show the density of identified chromatographic peaks per file along the retention time. See plotChromPeakImage() for details.

  • plotChromPeaks: indicate identified chromatographic peaks from one sample in the RT-m/z space. See plotChromPeaks() for details.

  • refineChromPeaks: refines identified chromatographic peaks in object. See refineChromPeaks() for details.

  • adjustedRtime: extract adjusted retention times. This is just an alias for rtime(object, adjusted = TRUE).

  • adjustRtime: performs retention time adjustment (alignment) of the data. See adjustRtime() for details.

  • applyAdjustedRtime: replaces the original (raw) retention times with the adjusted ones. See applyAdjustedRtime() for more information.

  • dropAdjustedRtime: drops alignment results (adjusted retention time) from the result object. This also reverts the retention times of identified chromatographic peaks if present in the result object. Note that any results from a correspondence analysis (i.e. feature definitions) will be dropped too (if the correspondence analysis was performed after the alignment). This can be overruled with keepAdjustedRtime = TRUE.

  • hasAdjustedRtime: whether alignment was performed on the object (i.e., the object contains alignment results).

  • plotAdjustedRtime: plot the alignment results; see plotAdjustedRtime() for more information.

  • dropFeatureDefinitions: removes any correspondence analysis results from object as well as any filled-in chromatographic peaks. By default (with parameter keepAdjustedRtime = FALSE) also all alignment results will be removed if alignment was performed after the correspondence analysis. This can be overruled with keepAdjustedRtime = TRUE.

  • featureArea: returns a matrix with columns "mzmin", "mzmax", "rtmin" and "rtmax" with the m/z and retention time range for each feature (row) in object. By default these represent the minimal m/z and retention times as well as maximal m/z and retention times for all chromatographic peaks assigned to that feature. Parameter features allows to extract these values for selected features only. Parameters mzmin, mzmax, rtmin and rtmax allow to define the function to calculate the reported "mzmin", "mzmax", "rtmin" and "rtmax" values.

  • featureChromatograms: extract ion chromatograms (EICs) for each feature in object. See featureChromatograms() for more details.

  • featureDefinitions: returns a data.frame with feature definitions or an empty data.frame if no correspondence analysis results are present. Parameters msLevel, mz, ppm and rt allow to define subsets of feature definitions that should be returned with the parameter type defining how these parameters should be used to subset the returned data.frame. See parameter descriptions for details.

  • featureSpectra: returns a Spectra() or List of Spectra with (MS1 or MS2) spectra associated to each feature. See featureSpectra() for more details and available parameters.

  • featuresSummary: calculate a simple summary on features. See featureSummary() for details.

  • groupChromPeaks: performs the correspondence analysis (i.e., grouping of chromatographic peaks into LC-MS features). See groupChromPeaks() for details.

  • hasFeatures: whether correspondence analysis results are presentin in object. The optional parameter msLevel allows to define the MS level(s) for which it should be determined if feature definitions are available.

  • overlappingFeatures: identify features that overlapping or close in m/z - rt dimension. See overlappingFeatures() for more information.

Extracting data and results from an XcmsExperiment

Preprocessing results can be extracted using the following functions:

  • chromPeaks: extract identified chromatographic peaks. See section on chromatographic peak detection for details.

  • featureDefinitions: extract the definition of features (chromatographic peaks grouped across samples). See section on correspondence analysis for details.

  • featureValues: extract a matrix of values for features from each sample (file). Rows are features, columns samples. Which value should be returned can be defined with parameter value, which can be any column of the chromPeaks matrix. By default (value = "into") the integrated chromatographic peak intensities are returned. With parameter msLevel it is possible to extract values for features from certain MS levels. During correspondence analysis, more than one chromatographic peak per sample can be assigned to the same feature (e.g. if they are very close in retention time). Parameter method allows to define the strategy to deal with such cases: method = "medret": report the value from the chromatographic peak with the apex position closest to the feautre's median retention time. method = "maxint": report the value from the chromatographic peak with the largest signal (parameter intensity allows to define the column in chromPeaks that should be selected; defaults to intensity = "into"). method = "sum"`: sum the values for all chromatographic peaks assigned to the feature in the same sample.

  • quantify: extract the correspondence analysis results as a SummarizedExperiment(). The feature values are used as assay in the returned SummarizedExperiment, rowData contains the featureDefinitions (without column "peakidx") and colData the sampleData of object. Additional parameters to the featureValues function (that is used to extract the feature value matrix) can be passed via ....

Visualization

  • plot: plot for each file the position of individual peaks in the m/z - retention time space (with color-coded intensity) and a base peak chromatogram. This function should ideally be called only on a data subset (i.e. after using filterRt and filterMz to restrict to a region of interest). Parameter msLevel allows to define from which MS level the plot should be created. If x is a XcmsExperiment with available identified chromatographic peaks, also the region defining the peaks are indicated with a rectangle. Parameter peakCol allows to define the color of the border for these rectangles.

  • plotAdjustedRtime: plot the alignment results; see plotAdjustedRtime() for more information.

  • plotChromPeakImage: show the density of identified chromatographic peaks per file along the retention time. See plotChromPeakImage() for details.

  • plotChromPeaks: indicate identified chromatographic peaks from one sample in the RT-m/z space. See plotChromPeaks() for details.

General functionality and functions for backward compatibility

  • uniqueMsLevels: returns the unique MS levels of the spectra in object.

The functions listed below ensure compatibility with the older XCMSnExp() xcms result object.

  • fileNames: returns the original data file names for the spectra data. Ideally, the dataOrigin or dataStorage spectra variables from the object's spectra should be used instead.

  • fromFile: returns the file (sample) index for each spectrum within object. Generally, subsetting by sample using the [ is the preferred way to get spectra from a specific sample.

  • polarity: returns the polarity information for each spectrum in object.

  • processHistory: returns a list with ProcessHistory process history objects that contain also the parameter object used for the different processings. Optional parameter type allows to query for specific processing steps.

  • rtime: extract retention times of the spectra from the MsExperiment or XcmsExperiment object. It is thus a shortcut for rtime(spectra(object)) which would be the preferred way to extract retention times from an MsExperiment. The rtime method for XcmsExperiment has an additional parameter adjusted which allows to define whether adjusted retention times (if present - adjusted = TRUE) or raw retention times (adjusted = FALSE) should be returned. By default adjusted retention times are returned if available.

Differences compared to the XCMSnExp() object

  • Subsetting by [ supports arbitrary ordering.

Author

Johannes Rainer

Examples


## Creating a MsExperiment object representing the data from an LC-MS
## experiment.
library(MsExperiment)

## Defining the raw data files
fls <- c(system.file('cdf/KO/ko15.CDF', package = "faahKO"),
         system.file('cdf/KO/ko16.CDF', package = "faahKO"),
         system.file('cdf/KO/ko18.CDF', package = "faahKO"))

## Defining a data frame with the sample characterization
df <- data.frame(mzML_file = basename(fls),
                sample = c("ko15", "ko16", "ko18"))
## Importing the data. This will initialize a `Spectra` object representing
## the raw data and assign these to the individual samples.
mse <- readMsExperiment(spectraFiles = fls, sampleData = df)

## Extract a total ion chromatogram and base peak chromatogram
## from the data
bpc <- chromatogram(mse, aggregationFun = "max")
tic <- chromatogram(mse)

## Plot them
par(mfrow = c(2, 1))
plot(bpc, main = "BPC")
plot(tic, main = "TIC")


## Extracting MS2 chromatographic data
##
## To show how MS2 chromatograms can be extracted we first load a DIA
## (SWATH) data set.
mse_dia <- readMsExperiment(system.file("TripleTOF-SWATH",
    "PestMix1_SWATH.mzML", package = "msdata"))

## Extracting MS2 chromatogram requires also to specify the isolation
## window from which to extract the data. Without that chromatograms
## will be empty:
chr_ms2 <- chromatogram(mse_dia, msLevel = 2L)
intensity(chr_ms2[[1L]])
#> numeric(0)

## First we list available isolation windows
table(isolationWindowTargetMz(spectra(mse_dia)))
#> 
#> 163.75 208.95 244.05 270.85  299.1  329.8 367.35 601.85 
#>   1000   1000   1000   1000   1000   1000   1000   1000 

## We can then extract the TIC of MS2 data for a specific isolation window
chr_ms2 <- chromatogram(mse_dia, msLevel = 2L,
    isolationWindowTargetMz = 244.05)
plot(chr_ms2)

####
## Chromatographic peak detection

## Perform peak detection on the data using the centWave algorith. Note
## that the parameters are chosen to reduce the run time of the example.
p <- CentWaveParam(noise = 10000, snthresh = 40, prefilter = c(3, 10000))
xmse <- findChromPeaks(mse, param = p)
xmse
#> Object of class XcmsExperiment 
#>  Spectra: MS1 (3834) 
#>  Experiment data: 3 sample(s)
#>  Sample data links:
#>   - spectra: 3 sample(s) to 3834 element(s).
#>  xcms results:
#>   - chromatographic peaks: 248 in MS level(s): 1 

## Have a quick look at the identified chromatographic peaks
head(chromPeaks(xmse))
#>          mz mzmin mzmax       rt    rtmin    rtmax       into       intb   maxo
#> CP001 453.2 453.2 453.2 2506.073 2501.378 2527.982  1007409.0  1007380.8  38152
#> CP002 302.0 302.0 302.0 2617.185 2595.275 2640.659   687146.6   671297.8  30552
#> CP003 344.0 344.0 344.0 2679.783 2646.919 2709.517  5210015.9  5135916.9 152320
#> CP004 430.1 430.1 430.1 2681.348 2639.094 2712.647  2395840.3  2299899.6  65752
#> CP005 366.0 366.0 366.0 2679.783 2642.224 2718.907  3365174.0  3279468.3  79928
#> CP006 343.0 343.0 343.0 2678.218 2637.529 2712.647 24147443.2 23703761.7 672064
#>          sn sample
#> CP001 38151      1
#> CP002    46      1
#> CP003    68      1
#> CP004    42      1
#> CP005    49      1
#> CP006    87      1

## Extract chromatographic peaks identified between 3000 and 3300 seconds
chromPeaks(xmse, rt = c(3000, 3300), type = "within")
#>          mz mzmin mzmax       rt    rtmin    rtmax        into        intb
#> CP021 453.2 453.2 453.2 3063.196 3035.027 3114.840  3001594.78  3001514.97
#> CP030 361.1 361.1 361.1 3147.704 3141.444 3153.964   367361.84   362278.28
#> CP031 340.3 340.3 340.3 3230.646 3225.952 3233.776    68671.25    68664.99
#> CP032 526.2 526.2 526.2 3243.166 3238.471 3246.296   374698.56   374692.30
#> CP033 313.1 313.1 313.1 3276.030 3254.121 3290.115  1199152.63  1162839.05
#> CP034 454.1 454.1 454.1 3276.030 3261.945 3294.809 12283448.95 12024951.77
#> CP108 307.1 307.1 307.1 3143.009 3121.099 3164.918  2191519.65  2063202.03
#> CP109 278.1 278.1 278.1 3196.217 3180.568 3213.432   657840.25   651048.36
#> CP110 526.1 526.1 526.1 3179.003 3150.833 3210.302 21334966.98 21147392.01
#> CP111 380.1 380.1 380.1 3152.398 3132.054 3174.308  2201517.91  2201474.09
#> CP112 380.1 380.1 380.1 3210.302 3210.302 3216.561   187031.44   187023.61
#> CP113 380.1 380.1 380.1 3202.477 3174.308 3216.561  1156780.73  1156736.91
#> CP114 380.1 380.1 380.1 3222.821 3216.561 3240.036   485818.25   485793.21
#> CP115 286.2 286.2 286.2 3258.815 3246.296 3280.725  1264118.64  1247233.70
#> CP116 308.1 308.1 308.1 3261.945 3241.601 3285.419  2066854.37  2025917.82
#> CP204 380.1 380.1 380.1 3150.835 3128.925 3171.179  2036803.99  1937655.24
#> CP205 286.2 286.2 286.2 3250.992 3233.777 3257.252   732016.23   722185.86
#> CP206 568.2 568.2 568.2 3207.173 3185.264 3232.212  3951832.25  3866503.74
#>         maxo    sn sample
#> CP021  53096 53095      1
#> CP030  49240    57      1
#> CP031  11391 11390      1
#> CP032  60800 60799      1
#> CP033  55392    50      1
#> CP034 554112    58      1
#> CP108 101400    43      2
#> CP109  30208    54      2
#> CP110 622144   251      2
#> CP111  95176 95175      2
#> CP112  25144 28791      2
#> CP113  28920 28919      2
#> CP114  23176 23535      2
#> CP115  67600    41      2
#> CP116  96272    78      2
#> CP204  82624    49      3
#> CP205  74824    68      3
#> CP206 164352   113      3

## Extract ion chromatograms (EIC) for the first two chromatographic
## peaks.
chrs <- chromatogram(xmse,
    mz = chromPeaks(xmse)[1:2, c("mzmin", "mzmax")],
    rt = chromPeaks(xmse)[1:2, c("rtmin", "rtmax")])
#> Processing chromatographic peaks

## An EIC for each sample and each of the two regions was extracted.
## Identified chromatographic peaks in the defined regions are extracted
## as well.
chrs
#> XChromatograms with 2 rows and 3 columns
#>                    1               2               3
#>      <XChromatogram> <XChromatogram> <XChromatogram>
#> [1,]        peaks: 1        peaks: 0        peaks: 0
#> [2,]        peaks: 1        peaks: 0        peaks: 0
#> phenoData with 3 variables
#> featureData with 4 variables
#> - - - xcms preprocessing - - -
#> Chromatographic peak detection:
#>  method: centWave 

## Plot the EICs for the second defined region
plot(chrs[2, ])


## Subsetting the data to the results (and data) for the second sample
a <- xmse[2]
nrow(chromPeaks(xmse))
#> [1] 248
nrow(chromPeaks(a))
#> [1] 100

## Filtering the result by retention time: keeping all spectra and
## chromatographic peaks within 3000 and 3500 seconds.
xmse_sub <- filterRt(xmse, rt = c(3000, 3500))
#> Filter spectra
xmse_sub
#> Object of class XcmsExperiment 
#>  Spectra: MS1 (960) 
#>  Experiment data: 3 sample(s)
#>  Sample data links:
#>   - spectra: 3 sample(s) to 960 element(s).
#>  xcms results:
#>   - chromatographic peaks: 79 in MS level(s): 1 
nrow(chromPeaks(xmse_sub))
#> [1] 79

## Perform an initial feature grouping to allow alignment using the
## peak groups method:
pdp <- PeakDensityParam(sampleGroups = rep(1, 3))
xmse <- groupChromPeaks(xmse, param = pdp)

## Perform alignment using the peak groups method.
pgp <- PeakGroupsParam(span = 0.4)
xmse <- adjustRtime(xmse, param = pgp)
#> Performing retention time correction using 19 peak groups.

## Visualizing the alignment results
plotAdjustedRtime(xmse)

## Performing the final correspondence analysis
xmse <- groupChromPeaks(xmse, param = pdp)

## Show the definition of the first 6 features
featureDefinitions(xmse) |> head()
#>      mzmed mzmin mzmax    rtmed    rtmin    rtmax npeaks 1      peakidx
#> FT01 279.0 279.0 279.0 2789.588 2787.430 2791.746      2 2      11, 199
#> FT02 286.2 286.2 286.2 3253.923 3245.811 3262.034      2 2     115, 205
#> FT03 300.2 300.2 300.2 3385.835 3384.068 3390.895      4 3 35, 125,....
#> FT04 301.0 301.0 301.0 2789.066 2787.430 2790.180      3 3  10, 97, 198
#> FT05 305.1 305.1 305.1 2927.922 2922.158 2933.686      2 2      14, 202
#> FT06 305.1 305.1 305.1 3000.543 2991.470 3009.616      2 2      15, 203
#>      ms_level
#> FT01        1
#> FT02        1
#> FT03        1
#> FT04        1
#> FT05        1
#> FT06        1

## Extract the feature values; show the results for the first 6 rows.
featureValues(xmse) |> head()
#>        ko15.CDF ko16.CDF   ko18.CDF
#> FT01 17140627.0       NA 16919266.9
#> FT02         NA  1264119   732016.2
#> FT03  4700903.2  5313736  5169558.2
#> FT04  3051847.8  1964444  2774885.3
#> FT05  1070389.9       NA  1983342.5
#> FT06   847473.1       NA  1003750.6

## The full results can also be extracted as a `SummarizedExperiment`
## that would eventually simplify subsequent analyses with other packages.
## Any additional parameters passed to the function are passed to the
## `featureValues` function that is called to generate the feature value
## matrix.
se <- quantify(xmse, method = "sum")

## EICs for all features can be extracted with the `featureChromatograms`
## function. Note that, depending on the data set, extracting this for
## all features might take some time. Below we extract EICs for the
## first 10 features by providing the feature IDs.
chrs <- featureChromatograms(xmse,
    features = rownames(featureDefinitions(xmse))[1:10])
chrs
#> XChromatograms with 10 rows and 3 columns
#>              ko15.CDF        ko16.CDF        ko18.CDF
#>       <XChromatogram> <XChromatogram> <XChromatogram>
#> [1,]         peaks: 1        peaks: 0        peaks: 1
#> [2,]         peaks: 0        peaks: 1        peaks: 1
#> ...              ...             ...             ... 
#> [9,]         peaks: 1        peaks: 2        peaks: 0
#> [10,]        peaks: 1        peaks: 1        peaks: 1
#> phenoData with 3 variables
#> featureData with 4 variables
#> - - - xcms preprocessing - - -
#> Chromatographic peak detection:
#>  method: centWave 
#> Correspondence:
#>  method: chromatographic peak density 
#>  10 feature(s) identified.

plot(chrs[3, ])