Title: | Quantitative Analysis of Mass Spectrometry Data |
---|---|
Description: | A complete analysis pipeline for matrix-assisted laser desorption/ionization-time-of-flight (MALDI-TOF) and other two-dimensional mass spectrometry data. In addition to commonly used plotting and processing methods it includes distinctive features, namely baseline subtraction methods such as morphological filters (TopHat) or the statistics-sensitive non-linear iterative peak-clipping algorithm (SNIP), peak alignment using warping functions, handling of replicated measurements as well as allowing spectra with different resolutions. |
Authors: | Sebastian Gibb [aut, cre] , Korbinian Strimmer [ths] , Sigurdur Smarason [ctb], Laurent Gatto [ctb] , Paolo Inglese [ctb] |
Maintainer: | Sebastian Gibb <[email protected]> |
License: | GPL (>= 3) |
Version: | 1.22.3 |
Built: | 2024-11-08 06:14:14 UTC |
Source: | https://github.com/sgibb/maldiquant |
MALDIquant provides a complete analysis pipeline for matrix-assisted laser desorption/ionization-time-of-flight (MALDI-TOF) and other two-dimensional mass spectrometry data.
In addition to commonly used plotting and processing methods it includes distinctive features, namely baseline subtraction methods such as morphological filters (TopHat) or the statistics-sensitive non-linear iterative peak-clipping algorithm (SNIP), peak alignment using warping functions, handling of replicated measurements as well as allowing spectra with different resolutions.
For a first overview see
vignette("MALDIquant-intro", package="MALDIquant")
and/or run
demo("MALDIquant")
.
Package: | MALDIquant |
License: | GPL (>= 3) |
URL: | https://strimmerlab.github.io/software/maldiquant/ |
Main classes:
MassPeaks
: Represents a peak list of a single
spectrum.
MassSpectrum
: Represents a single spectrum.
The accompanying website (see below) provides example R scripts to illustrate the functionality of this package, too.
Sebastian Gibb
Maintainer: Sebastian Gibb [email protected]
S. Gibb and K. Strimmer. 2012.
MALDIquant: a versatile R package for the analysis of mass spectrometry data.
Bioinformatics 28: 2270-2271.
doi:10.1093/bioinformatics/bts447
Website: https://strimmerlab.github.io/software/maldiquant/
Introduction:
vignette("MALDIquant-intro", package="MALDIquant")
.
Run demo files: demo("MALDIquant")
.
List all available manual pages: library(help="MALDIquant")
.
MALDIquant website: https://strimmerlab.github.io/software/maldiquant/.
more MALDIquant examples and complete analyses: https://github.com/sgibb/MALDIquantExamples/.
AbstractMassObject
is an abstract (means pure virtual)
class. It is the parent class of MassSpectrum
and
MassPeaks
.
It shouldn't create or handle by the user because it is for internal use only.
mass
:numeric
, mass or mass-to-charge ratio
intensity
:numeric
, intensities for measured
mass-to-charge ratios
metaData
:list
, some metadata to describe the
spectrum
signature(x = "AbstractMassObject", i = "numeric")
:
Extracts a range of an AbstractMassObject
object
and returns a new one.
signature(x = "AbstractMassObject")
:
Converts an AbstractMassObject
object to a
matrix with 2 columns (mass
, intensity
).
signature(object = "AbstractMassObject")
:
Accessor function for coordinates stored in object generated from imaging
mass spectrometry data.
signature(object = "AbstractMassObject",
value = "numeric|matrix")
Replacement function for coordinates used in imaging mass spectrometry
datasets.
signature(object = "AbstractMassObject")
:
Accessor function for slot intensity
.
signature(object = "AbstractMassObject",
value = "numeric")
Replacement function for slot intensity
.
signature(object = "AbstractMassObject")
:
Returns TRUE
if length of intensity
is 0 or all
intensity
values are 0.
signature(x = "AbstractMassObject")
:
Returns length of slot intensity
.
signature(x = "AbstractMassObject")
:
Extented function for adding AbstractMassObject
object as a
line to a specific plot.
See lines
for details.
signature(object = "AbstractMassObject")
:
Accessor function for slot mass
.
signature(object = "AbstractMassObject",
value = "numeric")
Replacement function for slot mass
.
signature(object = "AbstractMassObject")
:
Accessor function for slot mass
.
signature(object = "AbstractMassObject",
value = "numeric")
Replacement function for slot mass
.
signature(object = "AbstractMassObject")
:
Accessor function for slot metaData
.
signature(object = "AbstractMassObject")
:
Replacement function for slot metaData
.
signature(x = "AbstractMassObject", y = "missing")
:
Extented function for plotting an AbstractMassObject
object.
See plot,AbstractMassObject,missing-method
for
details.
signature(x = "AbstractMassObject")
:
Extented function for adding
AbstractMassObject
object as points
to a specific plot.
See points
for details.
signature(object = "AbstractMassObject",
range = "numeric")
: Trim an AbstractMassObject
object. See trim,AbstractMassObject,numeric-method
for details.
signature(object = "AbstractMassObject")
:
Transforms the intensities of an AbstractMassObject
object.
See transformIntensity,AbstractMassObject-method
for details.
Sebastian Gibb [email protected]
MassPeaks
,
MassSpectrum
,
plot,AbstractMassObject,missing-method
,
transformIntensity,AbstractMassObject-method
,
trim,AbstractMassObject,numeric-method
Website: https://strimmerlab.github.io/software/maldiquant/
## load package library("MALDIquant") ## create example spectrum s <- createMassSpectrum(mass=1:10, intensity=11:20, metaData=list(name="Example Spectrum")) ## get intensity intensity(s) ## get mass mass(s) ## get metaData metaData(s) ## replace metaData metaData(s) <- list(name="Spectrum") ## trim spectrum trim(s, c(2, 9)) ## select a range s[3:6]
## load package library("MALDIquant") ## create example spectrum s <- createMassSpectrum(mass=1:10, intensity=11:20, metaData=list(name="Example Spectrum")) ## get intensity intensity(s) ## get mass mass(s) ## get metaData metaData(s) ## replace metaData metaData(s) <- list(name="Spectrum") ## trim spectrum trim(s, c(2, 9)) ## select a range s[3:6]
This function aligns a list of MassSpectrum
objects
(spectra alignment is also known as warping/phase correction).
alignSpectra(spectra, halfWindowSize=20, noiseMethod="MAD", SNR=2, reference, tolerance=0.002, warpingMethod="lowess", allowNoMatches=FALSE, emptyNoMatches=FALSE, ...)
alignSpectra(spectra, halfWindowSize=20, noiseMethod="MAD", SNR=2, reference, tolerance=0.002, warpingMethod="lowess", allowNoMatches=FALSE, emptyNoMatches=FALSE, ...)
spectra |
|
halfWindowSize |
|
noiseMethod |
a noise estimation method; see
|
SNR |
single numeric value. |
reference |
|
tolerance |
|
warpingMethod |
used basic warping function; see
|
allowNoMatches |
|
emptyNoMatches |
|
... |
arguments to be passed to
|
alignSpectra
is a wrapper function around
detectPeaks
,
determineWarpingFunctions
and
warpMassSpectra
. Please call these functions
manually if you need finer control (e.g. plotting of warping functions).
Returns a list
of aligned MassSpectrum
objects.
Sebastian Gibb [email protected]
detectPeaks
,
determineWarpingFunctions
,
referencePeaks
,
warpMassSpectra
,
MassSpectrum
demo("warping")
Website: https://strimmerlab.github.io/software/maldiquant/
## load package library("MALDIquant") ## load example data data("fiedler2009subset", package="MALDIquant") ## running typical workflow ## transform intensities spectra <- transformIntensity(fiedler2009subset, method="sqrt") ## smooth spectra spectra <- smoothIntensity(spectra, method="MovingAverage") ## baseline correction spectra <- removeBaseline(spectra) ## align spectra spectra <- alignSpectra(spectra)
## load package library("MALDIquant") ## load example data data("fiedler2009subset", package="MALDIquant") ## running typical workflow ## transform intensities spectra <- transformIntensity(fiedler2009subset, method="sqrt") ## smooth spectra spectra <- smoothIntensity(spectra, method="MovingAverage") ## baseline correction spectra <- removeBaseline(spectra) ## align spectra spectra <- alignSpectra(spectra)
MassSpectrum
objects.
This function averages MassSpectrum
objects.
averageMassSpectra(l, labels, method=c("mean", "median", "sum"), ...)
averageMassSpectra(l, labels, method=c("mean", "median", "sum"), ...)
l |
|
labels |
|
method |
used aggregation function. |
... |
arguments to be passed to underlying functions (currently only
|
The mass of the averaged MassSpectrum
object
will be the mass of the first non-empty
MassSpectrum
object (of each group).
Returns a single (no labels
given) or a list
(labels
given) of averaged MassSpectrum
objects.
Sebastian Gibb [email protected]
Website: https://strimmerlab.github.io/software/maldiquant/
## load package library("MALDIquant") ## create four MassSpectrum objects and add them to a list s <- list(createMassSpectrum(mass=1:5, intensity=1:5), createMassSpectrum(mass=1:5, intensity=1:5), createMassSpectrum(mass=1:5, intensity=6:10), createMassSpectrum(mass=1:5, intensity=6:10)) ## average all four MassSpectrum objects into a single new one ## by sum their intensities ## (no labels, returns only one new MassSpectrum object) summedSpectra <- averageMassSpectra(s, method="sum") ## only average MassSpectrum objects in a group ## (e.g. useful for technical replicates) ## (two different labels, returns a list of two new MassPeaks objects) groups <- factor(c("a", "a", "b", "b"), levels=c("a", "b")) averagedSpectra <- averageMassSpectra(s, labels=groups, method="mean")
## load package library("MALDIquant") ## create four MassSpectrum objects and add them to a list s <- list(createMassSpectrum(mass=1:5, intensity=1:5), createMassSpectrum(mass=1:5, intensity=1:5), createMassSpectrum(mass=1:5, intensity=6:10), createMassSpectrum(mass=1:5, intensity=6:10)) ## average all four MassSpectrum objects into a single new one ## by sum their intensities ## (no labels, returns only one new MassSpectrum object) summedSpectra <- averageMassSpectra(s, method="sum") ## only average MassSpectrum objects in a group ## (e.g. useful for technical replicates) ## (two different labels, returns a list of two new MassPeaks objects) groups <- factor(c("a", "a", "b", "b"), levels=c("a", "b")) averagedSpectra <- averageMassSpectra(s, labels=groups, method="mean")
This function looks for similar peaks (mass) across
MassPeaks
objects and equalizes their mass.
binPeaks(l, method=c("strict", "relaxed", "reference"), tolerance=0.002)
binPeaks(l, method=c("strict", "relaxed", "reference"), tolerance=0.002)
l |
|
method |
bin creation rule. |
tolerance |
|
The algorithm is based on the following workflow:
Put all mass in a sorted vector.
Calculate differences between each neighbor.
Divide the mass vector at the largest gap (largest difference) and form a left and a right bin.
Rerun step 3 for the left and/or the right bin if they don't fulfill the following criteria:
All peaks in a bin are near to the mean (method == "strict"
or method == "relaxed"
)
(abs(mass-meanMass)/meanMass < tolerance
) or the
reference mass (method == "reference"
;
abs(mass-reference)/reference < tolerance
).
method == "strict"
: The bin doesn't contain two or more
peaks of the same sample.
method == "strict"
: The new peak positions (mass value) are the
mean mass of a bin. method == "relaxed"
: The new peak positions for the highest peaks of each
sample in a bin are generated by the mean mass of this peaks. The lower
peaks are not changed.
method == "reference"
: The new peak positions for the highest peaks of
each sample in a bin are generated by the mass of peaks of the first
MassPeaks
object. Lower peaks are not changed.
Returns a list
of mass adjusted
MassPeaks
objects.
Sebastian Gibb [email protected]
Website: https://strimmerlab.github.io/software/maldiquant/
## load package library("MALDIquant") ## create two MassPeaks objects p <- list(createMassPeaks(mass=seq(100, 500, 100), intensity=1:5), createMassPeaks(mass=c(seq(100.2, 300.2, 100), 395), intensity=1:4)) binnedPeaks <- binPeaks(p, tolerance=0.002) ## compare result iM1 <- intensityMatrix(p) iM2 <- intensityMatrix(binnedPeaks) all(dim(iM1) == c(2, 9)) # TRUE all(dim(iM2) == c(2, 6)) # TRUE show(iM2) ## increase tolerance binnedPeaks <- binPeaks(p, tolerance=0.1) iM3 <- intensityMatrix(binnedPeaks) all(dim(iM3) == c(2, 5)) # TRUE show(iM3) ## differences between "strict" and "relaxed" p <- c(createMassPeaks(mass=c(1, 1.01, 3), intensity=c(2, 1, 1)), createMassPeaks(mass=c(0.99, 3), intensity=rep(1, 2)), createMassPeaks(mass=c(1.02, 3), intensity=rep(1, 2))) intensityMatrix(binPeaks(p, method="strict", tolerance=0.05)) intensityMatrix(binPeaks(p, method="relaxed", tolerance=0.05)) ## use a reference ref <- createMassPeaks(mass=c(1, 3), intensity=rep(1, 2)) ## include the reference intensityMatrix(binPeaks(c(ref, p), method="reference", tolerance=0.05)) ## drop the reference intensityMatrix(binPeaks(c(ref, p), method="reference", tolerance=0.05)[-1])
## load package library("MALDIquant") ## create two MassPeaks objects p <- list(createMassPeaks(mass=seq(100, 500, 100), intensity=1:5), createMassPeaks(mass=c(seq(100.2, 300.2, 100), 395), intensity=1:4)) binnedPeaks <- binPeaks(p, tolerance=0.002) ## compare result iM1 <- intensityMatrix(p) iM2 <- intensityMatrix(binnedPeaks) all(dim(iM1) == c(2, 9)) # TRUE all(dim(iM2) == c(2, 6)) # TRUE show(iM2) ## increase tolerance binnedPeaks <- binPeaks(p, tolerance=0.1) iM3 <- intensityMatrix(binnedPeaks) all(dim(iM3) == c(2, 5)) # TRUE show(iM3) ## differences between "strict" and "relaxed" p <- c(createMassPeaks(mass=c(1, 1.01, 3), intensity=c(2, 1, 1)), createMassPeaks(mass=c(0.99, 3), intensity=rep(1, 2)), createMassPeaks(mass=c(1.02, 3), intensity=rep(1, 2))) intensityMatrix(binPeaks(p, method="strict", tolerance=0.05)) intensityMatrix(binPeaks(p, method="relaxed", tolerance=0.05)) ## use a reference ref <- createMassPeaks(mass=c(1, 3), intensity=rep(1, 2)) ## include the reference intensityMatrix(binPeaks(c(ref, p), method="reference", tolerance=0.05)) ## drop the reference intensityMatrix(binPeaks(c(ref, p), method="reference", tolerance=0.05)[-1])
This function calibrates (normalize) intensities of
MassSpectrum
objects.
## S4 method for signature 'MassSpectrum' calibrateIntensity(object, method=c("TIC", "PQN", "median"), range, ...) ## S4 method for signature 'list' calibrateIntensity(object, method=c("TIC", "PQN", "median"), range, ...)
## S4 method for signature 'MassSpectrum' calibrateIntensity(object, method=c("TIC", "PQN", "median"), range, ...) ## S4 method for signature 'list' calibrateIntensity(object, method=c("TIC", "PQN", "median"), range, ...)
object |
|
method |
the calibration method to be used. This should be one of
|
range |
|
... |
arguments to be passed to other functions. Currently only
|
A number of different calibration methods are provided:
"TIC"
:The TIC (Total Ion Current) of a
MassSpectrum
object is set to one. If range
is
given the TIC is only calculated for the intensities in the
specified mass range.
"PQN"
:The PQN (Probabilistic Quotient
Normalization) is described in Dieterle et al 2006.
calibrateIntensity
uses the following algorithm:
Calibrate all spectra using the "TIC"
calibration.
Calculate a median reference spectrum.
Calculate the quotients of all intensities of the spectra with those of the reference spectrum.
Calculate the median of these quotients for each spectrum.
Divide all intensities of each spectrum by its median of quotients.
"median"
:The median of intensities of a
MassSpectrum
object is set to one.
Returns a modified MassSpectrum
object with calibrated
intensities.
Sebastian Gibb [email protected]
F. Dieterle, A. Ross, G. Schlotterbeck, and Hans Senn. 2006. Probabilistic quotient normalization as robust method to account for dilution of complex biological mixtures. Application in 1H NMR metabonomics. Analytical Chemistry 78(13): 4281-4290.
Website: https://strimmerlab.github.io/software/maldiquant/
## load package library("MALDIquant") ## load example data data("fiedler2009subset", package="MALDIquant") ## baseline correction b <- removeBaseline(fiedler2009subset) ## calibrate intensity values calibrateIntensity(b, method="TIC") ## calibrate intensity values using TIC for a specific mass range calibrateIntensity(b, method="TIC", range=c(3000, 5000))
## load package library("MALDIquant") ## load example data data("fiedler2009subset", package="MALDIquant") ## baseline correction b <- removeBaseline(fiedler2009subset) ## calibrate intensity values calibrateIntensity(b, method="TIC") ## calibrate intensity values using TIC for a specific mass range calibrateIntensity(b, method="TIC", range=c(3000, 5000))
This function creates a MassPeaks
object. Normally it
shouldn't called by the user. Try
detectPeaks,MassSpectrum-method
instead.
createMassPeaks(mass, intensity, snr=rep.int(NA_real_, length(intensity)), metaData=list())
createMassPeaks(mass, intensity, snr=rep.int(NA_real_, length(intensity)), metaData=list())
mass |
|
intensity |
|
snr |
|
metaData |
|
Returns a MassPeaks
object.
Sebastian Gibb [email protected]
detectPeaks,MassSpectrum-method
,
MassPeaks
Website: https://strimmerlab.github.io/software/maldiquant/
## load package library("MALDIquant") ## create a MassPeaks object by default constructor s <- createMassPeaks(mass=1:100, intensity=rnorm(100)^2, metaData=list(name="example peaks")) ## show some details s
## load package library("MALDIquant") ## create a MassPeaks object by default constructor s <- createMassPeaks(mass=1:100, intensity=rnorm(100)^2, metaData=list(name="example peaks")) ## show some details s
This function creates a MassSpectrum
object.
createMassSpectrum(mass, intensity, metaData=list())
createMassSpectrum(mass, intensity, metaData=list())
mass |
|
intensity |
|
metaData |
|
Returns a MassSpectrum
object.
Sebastian Gibb [email protected]
Website: https://strimmerlab.github.io/software/maldiquant/
## load package library("MALDIquant") ## create a MassSpectrum object by default constructor s <- createMassSpectrum(mass=1:100, intensity=rnorm(100)^2, metaData=list(name="example spectrum")) ## show some details s
## load package library("MALDIquant") ## create a MassSpectrum object by default constructor s <- createMassSpectrum(mass=1:100, intensity=rnorm(100)^2, metaData=list(name="example spectrum")) ## show some details s
This method looks for peaks in mass spectrometry data
(represented by a MassSpectrum
object).
A peak is a local maximum above a user defined noise threshold.
## S4 method for signature 'MassSpectrum' detectPeaks(object, halfWindowSize=20, method=c("MAD", "SuperSmoother"), SNR=2, ...) ## S4 method for signature 'list' detectPeaks(object, ...)
## S4 method for signature 'MassSpectrum' detectPeaks(object, halfWindowSize=20, method=c("MAD", "SuperSmoother"), SNR=2, ...) ## S4 method for signature 'list' detectPeaks(object, ...)
object |
|
halfWindowSize |
|
method |
a noise estimation function; see
|
SNR |
single numeric value. |
... |
arguments to be passed to
|
Returns a MassPeaks
object.
Sebastian Gibb [email protected]
MassPeaks
,
MassSpectrum
,
estimateNoise,MassSpectrum-method
demo("peaks")
Website: https://strimmerlab.github.io/software/maldiquant/
## load package library("MALDIquant") ## load example data data("fiedler2009subset", package="MALDIquant") ## choose only the first mass spectrum s <- fiedler2009subset[[1]] ## transform intensities s <- transformIntensity(s, method="sqrt") ## smoothing spectrum s <- smoothIntensity(s, method="MovingAverage") ## remove baseline s <- removeBaseline(s) ## plot spectrum plot(s) ## call peak detection p <- detectPeaks(s) ## draw peaks on the plot points(p) ## label 10 highest peaks top10 <- intensity(p) %in% sort(intensity(p), decreasing=TRUE)[1:10] labelPeaks(p, index=top10)
## load package library("MALDIquant") ## load example data data("fiedler2009subset", package="MALDIquant") ## choose only the first mass spectrum s <- fiedler2009subset[[1]] ## transform intensities s <- transformIntensity(s, method="sqrt") ## smoothing spectrum s <- smoothIntensity(s, method="MovingAverage") ## remove baseline s <- removeBaseline(s) ## plot spectrum plot(s) ## call peak detection p <- detectPeaks(s) ## draw peaks on the plot points(p) ## label 10 highest peaks top10 <- intensity(p) %in% sort(intensity(p), decreasing=TRUE)[1:10] labelPeaks(p, index=top10)
This function determines a warping function for a list of
AbstractMassObject
objects
(warping is also known as phase correction/spectra alignment).
determineWarpingFunctions(l, reference, tolerance=0.002, method=c("lowess", "linear", "quadratic", "cubic"), allowNoMatches=FALSE, plot=FALSE, plotInteractive=FALSE, ...)
determineWarpingFunctions(l, reference, tolerance=0.002, method=c("lowess", "linear", "quadratic", "cubic"), allowNoMatches=FALSE, plot=FALSE, plotInteractive=FALSE, ...)
l |
|
reference |
|
tolerance |
|
method |
used basic warping function. |
allowNoMatches |
|
plot |
|
plotInteractive |
|
... |
arguments to be passed to |
warpingFunction
: determineWarpingFunctions
estimates a warping
function to overcome the difference between mass in reference
and in
the current sample. To calculate the differences each reference peak would
match with the highest sample peak in the nearer neighborhood
(defined by mass of reference peak*tolerance
). allowNoMatches
: If allowNoMatches
is TRUE
a warning instead
of an error is thrown if an MassPeaks
object could not match to the reference
. The returned list of warping
functions will contain NA
for this object (same index in the list).
plotInteractive
: If plot
is TRUE
a lot of output is
created (each sample in l
gets its own plot).
That's why an non-interactive devices is recommended:
## create a device pdf() ## calculate warping functions w <- determineWarpingFunctions(p, plot=TRUE) ## close device dev.off()
Returns a list
of individual warping functions. The attribute
nmatch
contains the number of matches of each
MassPeaks
element in l
against reference
.
Sebastian Gibb [email protected]
referencePeaks
,
warpMassPeaks
,
warpMassSpectra
,
MassPeaks
demo("warping")
Website: https://strimmerlab.github.io/software/maldiquant/
## load package library("MALDIquant") ## create a reference MassPeaks object r <- createMassPeaks(mass=1:5, intensity=1:5) ## create test samples p <- list(createMassPeaks(mass=((1:5)*1.01), intensity=1:5), createMassPeaks(mass=((1:5)*0.99), intensity=1:5)) ## create an interactive device with 2 rows par(mfrow=c(2, 1)) ## calculate warping function ## (using a linear function as basic warping function) ## and show warping plot w <- determineWarpingFunctions(p, tolerance=0.02, method="linear", plot=TRUE, plotInteractive=TRUE) par(mfrow=c(1, 1)) ## access number of matches attr(w, "nmatch") ## w contains the individual warping functions warpedPeaks <- warpMassPeaks(p, w) ## compare results all(mass(r) == mass(warpedPeaks[[1]])) # TRUE all(mass(r) == mass(warpedPeaks[[2]])) # TRUE ## realistic example ## load example data data("fiedler2009subset", package="MALDIquant") ## running typical workflow ## use only four spectra of the subset spectra <- fiedler2009subset[1:4] ## transform intensities spectra <- transformIntensity(spectra, method="sqrt") ## smooth spectra spectra <- smoothIntensity(spectra, method="MovingAverage") ## baseline correction spectra <- removeBaseline(spectra) ## detect peaks peaks <- detectPeaks(spectra) ## create an interactive device with 2 rows par(mfrow=c(4, 1)) ## calculate warping functions (using LOWESS based basic function [default]) w <- determineWarpingFunctions(peaks, plot=TRUE, plotInteractive=TRUE) par(mfrow=c(1, 1)) ## realistic example with user defined reference/calibration peaks ## use the workflow above for fiedler2009subset ## create reference peaks refPeaks <- createMassPeaks(mass=c(1207, 1264, 1351, 1466, 1616, 2769, 2932, 3191, 3262, 4091, 4209, 5904, 7762, 9285), intensity=rep(1, 14)) ## create an interactive device with 2 rows par(mfrow=c(4, 1)) ## calculate warping functions (using a quadratic function as basic function) w <- determineWarpingFunctions(peaks, reference=refPeaks, method="quadratic", plot=TRUE, plotInteractive=TRUE) par(mfrow=c(1, 1))
## load package library("MALDIquant") ## create a reference MassPeaks object r <- createMassPeaks(mass=1:5, intensity=1:5) ## create test samples p <- list(createMassPeaks(mass=((1:5)*1.01), intensity=1:5), createMassPeaks(mass=((1:5)*0.99), intensity=1:5)) ## create an interactive device with 2 rows par(mfrow=c(2, 1)) ## calculate warping function ## (using a linear function as basic warping function) ## and show warping plot w <- determineWarpingFunctions(p, tolerance=0.02, method="linear", plot=TRUE, plotInteractive=TRUE) par(mfrow=c(1, 1)) ## access number of matches attr(w, "nmatch") ## w contains the individual warping functions warpedPeaks <- warpMassPeaks(p, w) ## compare results all(mass(r) == mass(warpedPeaks[[1]])) # TRUE all(mass(r) == mass(warpedPeaks[[2]])) # TRUE ## realistic example ## load example data data("fiedler2009subset", package="MALDIquant") ## running typical workflow ## use only four spectra of the subset spectra <- fiedler2009subset[1:4] ## transform intensities spectra <- transformIntensity(spectra, method="sqrt") ## smooth spectra spectra <- smoothIntensity(spectra, method="MovingAverage") ## baseline correction spectra <- removeBaseline(spectra) ## detect peaks peaks <- detectPeaks(spectra) ## create an interactive device with 2 rows par(mfrow=c(4, 1)) ## calculate warping functions (using LOWESS based basic function [default]) w <- determineWarpingFunctions(peaks, plot=TRUE, plotInteractive=TRUE) par(mfrow=c(1, 1)) ## realistic example with user defined reference/calibration peaks ## use the workflow above for fiedler2009subset ## create reference peaks refPeaks <- createMassPeaks(mass=c(1207, 1264, 1351, 1466, 1616, 2769, 2932, 3191, 3262, 4091, 4209, 5904, 7762, 9285), intensity=rep(1, 14)) ## create an interactive device with 2 rows par(mfrow=c(4, 1)) ## calculate warping functions (using a quadratic function as basic function) w <- determineWarpingFunctions(peaks, reference=refPeaks, method="quadratic", plot=TRUE, plotInteractive=TRUE) par(mfrow=c(1, 1))
This method estimates the baseline of mass spectrometry data
(represented by a MassSpectrum
object).
## S4 method for signature 'MassSpectrum' estimateBaseline(object, method=c("SNIP", "TopHat", "ConvexHull", "median"), ...)
## S4 method for signature 'MassSpectrum' estimateBaseline(object, method=c("SNIP", "TopHat", "ConvexHull", "median"), ...)
object |
|
method |
used baseline estimation method, one of
|
... |
arguments to be passed to |
"SNIP"
:This baseline estimation is based on the Statistics-sensitive Non-linear
Iterative Peak-clipping algorithm (SNIP) described in Ryan et al 1988.
The algorithm based on the following equation:
It has two additional arguments namely iterations
and
decreasing
. iterations
controls the window size (k;
similar to halfWindowSize
in "TopHat"
, "Median"
) of
the algorithm. The resulting window reaches from
mass[cur_index-iterations]
to mass[cur_index+iterations]
.
decreasing
: In Morhac 2009 a decreasing clipping window is
suggested to get a smoother baseline. For decreasing = TRUE
(decreasing = FALSE
) k=iterations
is decreased
(increased) by one until zero (iterations
) is reached. The default
setting is decreasing = TRUE
.
"TopHat"
:This algorithm applies a moving minimum (erosion filter) and subsequently
a moving maximum (dilation filter) filter on the intensity values. The
implementation is based on van Herk 1996.
It has an additional halfWindowSize
argument determining the half
size of the moving window for the TopHat filter. The resulting window
reaches from mass[cur_index-halfWindowSize]
to
mass[cur_index+halfWindowSize]
.
"ConvexHull"
:The baseline estimation is based on a convex hull constructed below the spectrum.
"median"
:This baseline estimation uses a moving median. It is based on
runmed
. The additional argument halfWindowSize
corresponds to the k
argument in runmed
(k = 2 * halfWindowSize + 1
) and controls the half size of the
moving window. The resulting window reaches from
mass[cur_index-halfWindowSize]
to
mass[cur_index+halfWindowSize]
.
Returns a two column matrix (first column: mass, second column: intensity) of the estimated baseline.
Sebastian Gibb [email protected]
"SNIP"
:
C.G. Ryan, E. Clayton, W.L. Griffin, S.H. Sie, and D.R. Cousens. 1988.
Snip, a statistics-sensitive background treatment for the quantitative analysis
of pixe spectra in geoscience applications.
Nuclear Instruments and Methods in Physics Research Section B:
Beam Interactions with Materials and Atoms, 34(3): 396-402.
M. Morhac. 2009. An algorithm for determination of peak regions and baseline elimination in spectroscopic data. Nuclear Instruments and Methods in Physics Research Section A: Accelerators, Spectrometers, Detectors and Associated Equipment, 600(2), 478-487.
"TopHat"
:
M. van Herk. 1992.
A Fast Algorithm for Local Minimum and Maximum Filters on Rectangular and
Octagonal Kernels.
Pattern Recognition Letters 13.7: 517-521.
J. Y. Gil and M. Werman. 1996. Computing 2-Dimensional Min, Median and Max Filters. IEEE Transactions: 504-507.
"ConvexHull"
:
Andrew, A. M. 1979.
Another efficient algorithm for convex hulls in two dimensions.
Information Processing Letters, 9(5), 216-219.
MassSpectrum
,
removeBaseline,MassSpectrum-method
demo("baseline")
Website: https://strimmerlab.github.io/software/maldiquant/
## load package library("MALDIquant") ## load example data data("fiedler2009subset", package="MALDIquant") ## choose only the first mass spectrum s <- fiedler2009subset[[1]] ## SNIP plot(s) ## estimate baseline b <- estimateBaseline(s, method="SNIP", iterations=100) ## draw baseline on the plot lines(b, col="red") ## TopHat plot(s) ## estimate baseline (try different parameters) b1 <- estimateBaseline(s, method="TopHat", halfWindowSize=75) b2 <- estimateBaseline(s, method="TopHat", halfWindowSize=150) ## draw baselines on the plot lines(b1, col=2) lines(b2, col=3) ## draw legend legend(x="topright", lwd=1, legend=paste0("halfWindowSize=", c(75, 150)), col=c(2, 3)) ## ConvexHull plot(s) ## estimate baseline b <- estimateBaseline(s, method="ConvexHull") ## draw baseline on the plot lines(b, col="red") ## Median plot(s) ## estimate baseline b <- estimateBaseline(s, method="median") ## draw baseline on the plot lines(b, col="red")
## load package library("MALDIquant") ## load example data data("fiedler2009subset", package="MALDIquant") ## choose only the first mass spectrum s <- fiedler2009subset[[1]] ## SNIP plot(s) ## estimate baseline b <- estimateBaseline(s, method="SNIP", iterations=100) ## draw baseline on the plot lines(b, col="red") ## TopHat plot(s) ## estimate baseline (try different parameters) b1 <- estimateBaseline(s, method="TopHat", halfWindowSize=75) b2 <- estimateBaseline(s, method="TopHat", halfWindowSize=150) ## draw baselines on the plot lines(b1, col=2) lines(b2, col=3) ## draw legend legend(x="topright", lwd=1, legend=paste0("halfWindowSize=", c(75, 150)), col=c(2, 3)) ## ConvexHull plot(s) ## estimate baseline b <- estimateBaseline(s, method="ConvexHull") ## draw baseline on the plot lines(b, col="red") ## Median plot(s) ## estimate baseline b <- estimateBaseline(s, method="median") ## draw baseline on the plot lines(b, col="red")
This method estimates the noise of mass spectrometry data
(represented by a MassSpectrum
object).
## S4 method for signature 'MassSpectrum' estimateNoise(object, method=c("MAD", "SuperSmoother"), ...)
## S4 method for signature 'MassSpectrum' estimateNoise(object, method=c("MAD", "SuperSmoother"), ...)
object |
|
method |
used noise estimation method, one of
|
... |
arguments to be passed to |
"MAD"
:This function estimates the noise of mass spectrometry data
by calculating the median absolute deviation, see also
mad
.
"SuperSmoother"
:This function estimates the noise of mass spectrometry data using
Friedman's Super Smoother. Please refer
supsmu
for details and additional arguments.
Returns a two column matrix (first column: mass, second column: intensity) of the estimated noise.
Sebastian Gibb [email protected]
MassSpectrum
,
detectPeaks,MassSpectrum-method
,
mad
,
supsmu
Website: https://strimmerlab.github.io/software/maldiquant/
## load package library("MALDIquant") ## load example data data("fiedler2009subset", package="MALDIquant") ## choose only the first mass spectrum s <- fiedler2009subset[[1]] ## transform intensities s <- transformIntensity(s, method="sqrt") ## remove baseline s <- removeBaseline(s) ## plot spectrum plot(s) ## estimate noise nm <- estimateNoise(s, method="MAD") nss <- estimateNoise(s, method="SuperSmoother") ## draw noise on the plot lines(nm, col=2) lines(nss, col=4) ## draw legend legend(x="topright", lwd=1, legend=c("MAD", "SuperSmoother"), col=c(2, 4))
## load package library("MALDIquant") ## load example data data("fiedler2009subset", package="MALDIquant") ## choose only the first mass spectrum s <- fiedler2009subset[[1]] ## transform intensities s <- transformIntensity(s, method="sqrt") ## remove baseline s <- removeBaseline(s) ## plot spectrum plot(s) ## estimate noise nm <- estimateNoise(s, method="MAD") nss <- estimateNoise(s, method="SuperSmoother") ## draw noise on the plot lines(nm, col=2) lines(nss, col=4) ## draw legend legend(x="topright", lwd=1, legend=c("MAD", "SuperSmoother"), col=c(2, 4))
This dataset contains 16 example mass spectra.
It is used to demonstrate the usage of MALDIquant-package
.
data(fiedler2009subset)
data(fiedler2009subset)
A list containing 16 MassSpectrum-class
objects.
The dataset is a subset of data used in Fiedler et al 2009.
It contains spectra of 8 different patients (each one has 2 technical
replicates).
list_index | laboratory | patient_id | sex | age | type |
1 | Leipzig | LC77 | male | 37 | control |
2 | Leipzig | LC77 | male | 37 | control |
3 | Leipzig | LC213 | female | 51 | control |
4 | Leipzig | LC213 | female | 51 | control |
5 | Leipzig | LT178 | male | 58 | cancer |
6 | Leipzig | LT178 | male | 58 | cancer |
7 | Leipzig | LT157 | male | 60 | cancer |
8 | Leipzig | LT157 | male | 60 | cancer |
9 | Heidelberg | HC49 | male | 43 | control |
10 | Heidelberg | HC49 | male | 43 | control |
11 | Heidelberg | HC54 | female | 71 | control |
12 | Heidelberg | HC54 | female | 71 | control |
13 | Heidelberg | HT151 | male | 53 | cancer |
14 | Heidelberg | HT151 | male | 53 | cancer |
15 | Heidelberg | HT429 | female | 58 | cancer |
16 | Heidelberg | HT429 | female | 58 | cancer |
G.M. Fiedler, A.B. Leichtle, J. Kase, S. Baumann, U. Ceglarek, K. Felix,
T. Conrad, H. Witzigmann, A. Weimann, C. Schütte, J. Hauss, M. Büchler
and J. Thiery
“Serum Peptidome Profiling Revealed Platelet Factor 4 as a Potential
Discriminating Peptide Associated with Pancreatic Cancer”
Clinical Cancer Research, 11(15): 3812-3819, 2009
ISSN 1557-3265; doi:10.1158/1078-0432.CCR-08-2701
Website: https://strimmerlab.github.io/software/maldiquant/
This function removes infrequently occuring peaks in a list of
MassPeaks
objects.
filterPeaks(l, minFrequency, minNumber, labels, mergeWhitelists=FALSE)
filterPeaks(l, minFrequency, minNumber, labels, mergeWhitelists=FALSE)
l |
|
minFrequency |
|
minNumber |
|
labels |
|
mergeWhitelists |
|
For mergeWhitelists=FALSE
the filtering uses a separate peak whitelist
for each group specified by labels
, and is done independently in each
group. For mergeWhitelists=TRUE
the peak whitelists are combined,
which means that peaks that occur frequently in at least one group are also
kept in all other groups.
If both minFrequency
and minNumber
arguments are specified the
more stringent threshold is used.
Returns a list
of filtered
MassPeaks
objects.
Sebastian Gibb [email protected]
Website: https://strimmerlab.github.io/software/maldiquant/
## load package library("MALDIquant") ## create four MassPeaks objects and add them to the list p <- list(createMassPeaks(mass=1:2, intensity=1:2), createMassPeaks(mass=1:3, intensity=1:3), createMassPeaks(mass=1:4, intensity=1:4), createMassPeaks(mass=1:5, intensity=1:5)) ## only keep peaks which occur in all MassPeaks objects filteredPeaks <- filterPeaks(p, minFrequency=1) ## compare result intensities <- intensityMatrix(filteredPeaks) ## peaks at mass 3,4,5 are removed all(dim(intensities) == c(4, 2)) # TRUE all(intensities[,1] == 1) # TRUE all(intensities[,2] == 2) # TRUE ## only keep peaks which occur in all MassPeaks objects in a group ## (e.g. useful for technical replicates) groups <- factor(c("a", "a", "b", "b"), levels=c("a", "b")) filteredPeaks <- filterPeaks(p, minFrequency=1, labels=groups) ## peaks at mass 3 were removed in group "a" filteredPeaks[groups == "a"] ## peaks at mass 5 were removed in group "b" filteredPeaks[groups == "b"] ## only keep peaks which occur at least twice in a group groups <- factor(c("a", "a", "b", "b", "b"), levels=c("a", "b")) filteredPeaks <- filterPeaks(c(p, p[[3]]), minNumber=2, labels=groups) ## peaks at mass 3 were removed in group "a" filteredPeaks[groups == "a"] ## peaks at mass 5 were removed in group "b" filteredPeaks[groups == "b"] ## apply different minFrequency arguments to each group groups <- factor(c("a", "a", "b", "b", "b"), levels=c("a", "b")) filteredPeaks <- filterPeaks(c(p, p[[3]]), minFrequency=c(1, 2/3), labels=groups) intensityMatrix(filteredPeaks) # 1 2 3 4 #[1,] 1 2 NA NA #[2,] 1 2 NA NA #[3,] 1 2 3 4 #[4,] 1 2 3 4 #[4,] 1 2 3 4 ## demonstrate the use of mergeWhitelists groups <- factor(c("a", "a", "b", "b"), levels=c("a", "b")) ## default behaviour filteredPeaks <- filterPeaks(p, minNumber=2, labels=groups) intensityMatrix(filteredPeaks) # 1 2 3 4 #[1,] 1 2 NA NA #[2,] 1 2 NA NA #[3,] 1 2 3 4 #[4,] 1 2 3 4 ## use mergeWhitelists=TRUE to keep peaks of group "a" that match all filtering ## criteria in group "b" ## (please note that mass == 3 is not removed in the second MassPeaks object) filteredPeaks <- filterPeaks(p, minNumber=2, labels=groups, mergeWhitelists=TRUE) intensityMatrix(filteredPeaks) # 1 2 3 4 #[1,] 1 2 NA NA #[2,] 1 2 3 NA #[3,] 1 2 3 4 #[4,] 1 2 3 4
## load package library("MALDIquant") ## create four MassPeaks objects and add them to the list p <- list(createMassPeaks(mass=1:2, intensity=1:2), createMassPeaks(mass=1:3, intensity=1:3), createMassPeaks(mass=1:4, intensity=1:4), createMassPeaks(mass=1:5, intensity=1:5)) ## only keep peaks which occur in all MassPeaks objects filteredPeaks <- filterPeaks(p, minFrequency=1) ## compare result intensities <- intensityMatrix(filteredPeaks) ## peaks at mass 3,4,5 are removed all(dim(intensities) == c(4, 2)) # TRUE all(intensities[,1] == 1) # TRUE all(intensities[,2] == 2) # TRUE ## only keep peaks which occur in all MassPeaks objects in a group ## (e.g. useful for technical replicates) groups <- factor(c("a", "a", "b", "b"), levels=c("a", "b")) filteredPeaks <- filterPeaks(p, minFrequency=1, labels=groups) ## peaks at mass 3 were removed in group "a" filteredPeaks[groups == "a"] ## peaks at mass 5 were removed in group "b" filteredPeaks[groups == "b"] ## only keep peaks which occur at least twice in a group groups <- factor(c("a", "a", "b", "b", "b"), levels=c("a", "b")) filteredPeaks <- filterPeaks(c(p, p[[3]]), minNumber=2, labels=groups) ## peaks at mass 3 were removed in group "a" filteredPeaks[groups == "a"] ## peaks at mass 5 were removed in group "b" filteredPeaks[groups == "b"] ## apply different minFrequency arguments to each group groups <- factor(c("a", "a", "b", "b", "b"), levels=c("a", "b")) filteredPeaks <- filterPeaks(c(p, p[[3]]), minFrequency=c(1, 2/3), labels=groups) intensityMatrix(filteredPeaks) # 1 2 3 4 #[1,] 1 2 NA NA #[2,] 1 2 NA NA #[3,] 1 2 3 4 #[4,] 1 2 3 4 #[4,] 1 2 3 4 ## demonstrate the use of mergeWhitelists groups <- factor(c("a", "a", "b", "b"), levels=c("a", "b")) ## default behaviour filteredPeaks <- filterPeaks(p, minNumber=2, labels=groups) intensityMatrix(filteredPeaks) # 1 2 3 4 #[1,] 1 2 NA NA #[2,] 1 2 NA NA #[3,] 1 2 3 4 #[4,] 1 2 3 4 ## use mergeWhitelists=TRUE to keep peaks of group "a" that match all filtering ## criteria in group "b" ## (please note that mass == 3 is not removed in the second MassPeaks object) filteredPeaks <- filterPeaks(p, minNumber=2, labels=groups, mergeWhitelists=TRUE) intensityMatrix(filteredPeaks) # 1 2 3 4 #[1,] 1 2 NA NA #[2,] 1 2 3 NA #[3,] 1 2 3 4 #[4,] 1 2 3 4
These functions looks for empty AbstractMassObject
objects
in a list
.
findEmptyMassObjects(l) removeEmptyMassObjects(l)
findEmptyMassObjects(l) removeEmptyMassObjects(l)
l |
|
findEmptyMassObjects
: Returns a vector
of indices referring to empty
AbstractMassObject
objects.
removeEmptyMassObjects
: Returns a list
of
AbstractMassObject
objects but without empty ones.
Sebastian Gibb [email protected]
isEmpty,AbstractMassObject-method
,
AbstractMassObject
Website: https://strimmerlab.github.io/software/maldiquant/
## load package library("MALDIquant") ## create list peakList <- list() ## create two MassPeaks objects and add them to the list peakList[[1]] <- createMassPeaks(mass=1:100, intensity=1:100, metaData=list(name="example 1")) peakList[[2]] <- createMassPeaks(mass=1:100, intensity=1:100, metaData=list(name="example 2")) ## find empty objects (there should not be any one) findEmptyMassObjects(peakList) ## add an empty MassPeaks object to the list peakList[[3]] <- createMassPeaks(mass=double(), intensity=double(), metaData=list(name="empty MassPeaks object")) ## look for empty objects (isEmptyIdx == 3) (isEmptyIdx <- findEmptyMassObjects(peakList)) ## to remove all empty MassObjects from a list length(peakList) # 3 peakList <- removeEmptyMassObjects(peakList) length(peakList) # 2; WARNING: all indices could changed
## load package library("MALDIquant") ## create list peakList <- list() ## create two MassPeaks objects and add them to the list peakList[[1]] <- createMassPeaks(mass=1:100, intensity=1:100, metaData=list(name="example 1")) peakList[[2]] <- createMassPeaks(mass=1:100, intensity=1:100, metaData=list(name="example 2")) ## find empty objects (there should not be any one) findEmptyMassObjects(peakList) ## add an empty MassPeaks object to the list peakList[[3]] <- createMassPeaks(mass=double(), intensity=double(), metaData=list(name="empty MassPeaks object")) ## look for empty objects (isEmptyIdx == 3) (isEmptyIdx <- findEmptyMassObjects(peakList)) ## to remove all empty MassObjects from a list length(peakList) # 3 peakList <- removeEmptyMassObjects(peakList) length(peakList) # 2; WARNING: all indices could changed
This function converts a list
of
MassPeaks
objects into a matrix
.
intensityMatrix(peaks, spectra)
intensityMatrix(peaks, spectra)
peaks |
|
spectra |
|
peaks
have to be binned by binPeaks
before
calling intensityMatrix
.
Returns a matrix
containing intensities of all
MassPeaks
objects of peaks
and interpolated
intensity values for missing peaks if spectra
was given or NA
otherwise.
The matrix
has length(peaks)
rows
(one row for each sample) and length(unique(mass))
columns.
There is an additional attribute "mass"
that stores the mass values.
Sebastian Gibb [email protected]
binPeaks
,
MassPeaks
,
MassSpectrum
Website: https://strimmerlab.github.io/software/maldiquant/
## load package library("MALDIquant") ## create example MassPeaks objects p <- list(createMassPeaks(mass=1:4, intensity=11:14, metaData=list(name="test mass peaks 1")), createMassPeaks(mass=2:5, intensity=22:25, metaData=list(name="test mass peaks 2"))) ## converts MassPeaks objects into a matrix intensityMatrix(p) ## realistic example ## load example data data("fiedler2009subset", package="MALDIquant") ## transform intensities s <- transformIntensity(fiedler2009subset, method="sqrt") ## smoothing spectrum s <- smoothIntensity(s, method="MovingAverage") ## remove baseline s <- removeBaseline(s) ## call peak detection p <- detectPeaks(s) ## bin peaks p <- binPeaks(p) ## convert MassPeaks objects into a matrix with missing intensity ## values intensityMatrix(p) ## convert MassPeaks and MassSpectrum objects into a matrix without ## missing intensity values intensityMatrix(p, s)
## load package library("MALDIquant") ## create example MassPeaks objects p <- list(createMassPeaks(mass=1:4, intensity=11:14, metaData=list(name="test mass peaks 1")), createMassPeaks(mass=2:5, intensity=22:25, metaData=list(name="test mass peaks 2"))) ## converts MassPeaks objects into a matrix intensityMatrix(p) ## realistic example ## load example data data("fiedler2009subset", package="MALDIquant") ## transform intensities s <- transformIntensity(fiedler2009subset, method="sqrt") ## smoothing spectrum s <- smoothIntensity(s, method="MovingAverage") ## remove baseline s <- removeBaseline(s) ## call peak detection p <- detectPeaks(s) ## bin peaks p <- binPeaks(p) ## convert MassPeaks objects into a matrix with missing intensity ## values intensityMatrix(p) ## convert MassPeaks and MassSpectrum objects into a matrix without ## missing intensity values intensityMatrix(p, s)
These functions test for a MassSpectrum
or
MassPeaks
object.
isMassSpectrum(x) isMassPeaks(x)
isMassSpectrum(x) isMassPeaks(x)
x |
object to be tested. |
Returns TRUE
or FALSE
depending on whether its
argument is an MassSpectrum
or
MassPeaks
object.
Sebastian Gibb [email protected]
MassPeaks
,
MassSpectrum
,
AbstractMassObject
Website: https://strimmerlab.github.io/software/maldiquant/
## load package library("MALDIquant") ## create a MassPeaks object peaks <- createMassPeaks(mass=1:100, intensity=1:100, metaData=list(name="example 1")) ## test isMassPeaks(peaks) # returns TRUE isMassSpectrum(peaks) # returns FALSE isMassPeaks(double()) # returns FALSE
## load package library("MALDIquant") ## create a MassPeaks object peaks <- createMassPeaks(mass=1:100, intensity=1:100, metaData=list(name="example 1")) ## test isMassPeaks(peaks) # returns TRUE isMassSpectrum(peaks) # returns FALSE isMassPeaks(double()) # returns FALSE
These functions test a list
whether containing
MassSpectrum
or MassSpectrum
objects.
isMassSpectrumList(x) isMassPeaksList(x)
isMassSpectrumList(x) isMassPeaksList(x)
x |
object to be tested. |
Returns TRUE
or FALSE
depending on whether its
argument is a list
of MassSpectrum
or MassPeaks
objects.
Sebastian Gibb [email protected]
MassPeaks
,
MassSpectrum
,
AbstractMassObject
Website: https://strimmerlab.github.io/software/maldiquant/
## load package library("MALDIquant") ## create list p <- list() ## test list isMassPeaksList(p) # returns FALSE ## create two MassPeaks objects and add them to the list p <- createMassPeaks(mass=1:100, intensity=1:100, metaData=list(name="example 1")) p <- createMassPeaks(mass=1:100, intensity=1:100, metaData=list(name="example 2")) ## test list isMassPeaksList(p) # returns TRUE isMassSpectrumList(p) # returns FALSE
## load package library("MALDIquant") ## create list p <- list() ## test list isMassPeaksList(p) # returns FALSE ## create two MassPeaks objects and add them to the list p <- createMassPeaks(mass=1:100, intensity=1:100, metaData=list(name="example 1")) p <- createMassPeaks(mass=1:100, intensity=1:100, metaData=list(name="example 2")) ## test list isMassPeaksList(p) # returns TRUE isMassSpectrumList(p) # returns FALSE
labelPeaks
draws the corresponding mass values on top
of the peaks stored in a MassPeaks
object to a plot.
## S4 method for signature 'MassPeaks' labelPeaks(object, index, mass, labels, digits=3, underline=TRUE, verticalOffset=abs(diff(par("usr")[3:4]))*0.01, absoluteVerticalPos, adj=c(0.5, 0), cex=0.7, srt=0, avoidOverlap=FALSE, arrowLength=0, arrowLwd=0.5, arrowCol=1, ...)
## S4 method for signature 'MassPeaks' labelPeaks(object, index, mass, labels, digits=3, underline=TRUE, verticalOffset=abs(diff(par("usr")[3:4]))*0.01, absoluteVerticalPos, adj=c(0.5, 0), cex=0.7, srt=0, avoidOverlap=FALSE, arrowLength=0, arrowLwd=0.5, arrowCol=1, ...)
object |
|
index |
|
mass |
|
labels |
|
digits |
|
underline |
logical, underline peak values? |
verticalOffset |
|
absoluteVerticalPos |
|
adj |
|
cex |
|
srt |
|
avoidOverlap |
|
arrowLength , arrowLwd , arrowCol
|
arrow parameters, possible vectors.
|
... |
arguments to be passed to |
Please note that avoidOverlap = TRUE
is just supported for
srt %% 90 == 0
(means srt
has to be a multiple of 90 degree).
Sebastian Gibb
MassPeaks
,
plot,AbstractMassObject,missing-method
Website: https://strimmerlab.github.io/software/maldiquant/
## load package library("MALDIquant") ## create a MassPeaks object from scratch p <- createMassPeaks(mass=1:20, intensity=sample(x=100:10000, size=20), metaData=list(name="example")) ## plot peaks plot(p) ## label the first 5 peaks labelPeaks(p, index=1:5) ## label all peaks in mass range 15 to 20 labelPeaks(p, mass=15:20, underline=FALSE) ## label highest peaks (top 5) top5 <- intensity(p) %in% sort(intensity(p), decreasing=TRUE)[1:5] labelPeaks(p, index=top5, col="red") ## real example data("fiedler2009subset") ## a simplified preprocessing r <- removeBaseline(fiedler2009subset[[1]]) p <- detectPeaks(r) plot(p) ## label highest peaks (top 10) and avoid label overlap top10 <- sort(intensity(p), decreasing=TRUE, index.return=TRUE)$ix[1:10] labelPeaks(p, index=top10, avoidOverlap=TRUE, digits=1) ## use own labels and rotate by 90 degree plot(p) labelPeaks(p, index=top10, labels=paste("TOP", 1:10), underline=FALSE, srt=90, adj=c(0, 0.5), col=2)
## load package library("MALDIquant") ## create a MassPeaks object from scratch p <- createMassPeaks(mass=1:20, intensity=sample(x=100:10000, size=20), metaData=list(name="example")) ## plot peaks plot(p) ## label the first 5 peaks labelPeaks(p, index=1:5) ## label all peaks in mass range 15 to 20 labelPeaks(p, mass=15:20, underline=FALSE) ## label highest peaks (top 5) top5 <- intensity(p) %in% sort(intensity(p), decreasing=TRUE)[1:5] labelPeaks(p, index=top5, col="red") ## real example data("fiedler2009subset") ## a simplified preprocessing r <- removeBaseline(fiedler2009subset[[1]]) p <- detectPeaks(r) plot(p) ## label highest peaks (top 10) and avoid label overlap top10 <- sort(intensity(p), decreasing=TRUE, index.return=TRUE)$ix[1:10] labelPeaks(p, index=top10, avoidOverlap=TRUE, digits=1) ## use own labels and rotate by 90 degree plot(p) labelPeaks(p, index=top10, labels=paste("TOP", 1:10), underline=FALSE, srt=90, adj=c(0, 0.5), col=2)
MALDIquant
offers multi-core support using
mclapply
and mcmapply
. This
approach is limited to unix-based platforms.
Please note that not all functions benfit from parallelisation. Often the
overhead to create/copy objects outrun the time saving of parallel runs. This
is true for functions that are very fast to compute (e.g.
sqrt
-transformation). That's why the default value for the
mc.cores
argument in all functions is 1L
.
It depends on the size of the dataset which step (often only
removeBaseline
and
detectPeaks
) benefits from parallelisation.
In general it is faster to encapsulate the complete workflow into a function
and parallelise it using mclapply
instead of using the
mc.cores
argument of each method. The reason is the reduced overhead
for object management (only one split/combine is needed instead of doing these
operations in each function again and again).
The following functions/methods support the mc.cores
argument:
## Not run: ## load package library("MALDIquant") ## load example data data("fiedler2009subset", package="MALDIquant") ## run single-core baseline correction print(system.time( b1 <- removeBaseline(fiedler2009subset, method="SNIP") )) if(.Platform$OS.type == "unix") { ## run multi-core baseline correction print(system.time( b2 <- removeBaseline(fiedler2009subset, method="SNIP", mc.cores=2) )) stopifnot(all.equal(b1, b2)) } ## parallelise complete workflow workflow <- function(spectra, cores) { s <- transformIntensity(spectra, method="sqrt", mc.cores=cores) s <- smoothIntensity(s, method="SavitzkyGolay", halfWindowSize=10, mc.cores=cores) s <- removeBaseline(s, method="SNIP", iterations=100, mc.cores=cores) s <- calibrateIntensity(s, method="TIC", mc.cores=cores) detectPeaks(s, method="MAD", halfWindowSize=20, SNR=2, mc.cores=cores) } if(.Platform$OS.type == "unix") { ## parallelise the complete workflow is often faster because the overhead is ## reduced print(system.time( p1 <- unlist(parallel::mclapply(fiedler2009subset, function(x)workflow(list(x), cores=1), mc.cores=2), use.names=FALSE) )) print(system.time( p2 <- workflow(fiedler2009subset, cores=2) )) stopifnot(all.equal(p1, p2)) } ## End(Not run)
## Not run: ## load package library("MALDIquant") ## load example data data("fiedler2009subset", package="MALDIquant") ## run single-core baseline correction print(system.time( b1 <- removeBaseline(fiedler2009subset, method="SNIP") )) if(.Platform$OS.type == "unix") { ## run multi-core baseline correction print(system.time( b2 <- removeBaseline(fiedler2009subset, method="SNIP", mc.cores=2) )) stopifnot(all.equal(b1, b2)) } ## parallelise complete workflow workflow <- function(spectra, cores) { s <- transformIntensity(spectra, method="sqrt", mc.cores=cores) s <- smoothIntensity(s, method="SavitzkyGolay", halfWindowSize=10, mc.cores=cores) s <- removeBaseline(s, method="SNIP", iterations=100, mc.cores=cores) s <- calibrateIntensity(s, method="TIC", mc.cores=cores) detectPeaks(s, method="MAD", halfWindowSize=20, SNR=2, mc.cores=cores) } if(.Platform$OS.type == "unix") { ## parallelise the complete workflow is often faster because the overhead is ## reduced print(system.time( p1 <- unlist(parallel::mclapply(fiedler2009subset, function(x)workflow(list(x), cores=1), mc.cores=2), use.names=FALSE) )) print(system.time( p2 <- workflow(fiedler2009subset, cores=2) )) stopifnot(all.equal(p1, p2)) } ## End(Not run)
MassPeaks
represents extracted peaks of a single
spectrum of a MALDI-TOF mass spectrometry measurement.
createMassPeaks
: Creates a
MassPeaks
object.
Class AbstractMassObject
, directly.
snr
:vector
, signal-to-noise
ratio
signature(x = "MassPeaks")
:
Draws peak labels to plot.
See labelPeaks,MassPeaks-method
for
details.
signature(x = "MassPeaks")
:
Finds monoisotopic peaks in peak lists.
See monoisotopicPeaks,MassPeaks-method
for
details.
Sebastian Gibb [email protected]
createMassPeaks
,
detectPeaks,MassSpectrum-method
,
labelPeaks,MassPeaks-method
,
AbstractMassObject
Website: https://strimmerlab.github.io/software/maldiquant/
MassSpectrum
represents a single spectrum of a MALDI-TOF
mass spectrometry measurement. It provides an easy framework for doing
some preprocessing steps like peak detection, baseline correction and
much more.
createMassSpectrum
: Creates a
MassSpectrum
object.
Class AbstractMassObject
, directly.
signature(x = "MassSpectrum")
:
Calibrates the intensity of a
MassSpectrum
object.
See calibrateIntensity,MassSpectrum-method
for
details.
signature(x = "MassSpectrum")
:
Look for local maxima and estimate noise to extract peaks out of a
MassSpectrum
object.
See detectPeaks,MassSpectrum-method
for
details.
signature(x = "MassSpectrum")
:
Estimates the baseline of a
MassSpectrum
object.
See estimateBaseline,MassSpectrum-method
for
details.
signature(x = "MassSpectrum")
:
Estimates the noise of a
MassSpectrum
object.
See estimateNoise,MassSpectrum-method
for
details.
signature(object = "MassSpectrum")
:
Returns FALSE
if the frequency of mass values with irregular
intervals is greater than threshold
(because object
was measured in centroid mode or some intensity
values were filtered).
signature(x = "MassSpectrum")
:
Estimates and removes the baseline of a
MassSpectrum
object.
See removeBaseline,MassSpectrum-method
for
details.
signature(object = "MassSpectrum")
:
Smoothes the intensities of an MassSpectrum
object.
See smoothIntensity,MassSpectrum-method
for
details.
signature(object = "MassSpectrum")
:
Accessor function for Total Ion Current (TIC, area under the curve).
Sebastian Gibb [email protected]
createMassSpectrum
,
calibrateIntensity,MassSpectrum-method
,
detectPeaks,MassSpectrum-method
,
estimateBaseline,MassSpectrum-method
,
estimateNoise,MassSpectrum-method
,
removeBaseline,MassSpectrum-method
,
smoothIntensity,MassSpectrum-method
,
AbstractMassObject
Website: https://strimmerlab.github.io/software/maldiquant/
## load package library("MALDIquant") ## create a MassSpectrum object by default constructor s <- createMassSpectrum(mass=1:100, intensity=rnorm(100)^2, metaData=list(name="example")) ## show some details s ## plot spectrum plot(s) ## get TIC totalIonCurrent(s) ## modify intensity and metaData intensity(s)[1:50] <- 0 metaData(s) <- list(name="modified example") ## plot again plot(s)
## load package library("MALDIquant") ## create a MassSpectrum object by default constructor s <- createMassSpectrum(mass=1:100, intensity=rnorm(100)^2, metaData=list(name="example")) ## show some details s ## plot spectrum plot(s) ## get TIC totalIonCurrent(s) ## modify intensity and metaData intensity(s)[1:50] <- 0 metaData(s) <- list(name="modified example") ## plot again plot(s)
match.closest
returns a vector of the positions of (first) matches
its first arguments in its second. In contrast to the similar
match
it just accept numeric
arguments but
has an additional tolerance
argument that allows relaxed
matching.
match.closest(x, table, tolerance = Inf, nomatch = NA_integer_)
match.closest(x, table, tolerance = Inf, nomatch = NA_integer_)
x |
|
table |
|
tolerance |
|
nomatch |
|
An integer
vector of the same length as x
giving the
closest position in table
of the first match or nomatch
if
there is no match.
library("MALDIquant") match.closest(c(1.1, 1.4, 9.8), 1:10) # [1] 1 1 10 match.closest(c(1.1, 1.4, 9.8), 1:10, tolerance=0.25) # [1] 1 NA 10 match.closest(c(1.1, 1.4, 9.8), 1:10, tolerance=0.25, nomatch=0) # [1] 1 0 10 ## this function is most useful if you want to subset an intensityMatrix ## by a few (reference) peaks ## create an example intensityMatrix im <- matrix(1:10, nrow=2, dimnames=list(NULL, 1:5)) attr(im, "mass") <- 1:5 im # 1 2 3 4 5 # [1,] 1 3 5 7 9 # [2,] 2 4 6 8 10 # attr(,"mass") # [1] 1 2 3 4 5 ## reference peaks ref <- c(2.2, 4.8) im[, match.closest(ref, attr(im, "mass"), tolerance=0.25, nomatch=0)] # 2 5 # [1,] 3 9 # [2,] 4 10
library("MALDIquant") match.closest(c(1.1, 1.4, 9.8), 1:10) # [1] 1 1 10 match.closest(c(1.1, 1.4, 9.8), 1:10, tolerance=0.25) # [1] 1 NA 10 match.closest(c(1.1, 1.4, 9.8), 1:10, tolerance=0.25, nomatch=0) # [1] 1 0 10 ## this function is most useful if you want to subset an intensityMatrix ## by a few (reference) peaks ## create an example intensityMatrix im <- matrix(1:10, nrow=2, dimnames=list(NULL, 1:5)) attr(im, "mass") <- 1:5 im # 1 2 3 4 5 # [1,] 1 3 5 7 9 # [2,] 2 4 6 8 10 # attr(,"mass") # [1] 1 2 3 4 5 ## reference peaks ref <- c(2.2, 4.8) im[, match.closest(ref, attr(im, "mass"), tolerance=0.25, nomatch=0)] # 2 5 # [1,] 3 9 # [2,] 4 10
MassPeaks
objects.
This function merges MassPeaks
objects.
mergeMassPeaks(l, labels, method=c("mean", "median", "sum"), ignore.na=TRUE, ...)
mergeMassPeaks(l, labels, method=c("mean", "median", "sum"), ignore.na=TRUE, ...)
l |
|
labels |
|
method |
used merge method. |
ignore.na |
Should |
... |
arguments to be passed to underlying functions (currently only
|
Returns a single (no labels
given) or a list
(labels
given) of merged MassPeaks
objects.
Sebastian Gibb [email protected]
Website: https://strimmerlab.github.io/software/maldiquant/
## load package library("MALDIquant") ## create four MassPeaks objects and add them to the list p <- list(createMassPeaks(mass=1:2, intensity=1:2), createMassPeaks(mass=1:3, intensity=1:3), createMassPeaks(mass=1:4, intensity=1:4), createMassPeaks(mass=1:5, intensity=1:5)) ## merge all four MassPeaks objects into a single new one ## by sum their intensities ## (no labels, returns only one new MassPeaks object) mergedPeaks <- mergeMassPeaks(p, method="sum") ## only merge MassPeaks objects in a group ## (two different labels, returns a list of two new MassPeaks objects) groups <- factor(c("a", "a", "b", "b"), levels=c("a", "b")) mergedPeaks <- mergeMassPeaks(p, labels=groups, method="mean") ## the same, but treat NA as zero mergedPeaks <- mergeMassPeaks(p, labels=groups, method="mean", ignore.na=FALSE)
## load package library("MALDIquant") ## create four MassPeaks objects and add them to the list p <- list(createMassPeaks(mass=1:2, intensity=1:2), createMassPeaks(mass=1:3, intensity=1:3), createMassPeaks(mass=1:4, intensity=1:4), createMassPeaks(mass=1:5, intensity=1:5)) ## merge all four MassPeaks objects into a single new one ## by sum their intensities ## (no labels, returns only one new MassPeaks object) mergedPeaks <- mergeMassPeaks(p, method="sum") ## only merge MassPeaks objects in a group ## (two different labels, returns a list of two new MassPeaks objects) groups <- factor(c("a", "a", "b", "b"), levels=c("a", "b")) mergedPeaks <- mergeMassPeaks(p, labels=groups, method="mean") ## the same, but treat NA as zero mergedPeaks <- mergeMassPeaks(p, labels=groups, method="mean", ignore.na=FALSE)
This method looks for monoisotopic peaks in peak list data
(represented by a MassPeaks
objects).
It is based on the poisson model for isotopic patterns described in Breen et al
2000.
## S4 method for signature 'MassPeaks' monoisotopicPeaks(object, minCor=0.95, tolerance=1e-4, distance=1.00235, size=3L:10L) ## S4 method for signature 'list' monoisotopicPeaks(object, ...)
## S4 method for signature 'MassPeaks' monoisotopicPeaks(object, minCor=0.95, tolerance=1e-4, distance=1.00235, size=3L:10L) ## S4 method for signature 'list' monoisotopicPeaks(object, ...)
object |
|
minCor |
|
tolerance |
|
distance |
|
size |
|
... |
arguments to be passed to
|
Returns a MassPeaks
object with monoisotopic peaks only.
Sebastian Gibb [email protected]
K. Park, J.Y. Yoon, S. Lee, E. Paek, H. Park, H.J. Jung, and S.W. Lee. 2008. Isotopic peak intensity ratio based algorithm for determination of isotopic clusters and monoisotopic masses of polypeptides from high-resolution mass spectrometric data. Analytical Chemistry, 80: 7294-7303.
E.J. Breen, F.G. Hopwood, K.L. Williams, and M.R. Wilkins. 2000. Automatic poisson peak harvesting for high throughput protein identification. Electrophoresis 21: 2243-2251.
MassPeaks
,
detectPeaks,MassSpectrum-method
Website: https://strimmerlab.github.io/software/maldiquant/
## load package library("MALDIquant") ## create example peaks p <- createMassPeaks(mass=995:1005, intensity=c(100, 10, 30, 10, 40, # noise 550, 330, 110, 10, # isotopic pattern 5, 15)) # more noise m <- monoisotopicPeaks(p) as.matrix(m) ## plot the peaks and mark the monoisotopic one plot(p) points(m, col=2, pch=4)
## load package library("MALDIquant") ## create example peaks p <- createMassPeaks(mass=995:1005, intensity=c(100, 10, 30, 10, 40, # noise 550, 330, 110, 10, # isotopic pattern 5, 15)) # more noise m <- monoisotopicPeaks(p) as.matrix(m) ## plot the peaks and mark the monoisotopic one plot(p) points(m, col=2, pch=4)
This function turns a mass spectrometry imaging dataset represented by a
list
of AbstractMassObject
objects into an
intensityMatrix
for each slice (stored in an
array
).
msiSlices(x, center, tolerance, method=c("sum", "mean", "median"), adjust=TRUE)
msiSlices(x, center, tolerance, method=c("sum", "mean", "median"), adjust=TRUE)
x |
a |
center |
|
tolerance |
|
method |
used aggregation function. |
adjust |
|
Each MassSpectrum
/MassPeaks
object in
x
must contain a list
named imaging
with an element
pos
that stores the x
and y
value of the spectrum, e.g.:
> metaData(spectra[[1]])$imaging$pos x y 1 5
Returns an array
of three dimensions. The first and second
dimensions contains the x and y coordinates of the image. The third dimension
represents the index of the center
of each slice. There are two
additional attributes, namely "center"
and "tolerance"
which store
the original center
and tolerance
information.
Sebastian Gibb [email protected]
AbstractMassObject
,
MassSpectrum
,
MassPeaks
,
coordinates,AbstractMassObject-method
,
plotMsiSlice,list-method
Please find real examples on:
Website: https://strimmerlab.github.io/software/maldiquant/
Vignette: https://github.com/sgibb/MALDIquantExamples/raw/master/inst/doc/nyakas2013.pdf
Shiny: https://github.com/sgibb/ims-shiny/
## load package library("MALDIquant") ## load example data data("fiedler2009subset", package="MALDIquant") ## please note: this is NOT a MSI data set ## we just add some coordinates for demonstration coordinates(fiedler2009subset) <- cbind(x=rep(1:4, 2), y=rep(1:2, each=4)) slices <- msiSlices(fiedler2009subset, center=c(5864.49, 8936.97), tolerance=0.25) slices
## load package library("MALDIquant") ## load example data data("fiedler2009subset", package="MALDIquant") ## please note: this is NOT a MSI data set ## we just add some coordinates for demonstration coordinates(fiedler2009subset) <- cbind(x=rep(1:4, 2), y=rep(1:2, each=4)) slices <- msiSlices(fiedler2009subset, center=c(5864.49, 8936.97), tolerance=0.25) slices
This is an overloaded method to allow plotting of an
AbstractMassObject
object.
## S4 method for signature 'AbstractMassObject,missing' plot(x, col="black", xlab=expression(italic(m/z)), ylab="intensity", type=ifelse(isMassPeaks(x), "h", "l"), xlim=c(ifelse(length(x@mass), min(x@mass, na.rm=TRUE), 0), ifelse(length(x@mass), max(x@mass, na.rm=TRUE), 1)), ylim=c(0, ifelse(length(x@intensity), max(x@intensity, na.rm=TRUE), 1)), main=x@metaData$name, sub=x@metaData$file, cex.sub=0.75, col.sub="#808080", ...)
## S4 method for signature 'AbstractMassObject,missing' plot(x, col="black", xlab=expression(italic(m/z)), ylab="intensity", type=ifelse(isMassPeaks(x), "h", "l"), xlim=c(ifelse(length(x@mass), min(x@mass, na.rm=TRUE), 0), ifelse(length(x@mass), max(x@mass, na.rm=TRUE), 1)), ylim=c(0, ifelse(length(x@intensity), max(x@intensity, na.rm=TRUE), 1)), main=x@metaData$name, sub=x@metaData$file, cex.sub=0.75, col.sub="#808080", ...)
x |
|
col |
line colour, see |
xlab |
title for the x-axis, see |
ylab |
title for the y-axis, see |
type |
type of plot: see |
xlim |
the x limits (x1, x2) of the plot, see
|
ylim |
the y limits (y1, y2) of the plot, see
|
main |
title for the plot, see |
sub |
sub title for the plot, see |
cex.sub |
sub title font size, see |
col.sub |
sub title color, see |
... |
arguments to be passed to |
Sebastian Gibb [email protected]
Website: https://strimmerlab.github.io/software/maldiquant/
## load package library("MALDIquant") ## create a MassSpectrum object by default constructor s <- createMassSpectrum(mass=1:100, intensity=rnorm(100)^2, metaData=list(name="example")) ## show some details s ## plot spectrum plot(s)
## load package library("MALDIquant") ## create a MassSpectrum object by default constructor s <- createMassSpectrum(mass=1:100, intensity=rnorm(100)^2, metaData=list(name="example")) ## show some details s ## plot spectrum plot(s)
This function allows to plot a slice of a mass spectrometry imaging dataset
represented by a list
of AbstractMassObject
objects
or an array
or a matrix
.
## S4 method for signature 'list' plotMsiSlice(x, center, tolerance, colRamp=colorRamp(c("black", "blue", "green", "yellow", "red")), interpolate=FALSE, legend=TRUE, alignLabels=FALSE, combine=FALSE, ...) ## S4 method for signature 'array' plotMsiSlice(x, colRamp=colorRamp(c("black", "blue", "green", "yellow", "red")), interpolate=FALSE, legend=TRUE, alignLabels=FALSE, combine=FALSE, plotInteractive=FALSE, ...) ## S4 method for signature 'matrix' plotMsiSlice(x, colRamp=colorRamp(c("black", "blue", "green", "yellow", "red")), interpolate=FALSE, scale=TRUE, legend=scale, ...)
## S4 method for signature 'list' plotMsiSlice(x, center, tolerance, colRamp=colorRamp(c("black", "blue", "green", "yellow", "red")), interpolate=FALSE, legend=TRUE, alignLabels=FALSE, combine=FALSE, ...) ## S4 method for signature 'array' plotMsiSlice(x, colRamp=colorRamp(c("black", "blue", "green", "yellow", "red")), interpolate=FALSE, legend=TRUE, alignLabels=FALSE, combine=FALSE, plotInteractive=FALSE, ...) ## S4 method for signature 'matrix' plotMsiSlice(x, colRamp=colorRamp(c("black", "blue", "green", "yellow", "red")), interpolate=FALSE, scale=TRUE, legend=scale, ...)
x |
The mass spectrometry imaging dataset. It could be a |
center |
|
tolerance |
|
colRamp |
colours as |
interpolate |
|
scale |
|
legend |
|
alignLabels |
|
combine |
|
plotInteractive |
|
... |
arguments to be passed to |
Each MassSpectrum
/MassPeaks
object in
x
must contain a list
named imaging
with an element
pos
that stores the x
and y
value of the spectrum, e.g.:
> metaData(spectra[[1]])$imaging$pos x y 1 5
Sebastian Gibb [email protected]
AbstractMassObject
,
MassSpectrum
,
MassPeaks
,
coordinates,AbstractMassObject-method
,
msiSlices
,
plot,MassSpectrum,missing-method
Please find real examples on:
Website: https://strimmerlab.github.io/software/maldiquant/
Vignette: https://github.com/sgibb/MALDIquantExamples/raw/master/inst/doc/nyakas2013.pdf
Shiny: https://github.com/sgibb/ims-shiny/
## load package library("MALDIquant") ## load example data data("fiedler2009subset", package="MALDIquant") ## please note: this is NOT a MSI data set ## we just add some coordinates for demonstration coordinates(fiedler2009subset) <- cbind(x=rep(1:4, 2), y=rep(1:2, each=4)) plotMsiSlice(fiedler2009subset, center=8936.97, tolerance=0.25) plotMsiSlice(fiedler2009subset, center=c(5864.49, 8936.97), tolerance=0.25, combine=TRUE, colRamp=list(colorRamp(c("#000000", "#FF00FF")), colorRamp(c("#000000", "#00FF00"))))
## load package library("MALDIquant") ## load example data data("fiedler2009subset", package="MALDIquant") ## please note: this is NOT a MSI data set ## we just add some coordinates for demonstration coordinates(fiedler2009subset) <- cbind(x=rep(1:4, 2), y=rep(1:2, each=4)) plotMsiSlice(fiedler2009subset, center=8936.97, tolerance=0.25) plotMsiSlice(fiedler2009subset, center=c(5864.49, 8936.97), tolerance=0.25, combine=TRUE, colRamp=list(colorRamp(c("#000000", "#FF00FF")), colorRamp(c("#000000", "#00FF00"))))
MassPeaks
object.
This function creates a reference MassPeaks
object (also
called Anchor Peaks) from a list of MassPeaks
objects.
Generally it is a combination of binPeaks
and
filterPeaks
referencePeaks(l, method=c("strict", "relaxed"), minFrequency=0.9, tolerance=0.002)
referencePeaks(l, method=c("strict", "relaxed"), minFrequency=0.9, tolerance=0.002)
l |
|
method |
bin creation rule (see |
minFrequency |
|
tolerance |
|
Returns a new MassPeaks
objects.
The intensity
slot of the returned
MassPeaks
represents the frequency of this mass position
in all samples.
Sebastian Gibb [email protected]
binPeaks
,
filterPeaks
,
MassPeaks
Website: https://strimmerlab.github.io/software/maldiquant/
## load package library("MALDIquant") ## create four MassPeaks objects and add them to the list p<- list(createMassPeaks(mass=1:2, intensity=1:2), createMassPeaks(mass=1:3, intensity=1:3), createMassPeaks(mass=1:4, intensity=1:4), createMassPeaks(mass=1:5, intensity=1:5)) ## only use peaks which occur in all MassPeaks objects as reference peaks refPeaks <- referencePeaks(p, minFrequency=1) mass(refPeaks) # 1:2 intensity(refPeaks) # c(1, 1)
## load package library("MALDIquant") ## create four MassPeaks objects and add them to the list p<- list(createMassPeaks(mass=1:2, intensity=1:2), createMassPeaks(mass=1:3, intensity=1:3), createMassPeaks(mass=1:4, intensity=1:4), createMassPeaks(mass=1:5, intensity=1:5)) ## only use peaks which occur in all MassPeaks objects as reference peaks refPeaks <- referencePeaks(p, minFrequency=1) mass(refPeaks) # 1:2 intensity(refPeaks) # c(1, 1)
This method removes the baseline of mass spectrometry data
(represented by a MassSpectrum
object).
The intensity of the mass spectrometry data would be reduced by
baseline
.
## S4 method for signature 'MassSpectrum' removeBaseline(object, method=c("SNIP", "TopHat", "ConvexHull", "median"), ...) ## S4 method for signature 'list' removeBaseline(object, ...)
## S4 method for signature 'MassSpectrum' removeBaseline(object, method=c("SNIP", "TopHat", "ConvexHull", "median"), ...) ## S4 method for signature 'list' removeBaseline(object, ...)
object |
|
method |
used baseline estimation method, one of
|
... |
arguments to be passed to
|
Returns a modified MassSpectrum
object with reduced
intensities.
Sebastian Gibb [email protected]
MassSpectrum
,
estimateBaseline,MassSpectrum-method
demo("baseline")
Website: https://strimmerlab.github.io/software/maldiquant/
## load package library("MALDIquant") ## load example data data("fiedler2009subset", package="MALDIquant") ## choose only the first mass spectrum s <- fiedler2009subset[[1]] ## plot spectrum plot(s) ## subtract baseline b <- removeBaseline(s, method="SNIP") ## draw modified spectrum on the plot lines(b, col="blue")
## load package library("MALDIquant") ## load example data data("fiedler2009subset", package="MALDIquant") ## choose only the first mass spectrum s <- fiedler2009subset[[1]] ## plot spectrum plot(s) ## subtract baseline b <- removeBaseline(s, method="SNIP") ## draw modified spectrum on the plot lines(b, col="blue")
This method smoothes the intensity values of a
MassSpectrum
object.
## S4 method for signature 'MassSpectrum' smoothIntensity(object, method=c("SavitzkyGolay", "MovingAverage"), halfWindowSize, ...)
## S4 method for signature 'MassSpectrum' smoothIntensity(object, method=c("SavitzkyGolay", "MovingAverage"), halfWindowSize, ...)
object |
|
method |
used smoothing method, one of |
halfWindowSize |
half window size. The resulting
window reaches from |
... |
arguments to be passed to |
halfWindowSize
: Depends on the selected method
.
For the SavitzkyGolay
the halfWindowSize
should be smaller than
FWHM of the peaks (full width at half maximum; please find details in
Bromba and Ziegler 1981).
In general the halfWindowSize
for the MovingAverage
has to be
much smaller than for SavitzkyGolay
to conserve the peak shape.
Sebastian Gibb [email protected]
Weighted moving average: Sigurdur Smarason
A. Savitzky and M. J. Golay. 1964. Smoothing and differentiation of data by simplified least squares procedures. Analytical chemistry, 36(8), 1627-1639.
M. U. Bromba and H. Ziegler. 1981. Application hints for Savitzky-Golay digital smoothing filters. Analytical Chemistry, 53(11), 1583-1586.
Website: https://strimmerlab.github.io/software/maldiquant/
## load package library("MALDIquant") ## load example data data("fiedler2009subset", package="MALDIquant") ## smooth spectra s <- smoothIntensity(fiedler2009subset, method="MovingAverage", halfWindowSize=2) ## or s <- smoothIntensity(fiedler2009subset, method="MovingAverage", halfWindowSize=2, weighted=TRUE) ## or s <- smoothIntensity(fiedler2009subset, method="SavitzkyGolay", halfWindowSize=10)
## load package library("MALDIquant") ## load example data data("fiedler2009subset", package="MALDIquant") ## smooth spectra s <- smoothIntensity(fiedler2009subset, method="MovingAverage", halfWindowSize=2) ## or s <- smoothIntensity(fiedler2009subset, method="MovingAverage", halfWindowSize=2, weighted=TRUE) ## or s <- smoothIntensity(fiedler2009subset, method="SavitzkyGolay", halfWindowSize=10)
This method performs a transformation (e.g. sqrt-transformation) on the
intensities of an AbstractMassObject
object.
## S4 method for signature 'AbstractMassObject' transformIntensity(object, method=c("sqrt", "log", "log2", "log10")) ## S4 method for signature 'list' transformIntensity(object, ...)
## S4 method for signature 'AbstractMassObject' transformIntensity(object, method=c("sqrt", "log", "log2", "log10")) ## S4 method for signature 'list' transformIntensity(object, ...)
object |
|
method |
used transformation method. |
... |
arguments to be passed to underlying functions. If
|
Sebastian Gibb [email protected]
AbstractMassObject
,
MassSpectrum
Website: https://strimmerlab.github.io/software/maldiquant/
## load package library("MALDIquant") ## load example data data("fiedler2009subset", package="MALDIquant") ## choose only the first mass spectrum s <- fiedler2009subset[[1]] ## transform spectrum t <- transformIntensity(s, method="sqrt") ## plot spectrum par(mfrow=c(2, 1)) plot(s, main="raw spectrum") plot(t, main="transformed spectrum") par(mfrow=c(1, 1))
## load package library("MALDIquant") ## load example data data("fiedler2009subset", package="MALDIquant") ## choose only the first mass spectrum s <- fiedler2009subset[[1]] ## transform spectrum t <- transformIntensity(s, method="sqrt") ## plot spectrum par(mfrow=c(2, 1)) plot(s, main="raw spectrum") plot(t, main="transformed spectrum") par(mfrow=c(1, 1))
This method trims an AbstractMassObject
object.
That is useful if some mass ranges should be excluded from further analysis.
## S4 method for signature 'AbstractMassObject,numeric' trim(object, range) ## S4 method for signature 'list,numeric' trim(object, range, ...) ## S4 method for signature 'list,missing' trim(object, range, ...)
## S4 method for signature 'AbstractMassObject,numeric' trim(object, range) ## S4 method for signature 'list,numeric' trim(object, range, ...) ## S4 method for signature 'list,missing' trim(object, range, ...)
object |
|
range |
|
... |
arguments to be passed to underlying functions (currently only
|
Sebastian Gibb [email protected]
AbstractMassObject
,
MassPeaks
,
MassSpectrum
Website: https://strimmerlab.github.io/software/maldiquant/
## load package library("MALDIquant") ## load example data data("fiedler2009subset", package="MALDIquant") ## select only one spectrum s <- fiedler2009subset[[1]] ## remove all mass lower 3000 trim(s, range=c(3000, Inf)) ## remove all mass higher 8000 trim(s, range=c(0, 8000)) ## remove all mass lower 3000 and higher 8000 trim(s, range=c(3000, 8000)) ## choose largest overlapping mass range for all spectra trim(fiedler2009subset)
## load package library("MALDIquant") ## load example data data("fiedler2009subset", package="MALDIquant") ## select only one spectrum s <- fiedler2009subset[[1]] ## remove all mass lower 3000 trim(s, range=c(3000, Inf)) ## remove all mass higher 8000 trim(s, range=c(0, 8000)) ## remove all mass lower 3000 and higher 8000 trim(s, range=c(3000, 8000)) ## choose largest overlapping mass range for all spectra trim(fiedler2009subset)
These functions run warping functions on
AbstractMassObject
objects
(warping is also known as phase correction).
warpMassPeaks(l, w, emptyNoMatches=FALSE) warpMassSpectra(l, w, emptyNoMatches=FALSE)
warpMassPeaks(l, w, emptyNoMatches=FALSE) warpMassSpectra(l, w, emptyNoMatches=FALSE)
l |
|
w |
a |
emptyNoMatches |
|
The warping function w
is called in the following way:
Returns a list
of warped MassPeaks
or
MassSpectrum
objects.
Sebastian Gibb [email protected]
determineWarpingFunctions
,
MassPeaks
,
MassSpectrum
Website: https://strimmerlab.github.io/software/maldiquant/
## load package library("MALDIquant") ## create a MassPeaks object p <- createMassPeaks(mass=1:5, intensity=1:5) ## stupid warping function for demonstration ## (please use determineWarpingFunctions in real life applications) simpleWarp <- function(x) { return(1) } ## run warping function w <- warpMassPeaks(list(p), list(simpleWarp))[[1]] ## compare results all(mass(w) == mass(p)+1) # TRUE ## no warping (MassPeaks object is not changed) warpMassPeaks(list(p), list(NA)) ## no warping (intensity values of MassPeaks object are set to zero) warpMassPeaks(list(p), list(NA), emptyNoMatches=TRUE)
## load package library("MALDIquant") ## create a MassPeaks object p <- createMassPeaks(mass=1:5, intensity=1:5) ## stupid warping function for demonstration ## (please use determineWarpingFunctions in real life applications) simpleWarp <- function(x) { return(1) } ## run warping function w <- warpMassPeaks(list(p), list(simpleWarp))[[1]] ## compare results all(mass(w) == mass(p)+1) # TRUE ## no warping (MassPeaks object is not changed) warpMassPeaks(list(p), list(NA)) ## no warping (intensity values of MassPeaks object are set to zero) warpMassPeaks(list(p), list(NA), emptyNoMatches=TRUE)