Title: | Mass Spectrometry Quantitation |
---|---|
Description: | This package quantitates raw mass spectrometry data contained in Spectra objects into QFeaures instances. |
Authors: | Laurent Gatto [aut, cre] |
Maintainer: | Laurent Gatto <[email protected]> |
License: | Artistic-2.0 |
Version: | 0.0.2 |
Built: | 2024-10-27 04:20:45 UTC |
Source: | https://github.com/rformassspectrometry/MsQuantitation |
Function that converses a Spectra
instance into an object of
class QFeatures
by quantifying reporter ions of the MS2 spectra.
quantify_labelled_ms2(x, reporters, ...)
quantify_labelled_ms2(x, reporters, ...)
x |
An instance of class Spectra with MS2 spectra. |
reporters |
An instance of class ReporterIons. |
... |
Additional parameters passed to |
A instance of class QFeatures with as many assays as
there where acquisitions in x
.
Laurent Gatto
Function that converses a Spectra
instance into an object of
class QFeatures
by quantifying reporter ions of the MS3
spectra. The assays rowData are pulled from the respective MS2
precursor scans.
quantify_labelled_ms3(x, reporters, ...)
quantify_labelled_ms3(x, reporters, ...)
x |
An instance of class Spectra with MS2 and MS3 spectra. |
reporters |
An instance of class ReporterIons. |
... |
Additional parameters passed to |
A instance of class QFeatures with as many assays as
there where acquisitions in x
.
Laurent Gatto
The quantify()
method is the entry point for quantitation of raw
mass spectrometry data. The raw data is contained in a Spectra()
object while the details and parameters of the quantitation method
are defined in a dedicated QuantParam()
.
## S4 method for signature 'Spectra' quantify(object, param, ...)
## S4 method for signature 'Spectra' quantify(object, param, ...)
object |
An instance of class |
param |
An instance of class |
... |
additional parameters controlled parallel processing of
the data. See |
An instance of QFeatures()
with as many assays as there
where acquisitions (files) in object
.
## ---------------------------------------- ## Labelled MS2 quantitation ## ---------------------------------------- ## Test data from the msdata package f <- msdata::proteomics(pattern = "01.mzML.gz", full.names = TRUE) rw <- Spectra(f) ## Quantitation parameters p <- QuantParam(msLevel = 2L, label = TRUE, params = list(reporters = TMT6)) p quantify(rw, p) ## Simulate data from 2 files rw <- filterMsLevel(rw, 2L) rw <- setBackend(rw, MsBackendDataFrame()) rw$dataOrigin <- sample(c("file1", "file2"), length(rw), replace = TRUE) quantify(rw, p) ## ---------------------------------------- ## Labelled MS3 quantitation ## ---------------------------------------- ## Test data from the msdata pakage basename(f <- msdata::proteomics(pattern = "MS3TMT11", full.names = TRUE)) x <- Spectra(f) ## Quantitation parameters p <- QuantParam(msLevel = 3L, label = TRUE, params = list(reporters = TMT11)) p quantify(x, p)
## ---------------------------------------- ## Labelled MS2 quantitation ## ---------------------------------------- ## Test data from the msdata package f <- msdata::proteomics(pattern = "01.mzML.gz", full.names = TRUE) rw <- Spectra(f) ## Quantitation parameters p <- QuantParam(msLevel = 2L, label = TRUE, params = list(reporters = TMT6)) p quantify(rw, p) ## Simulate data from 2 files rw <- filterMsLevel(rw, 2L) rw <- setBackend(rw, MsBackendDataFrame()) rw$dataOrigin <- sample(c("file1", "file2"), length(rw), replace = TRUE) quantify(rw, p) ## ---------------------------------------- ## Labelled MS3 quantitation ## ---------------------------------------- ## Test data from the msdata pakage basename(f <- msdata::proteomics(pattern = "MS3TMT11", full.names = TRUE)) x <- Spectra(f) ## Quantitation parameters p <- QuantParam(msLevel = 3L, label = TRUE, params = list(reporters = TMT11)) p quantify(x, p)
Quantitation of mass spectrometry data is implemented through the
quantify()
method, that converts raw data objects such as
Spectra()
into QFeatures()
objects that store and hande
quantitative data. Quantitation methods for MS data a numerous and
varied, and the exact method to be applied on the raw data is
paramterised by a QuantParam
instance, i.e. a dedicated class
that stores all the parameters needed to perform a specific
quantitation method.
MS quantitation can be performed at the MS1 or MS2 (and 3)
levels. Additionally, quantitation can be label-free or
labelled. These two states are defined by the msLevel
(an
integer) and label
(a logical) parameters.
Additional details about the quantitation method are defined as a
named list, as described in the 'Details' section. The example
below illustrates how to defined labelled (label = TRUE
) MS2
(msLevel = 2L) quantitation using TMT 10-plex isobaric tagging (defined by the
TMT10' ReporterIons()
instance).
## S4 method for signature 'QuantParam' isEmpty(x) QuantParam(msLevel = NA_integer_, label = NA, params = list())
## S4 method for signature 'QuantParam' isEmpty(x) QuantParam(msLevel = NA_integer_, label = NA, params = list())
x |
A instance of class |
msLevel |
|
label |
|
params |
|
Labelled MS2 Quantitation:
The params
need only to contain a single mandatory element, a
ReporterIons()
instance, named reporters
.
Labelled MS3 Quantitation:
The params
need to contain a mandatory element, a
ReporterIons()
instance, named reporters
.
Two optional character strings, named acquisitionNum
and
precScanNum
, containing the respective column names of the
acquition number and precursor scan acquisition numbers in the
raw data to be quantified. The default values are
"acquisitionNum"
and "precScanNum"
respectively, which match
the names in standard Spectra
objects.
msLevel
integer(1)
indicating the MS level the quantitation
will be performed on.
label
params
named list()
of additional parameters.
Laurent Gatto
## default (empty) parameters QuantParam() isEmpty(QuantParam()) ## MS2 quantitation using TMT10 plex QuantParam(msLevel = 2L, label = TRUE, params = list(reporters = TMT10)) ## MS3 quantitation using TMT16 plex and non-standard acquisition ## number and precursor scan acqusition number QuantParam(msLevel = 2L, label = TRUE, params = list(reporters = TMT10, acqusitionNum = "aquisition_number", preScanNum = "prec_scan_number"))
## default (empty) parameters QuantParam() isEmpty(QuantParam()) ## MS2 quantitation using TMT10 plex QuantParam(msLevel = 2L, label = TRUE, params = list(reporters = TMT10)) ## MS3 quantitation using TMT16 plex and non-standard acquisition ## number and precursor scan acqusition number QuantParam(msLevel = 2L, label = TRUE, params = list(reporters = TMT10, acqusitionNum = "aquisition_number", preScanNum = "prec_scan_number"))
The ReporterIons
class defines a set of isobaric reporter ions
that are used for quantification of labelled MS2 data, e.g. iTRAQ
(isobaric tag for relative and absolute quantitation) or TMT
(tandem mass tags).
Many commercial reporter ions are readily available. Custom
instances can be created with the ReporterIons()
function.
ReporterIons(name, reporterNames, mz, width) ## S4 method for signature 'ReporterIons' show(object) ## S4 method for signature 'ReporterIons,ANY,ANY,ANY' x[i = "numeric", j = "missing", drop = "missing"] ## S4 method for signature 'ReporterIons' length(x) ## S4 method for signature 'ReporterIons' mz(object) ## S4 method for signature 'ReporterIons' width(x) reporterNames(x) reporterNames(x) <- value ## S4 method for signature 'ReporterIons' names(x) iTRAQ4 iTRAQ5 iTRAQ8 iTRAQ9 TMT10HCD TMT10 TMT10ETD TMT11HCD TMT11 TMT16HCD TMT16 TMT6 TMT6b TMT7 TMT7b
ReporterIons(name, reporterNames, mz, width) ## S4 method for signature 'ReporterIons' show(object) ## S4 method for signature 'ReporterIons,ANY,ANY,ANY' x[i = "numeric", j = "missing", drop = "missing"] ## S4 method for signature 'ReporterIons' length(x) ## S4 method for signature 'ReporterIons' mz(object) ## S4 method for signature 'ReporterIons' width(x) reporterNames(x) reporterNames(x) <- value ## S4 method for signature 'ReporterIons' names(x) iTRAQ4 iTRAQ5 iTRAQ8 iTRAQ9 TMT10HCD TMT10 TMT10ETD TMT11HCD TMT11 TMT16HCD TMT16 TMT6 TMT6b TMT7 TMT7b
name |
Parameter to set the |
reporterNames |
Parameter to set the |
mz |
Parameter to set the |
width |
Parameter to set the |
object |
An instance of class |
x |
An instance of class |
i |
|
j |
ignored. |
drop |
ignored. |
value |
A value for replacement. |
An object of class ReporterIons
of length 4.
An object of class ReporterIons
of length 5.
An object of class ReporterIons
of length 8.
An object of class ReporterIons
of length 9.
An object of class ReporterIons
of length 10.
An object of class ReporterIons
of length 10.
An object of class ReporterIons
of length 10.
An object of class ReporterIons
of length 11.
An object of class ReporterIons
of length 11.
An object of class ReporterIons
of length 16.
An object of class ReporterIons
of length 16.
An object of class ReporterIons
of length 6.
An object of class ReporterIons
of length 6.
An object of class ReporterIons
of length 7.
An object of class ReporterIons
of length 7.
name
character(1)
naming the instance.
reporterNames
character()
of length equal to the number of
reporter ions in the instance. Used to uniquely name each ion.
mz
numeric()
with the reporter ions' m/z values.
width
numeric(1)
defining the region where the reporter
ion can be expected. This region is calculated as mz
+/-
width
for earch reporter ion.
Laurent Gatto
Ross PL, Huang YN, Marchese JN, Williamson B, Parker K, Hattan S, Khainovski N, Pillai S, Dey S, Daniels S, Purkayastha S, Juhasz P, Martin S, Bartlet-Jones M, He F, Jacobson A, Pappin DJ. Multiplexed protein quantitation in Saccharomyces cerevisiae using amine-reactive isobaric tagging reagents. Mol Cell Proteomics, 2004 Dec;3(12):1154-69. Epub 2004 Sep 22. PubMed PMID: 15385600.
Thompson A, Schäfer J, Kuhn K, Kienle S, Schwarz J, Schmidt G, Neumann T, Johnstone R, Mohammed AK, Hamon C. Tandem mass tags: a novel quantification strategy for comparative analysis of complex protein mixtures by MS/MS. Anal Chem. 2003 Apr 15;75(8):1895-904. Erratum in: Anal Chem. 2006 Jun 15;78(12):4235. Mohammed, A Karim A and Anal Chem. 2003 Sep 15;75(18):4942. Johnstone, R. PubMed PMID: 12713048.
## Manual construction of a ReporterIons instance: ReporterIons(name = "iTRAQ4", reporterNames = c("iTRAQ4.114", "iTRAQ4.115", "iTRAQ4.116", "iTRAQ4.117"), mz = c(114.1, 115.1, 116.1, 117.1), width = 0.05) ## Some pre-defined reporter ions constructors: TMT10 iTRAQ4 TMT10[1:5] reporterNames(TMT10) names(TMT10) mz(TMT10) width(TMT10) as(TMT10, "DataFrame")
## Manual construction of a ReporterIons instance: ReporterIons(name = "iTRAQ4", reporterNames = c("iTRAQ4.114", "iTRAQ4.115", "iTRAQ4.116", "iTRAQ4.117"), mz = c(114.1, 115.1, 116.1, 117.1), width = 0.05) ## Some pre-defined reporter ions constructors: TMT10 iTRAQ4 TMT10[1:5] reporterNames(TMT10) names(TMT10) mz(TMT10) width(TMT10) as(TMT10, "DataFrame")