View source: R/XcmsExperiment-functions.R
filterFeatureDefinitions | R Documentation |
xcms
Result ObjectThe 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.
filterFeatureDefinitions(object, ...)
## S4 method for signature 'MsExperiment'
filterRt(object, rt = numeric(), ...)
## S4 method for signature 'MsExperiment'
filterMzRange(object, mz = numeric(), msLevel. = uniqueMsLevels(object))
## S4 method for signature 'MsExperiment'
filterMz(object, mz = numeric(), msLevel. = uniqueMsLevels(object))
## S4 method for signature 'MsExperiment'
filterMsLevel(object, msLevel. = uniqueMsLevels(object))
## S4 method for signature 'MsExperiment'
uniqueMsLevels(object)
## S4 method for signature 'MsExperiment'
filterFile(object, file = integer(), ...)
## S4 method for signature 'MsExperiment'
rtime(object)
## S4 method for signature 'MsExperiment'
fromFile(object)
## S4 method for signature 'MsExperiment'
fileNames(object)
## S4 method for signature 'MsExperiment'
polarity(object)
## S4 method for signature 'MsExperiment'
filterIsolationWindow(object, mz = numeric())
## S4 method for signature '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()
)
featureArea(
object,
mzmin = min,
mzmax = max,
rtmin = min,
rtmax = max,
features = character()
)
## S4 method for signature 'MsExperiment,missing'
plot(x, y, msLevel = 1L, peakCol = "#ff000060", ...)
## S4 method for signature 'XcmsExperiment,ANY,ANY,ANY'
x[i, j, ..., drop = TRUE]
## S4 method for signature 'XcmsExperiment'
filterIsolationWindow(object, mz = numeric())
## S4 method for signature 'XcmsExperiment'
filterRt(object, rt, msLevel.)
## S4 method for signature 'XcmsExperiment'
filterMzRange(object, mz = numeric(), msLevel. = uniqueMsLevels(object))
## S4 method for signature 'XcmsExperiment'
filterMsLevel(object, msLevel. = uniqueMsLevels(object))
## S4 method for signature 'XcmsExperiment'
hasChromPeaks(object, msLevel = integer())
## S4 method for signature 'XcmsExperiment'
dropChromPeaks(object, keepAdjustedRtime = FALSE)
## S4 replacement method for signature 'XcmsExperiment'
chromPeaks(object) <- value
## S4 method for signature 'XcmsExperiment'
chromPeaks(
object,
rt = numeric(),
mz = numeric(),
ppm = 0,
msLevel = integer(),
type = c("any", "within", "apex_within"),
isFilledColumn = FALSE
)
## S4 replacement method for signature 'XcmsExperiment'
chromPeakData(object) <- value
## S4 method for signature 'XcmsExperiment'
chromPeakData(
object,
msLevel = integer(),
return.type = c("DataFrame", "data.frame")
)
## S4 method for signature 'XcmsExperiment'
filterChromPeaks(
object,
keep = rep(TRUE, nrow(.chromPeaks(object))),
method = "keep",
...
)
## S4 method for signature 'XcmsExperiment'
dropAdjustedRtime(object)
## S4 method for signature 'MsExperiment'
hasAdjustedRtime(object)
## S4 method for signature 'XcmsExperiment'
rtime(object, adjusted = hasAdjustedRtime(object))
## S4 method for signature 'XcmsExperiment'
adjustedRtime(object)
## S4 method for signature 'XcmsExperiment'
hasFeatures(object, msLevel = integer())
## S4 replacement method for signature 'XcmsExperiment'
featureDefinitions(object) <- value
## S4 method for signature 'XcmsExperiment'
featureDefinitions(
object,
mz = numeric(),
rt = numeric(),
ppm = 0,
type = c("any", "within", "apex_within"),
msLevel = integer()
)
## S4 method for signature 'XcmsExperiment'
dropFeatureDefinitions(object, keepAdjustedRtime = FALSE)
## S4 method for signature 'XcmsExperiment'
filterFeatureDefinitions(object, features = integer())
## S4 method for signature 'XcmsExperiment'
hasFilledChromPeaks(object)
## S4 method for signature 'XcmsExperiment'
dropFilledChromPeaks(object)
## S4 method for signature 'XcmsExperiment'
quantify(object, ...)
## S4 method for signature 'XcmsExperiment'
featureValues(
object,
method = c("medret", "maxint", "sum"),
value = "into",
intensity = "into",
filled = TRUE,
missing = NA_real_,
msLevel = integer()
)
## S4 method for signature '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 signature 'XcmsExperiment'
processHistory(object, type)
## S4 method for signature 'XcmsExperiment'
filterFile(
object,
file,
keepAdjustedRtime = hasAdjustedRtime(object),
keepFeatures = FALSE,
...
)
object |
An |
... |
Additional optional parameters. For |
rt |
For |
mz |
For |
msLevel. |
For |
file |
For |
aggregationFun |
For |
msLevel |
|
isolationWindowTargetMz |
For |
chunkSize |
For |
return.type |
For |
BPPARAM |
For |
mzmin |
For |
mzmax |
For |
rtmin |
For |
rtmax |
For |
features |
For |
x |
An |
y |
For |
peakCol |
For |
i |
For |
j |
For |
drop |
For |
keepAdjustedRtime |
|
value |
For |
ppm |
For |
type |
For |
isFilledColumn |
For |
keep |
For |
method |
For |
adjusted |
For |
intensity |
For |
filled |
For |
missing |
For |
include |
For |
chromPeaks |
For |
keepFeatures |
for most subsetting functions ( |
[
: 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. By default, with parameter
isolationWindowTargetMz = NULL
or isolationWindowTargetMz = NA_real_
,
data from all MS2 spectra will be considered in the chromatogram
extraction. If MS2 data was generated within different m/z isolation
windows (such as e.g. with Scies SWATH data), the parameter
isolationWindowTargetMz
should be used to ensure signal is only extracted
from the respective isolation window. The isolationWindowTargetMz()
function on the Spectra
object can be used to inspect/list available
isolation windows of a data set. See also the xcms LC-MS/MS vignette for
examples and details.
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.
plotPrecursorIons
: general visualization of precursor ions of
LC-MS/MS data. See plotPrecursorIons()
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.
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 ...
.
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.
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.
XCMSnExp()
object Subsetting by [
supports arbitrary ordering.
Johannes Rainer
## 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]])
## First we list available isolation windows
table(isolationWindowTargetMz(spectra(mse_dia)))
## 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
## Have a quick look at the identified chromatographic peaks
head(chromPeaks(xmse))
## Extract chromatographic peaks identified between 3000 and 3300 seconds
chromPeaks(xmse, rt = c(3000, 3300), type = "within")
## 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")])
## 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
## 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))
nrow(chromPeaks(a))
## 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))
xmse_sub
nrow(chromPeaks(xmse_sub))
## 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)
## 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()
## Extract the feature values; show the results for the first 6 rows.
featureValues(xmse) |> head()
## 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
plot(chrs[3, ])
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.