Title: | Isotope Pattern, Profile and Centroid Calculation for Mass Spectrometry |
---|---|
Description: | Fast and very memory-efficient calculation of isotope patterns, subsequent convolution to theoretical envelopes (profiles) plus valley detection and centroidization or intensoid calculation. Batch processing, resolution interpolation, wrapper, adduct calculations and molecular formula parsing. Loos, M., Gerber, C., Corona, F., Hollender, J., Singer, H. (2015) <doi:10.1021/acs.analchem.5b00941>. |
Authors: | Martin Loos, Christian Gerber |
Maintainer: | Martin Loos <[email protected]> |
License: | GPL-2 |
Version: | 2.7 |
Built: | 2024-11-03 05:35:37 UTC |
Source: | https://github.com/blosloos/envipat |
Fast and memory-efficient calculation of isotope patterns (fine structures) for up to very large molecules, based on three different algorithms. Subsequent convolution of isotope patterns with a peak shape function to theoretical envelopes (profiles). Based on envelopes, valley detection and centroidization/intensoid calculation. Allows for batch processing of chemical formulas and interpolation of measurement resolutions. Includes a wrapper combining all of the above functionalities.
Furthermore, includes (1) a check for consistency of chemical formulas, (2) a check for molecules with overlapping isotope patterns, (3) a list of all stable isotopes, (4) a list of different resolution data sets for Thermo Orbitrap and QExactive high-resolution mass spectrometers and (5) a list of adducts formed during electorspray ionization (ESI).
A web-based GUI for enviPat is freely available under https://www.envipat.eawag.ch.
Package: | enviPat |
Type: | Package |
Version: | 1.0 |
Date: | 2013-03-05 |
License: | GPL-2 |
LazyLoad: | yes |
Martin Loos, Christian Gerber
Maintainer: Martin Loos <[email protected]>
Loos, M., Gerber, C., Corona, F., Hollender, J., Singer, H. (2015). Accelerated isotope fine structure calculation using pruned transition trees, Analytical Chemistry 87(11), 5738-5744.
https://pubs.acs.org/doi/abs/10.1021/acs.analchem.5b00941
check_chemform
getR
isopattern
envelope
vdetect
isowrap
check_several
isotopes
resolution_list
chemforms
adducts
check_ded
mergeform
subform
multiform
List of common adducts observed for ESI-MS measurements in soft positive and negative ionization modes.
data(adducts)
data(adducts)
A data frame with 47 observations on the following 6 variables.
Name
Adduct name
calc
Equation for calculating adduct m/z from uncharged non-adduct molecular mass M (m/z = M/z + X)
Charge
z
Mult
1/z
Mass
X
Ion_mode
Ionization mode (positive or negative)
Formula_add
Adduct chemical formula to be added
Formula_ded
Adduct chemical formula to be subtracted
Multi
Factor to multiply chemical formula with
The correct way to calculate the isotopic pattern of a specific adduct is the following.
First, multiply the chemical formula of the molecule by the times it appears in the final adduct; multiform
.
Second, add the chemical formula of any adduct to that of the molecule; mergeform
.
Third, subtract the chemical formula of any deduct from that of the molecule; check_ded
& subform
.
Finally, calculate the isotopic fine structure using the correct charge
argument in isopattern
.
Chemical formulas must conform to what is described in check_chemform
.
https://fiehnlab.ucdavis.edu/staff/kind/Metabolomics/MS-Adduct-Calculator/
Huang N., Siegel M.M., Kruppa G.H., Laukien F.H., J. Am. Soc. Mass. Spectrom. 1999, 10. Automation of a Fourier transform ion cyclotron resonance mass spectrometer for acquisition, analysis, and e-mailing of high-resolution exact-mass electrospray ionization mass spectral data
multiform
mergeform
check_ded
subform
# example of M+H adduct batch calculation data(adducts) data(isotopes) data(chemforms) # (1) check formulas for consistency - recommended checked_chemforms <- check_chemform(isotopes, chemforms) # (2) multiply, see column 4 of adducts chemforms <-multiform(checked_chemforms[,2],1) # (3) add adduct - see column 7 of adducts chemforms<-mergeform(chemforms,"H1") # (4) calculate fine structure patterns <- isopattern(isotopes, chemforms)
# example of M+H adduct batch calculation data(adducts) data(isotopes) data(chemforms) # (1) check formulas for consistency - recommended checked_chemforms <- check_chemform(isotopes, chemforms) # (2) multiply, see column 4 of adducts chemforms <-multiform(checked_chemforms[,2],1) # (3) add adduct - see column 7 of adducts chemforms<-mergeform(chemforms,"H1") # (4) calculate fine structure patterns <- isopattern(isotopes, chemforms)
Checks chemical formulas (=a vector of character strings) for consistency with usage in isopattern
;
calculation of the molecular mass.
check_chemform(isotopes,chemforms,get_sorted=FALSE,get_list=FALSE)
check_chemform(isotopes,chemforms,get_sorted=FALSE,get_list=FALSE)
isotopes |
|
chemforms |
Vector of character strings with chemical formulas |
get_sorted |
Should elements in each formula be sorted according to their order in |
get_list |
Return list with vectors of elementwise atom counts contained in each chemical formula? |
Default checks if (1) a chemical formula contains only letters, numbers and square or round brackets, (2) elements can be found
in isotopes
and (3) letters and round brackets are all followed by a number of counts.
Where (3) are missing, they are set to 1.
(2) must consist of an upper case letter, possibly followed by lower case letters; to refer to individual isotopes (e.g., from isotope labelling of a molecule, e.g., N5 vs. [15]N2N3), square brackets may precede the capital letter. Any other symbols which may be part of a chemical formula (e.g., charges (+), dashes, asterisks, ...) are not permissible.
The molecular mass will be calculated from isotope masses and abundances listed in isotopes
.
Dataframe with 3 columns for get_list=FALSE
:
warning |
Correct chemical formula, |
new_formula |
Chemical formula |
monoisotopic_mass |
Monoisotopic mass |
Or list containing vector of elements for get_list=TRUE
.
Highly recommended for usage with isopattern
Martin Loos, Christian Gerber
# Check package data set of chemical formulas ############# data(chemforms); data(isotopes); checked<-check_chemform(isotopes,chemforms); checked; # Check for some senseless molecular formulas ############# chemforms<-c("C900Cl4H49","O82394","C8G500Zn9","Br1","6DBr9889"); data(isotopes); checked<-check_chemform(isotopes,chemforms); checked; # Molecular mass with and without isotope labelling ####### chemforms<-c("C10H5N4O5","[13]C2C8D2H3[15]N2N2[18]O2O3"); data(isotopes); checked<-check_chemform(isotopes,chemforms); checked;
# Check package data set of chemical formulas ############# data(chemforms); data(isotopes); checked<-check_chemform(isotopes,chemforms); checked; # Check for some senseless molecular formulas ############# chemforms<-c("C900Cl4H49","O82394","C8G500Zn9","Br1","6DBr9889"); data(isotopes); checked<-check_chemform(isotopes,chemforms); checked; # Molecular mass with and without isotope labelling ####### chemforms<-c("C10H5N4O5","[13]C2C8D2H3[15]N2N2[18]O2O3"); data(isotopes); checked<-check_chemform(isotopes,chemforms); checked;
Check if a chemical formula is contained in another chemical formula
check_ded(formulas, deduct)
check_ded(formulas, deduct)
formulas |
Vector with the containing chemical formula(s) |
deduct |
Chemical formula to be contained ("deduct") |
Returns a vector with length of input formulas
, with TRUE
if deduct is not contained and FALSE
otherwise.
Might be used used prior to subtracting a "deduct" chemical formula from that of a molecule when including adducts in the calculation of
isotopic patterns. Chemical formulas must conform to what is described in check_chemform
.
Martin Loos
formulas<-c("C8H4Cl2","C10H16O2","C3H10") deduct<-c("C4H10") check_ded(formulas, deduct)
formulas<-c("C8H4Cl2","C10H16O2","C3H10") deduct<-c("C4H10") check_ded(formulas, deduct)
Check for molecules overlapping in m/z,
based on isotope fine structures from isopattern
or on centroids/intensoids from envelope
.
check_several(pattern, dmz, ppm = TRUE)
check_several(pattern, dmz, ppm = TRUE)
pattern |
Output from |
dmz |
m/z window. In combination with |
ppm |
Should m/z window be set in ppm ( |
Overlaps in m/z among molecules are screened for within the m/z tolerance defined by the arguments dmz
and ppm
.
Dataframe with 4 columns, with number of rows equal to the length of argument pattern
compound |
Chemical formula of the compound |
warning |
Overlap detected? |
to? |
If overlap: with wich other compound(s)? Refers to row number, recycled for peak_number. |
peak_number |
If overlap: with which peak(s) of the other compound(s)? Refers to peak number. |
Martin Loos, Christian Gerber
data(isotopes) data(chemforms) pattern<-isopattern( isotopes, chemforms, threshold=0.1, plotit=TRUE, charge=FALSE, emass=0.00054858, algo=1 ) check_several(pattern,dmz=0.001,ppm=FALSE)
data(isotopes) data(chemforms) pattern<-isopattern( isotopes, chemforms, threshold=0.1, plotit=TRUE, charge=FALSE, emass=0.00054858, algo=1 ) check_several(pattern,dmz=0.001,ppm=FALSE)
Vector with character strings of exemplary chemical formulas (pesticides, pharmaceuticals)
data(chemforms)
data(chemforms)
Vector with character strings
data(chemforms) chemforms
data(chemforms) chemforms
Convolutes an isotope pattern from isopattern
with a peak shape function (Gaussian or Cauchy-Lorentz function)
to its theoretical envelope (profile), at a given measurement resolution.
The envelope is represented by sticks, i.e. measurement abundances at discrete m/z intervals.
envelope(pattern, ppm = FALSE, dmz = "get", frac = 1/4, env = "Gaussian", resolution = 5e+05, plotit = FALSE, verbose = TRUE)
envelope(pattern, ppm = FALSE, dmz = "get", frac = 1/4, env = "Gaussian", resolution = 5e+05, plotit = FALSE, verbose = TRUE)
pattern |
List of isotope pattern(s) as generated by |
ppm |
Should stick discretization be set in ppm ( |
dmz |
Stick discretization. Set to |
frac |
Used if |
env |
Peak shape function; either |
resolution |
Single resolution value or vector of resolutions with length equal to the number of entries in list |
plotit |
Should results be plotted, |
verbose |
Verbose, |
The theoretical profiles are represented by sticks, i.e. abundances at discrete m/z intervals.
While the profile width is set by argument resolution
, the mass discretization between adjacent sticks can be set in two different ways.
On the one hand, discretization can be given as a numerical value, either in ppm or absolute m/z.
To do so, set argument dmz
to a numerical value and specify with argument ppm
if this value is stating the discretization in ppm or
as absolute m/z.
On the other hand, discretization can be derived from the measurement resolution (R) set by
argument resolution
. To do so, set dmz
to "get"
, which leads to argument ppm
being ignored.
In this case, the stick discretization is retrieved from (dm/z)*frac
, with (dm/z) = (m/z)/R = peak width at half maximum.
List with length equal to length of list pattern
, with equal names of list entries.
Each entry in that list contains the sticks of the envelope in two columns:
m/z |
Stick m/z |
abundance |
Stick abundance |
The resolution R is defined as R=(m/z)/(dm/z), with dm/z = peak width at half maximum, cp. resolution_list
.
Martin Loos, Christian Gerber
Li, L., Kresh, J., Karabacak, N., Cobb, J., Agar, J. and Hong, P. (2008). A Hierarchical Algorithm for Calculating the Isotopic Fine Structures of Molecules. Journal of the American Society for Mass Spectrometry, 19, 1867–1874.
############################ # batch of chemforms ####### data(isotopes) data(chemforms) chemforms<-chemforms[1:5] pattern<-isopattern( isotopes, chemforms, threshold=0.1, plotit=TRUE, charge=FALSE, emass=0.00054858, algo=2 ) profiles<-envelope( pattern, ppm=FALSE, dmz=0.0001, frac=1/4, env="Gaussian", resolution=1E6, plotit=TRUE ) ############################
############################ # batch of chemforms ####### data(isotopes) data(chemforms) chemforms<-chemforms[1:5] pattern<-isopattern( isotopes, chemforms, threshold=0.1, plotit=TRUE, charge=FALSE, emass=0.00054858, algo=2 ) profiles<-envelope( pattern, ppm=FALSE, dmz=0.0001, frac=1/4, env="Gaussian", resolution=1E6, plotit=TRUE ) ############################
Given a set of MS measurement resolutions (R) as a function of measurement mass (m/z), getR interpolates R for any given molecular mass(es)
calculated by check_chemform
using smooth.spline
.
getR(checked, resmass, nknots = 13, spar = 0.1, plotit = TRUE)
getR(checked, resmass, nknots = 13, spar = 0.1, plotit = TRUE)
checked |
Dataframe produced by |
resmass |
Dataframe with two columns, resolution and mass; such as the list entries in |
nknots |
Integer number of knots to use for the smoothing spline. Default = 6. See also |
spar |
Smoothing parameter, (0,1]. See also |
plotit |
Plot results, |
Vector with resolutions.
check_chemform
gives molecular masses (m/z) for z=+/-1 only.
If z>1 or z<-1 is required, molecular mass entries in argument checked have to be divided accordingly to be consistent.
Martin Loos, Christian Gerber
smooth.spline
check_chemform
resolution_list
data(resolution_list) resmass<-resolution_list[[4]] data(isotopes) data(chemforms) checked<-check_chemform(isotopes,chemforms) resolution<-getR(checked,resmass,nknots=13,spar=0.1,plotit=TRUE) # same for z=-2: checked<-check_chemform(isotopes,chemforms) checked[,3]<-(checked[,3]/abs(-2)) resolution<-getR(checked,resmass,nknots=13,spar=0.1,plotit=TRUE)
data(resolution_list) resmass<-resolution_list[[4]] data(isotopes) data(chemforms) checked<-check_chemform(isotopes,chemforms) resolution<-getR(checked,resmass,nknots=13,spar=0.1,plotit=TRUE) # same for z=-2: checked<-check_chemform(isotopes,chemforms) checked[,3]<-(checked[,3]/abs(-2)) resolution<-getR(checked,resmass,nknots=13,spar=0.1,plotit=TRUE)
The function calculates the isotopologues ("isotope fine structure") of a given chemical formula or a set of chemical formulas (batch calculation) with fast and memory efficient transition tree algorithms, which can handle relative pruning thresholds. Returns accurate masses, probabilities and isotopic compositions of individual isotopologues. The isotopes of elements can be defined by the user.
isopattern(isotopes, chemforms, threshold = 0.001, charge = FALSE, emass = 0.00054857990924, plotit = FALSE, algo=1, rel_to = 0, verbose = TRUE, return_iso_calc_amount = FALSE)
isopattern(isotopes, chemforms, threshold = 0.001, charge = FALSE, emass = 0.00054857990924, plotit = FALSE, algo=1, rel_to = 0, verbose = TRUE, return_iso_calc_amount = FALSE)
isotopes |
Dataframe listing all relevant isotopes, such as |
chemforms |
Vector with character strings of chemical formulas, such as data set |
threshold |
Probability below which isotope peaks can be omitted, as specified by argument |
charge |
z in m/z. Either a single integer or a vector of integers with length equal to that of argument |
emass |
Electrone mass; only relevant if |
plotit |
Should results be plotted, |
algo |
Which algorithm to use? Type |
rel_to |
Probability definition, numeric |
verbose |
Verbose, |
return_iso_calc_amount |
Ignore; number of intermediate isotopologues. |
Isotope pattern calculation can be done by chosing one of two algorithms, set by argument algo
. Both algorithms use
transition tree updates to derive the exact mass and probability of a new isotopologue from existing ones, by steps of single isotope replacements.
These transition tree approaches are memory-efficient and fast for a wide range of molecular formulas and are able to reproduce the isotope fine structure
of molecules. The latter must often be pruned during calculation, c.p. argument rel_to
.
algo=1
grows transition trees within element-wise sub-molecules, whereas algo==2
grows them in larger sub-molecules of two elements, if available.
The latter approach can be slightly more efficient for very large or very complex molecules. The sub-isotopologues within sub-molecules are finally combined to
the isotopologuees of the full molecule. In contrast, intermediate counts of sub-isotopologues instead of fine structures are returned for return_iso_calc_amount==TRUE
rel_to
offers 5 possibilities of how probabilities are defined and pruned, each affecting the threshold
argument differently.
Default option rel_to=0
prunes and returns probabilities relative to the most intense isotope peak;
threshold
states a percentage of the intensity of this latter peak.
Similarly, option rel_to=1
normalizes relative to the peak consisting of the most abundant isotopes for each element, which
is often the monoisotopic one.
Option rel_to=2
prunes and returns absolute probabilities ; threshold
is not a percentage but an abolute cutoff.
Options rel_to=3
and rel_to=4
prune relative to the most intense and "monoisotopic" peak, respectively.
Although threshold
is a percentage, both options return absolut probabilities .
List with length equal to length of vector chemforms
; names of entries in list = chemical formula in chemform.
Each entry in that list contains information on individual isotopologues (rows) with columns:
m/z |
First column; m/z of an isotope peak. |
abundance |
Second column; abundance of an isotope peak. Probabilities are set relative to the most abundant peak of the isotope pattern. |
12C , 13C , 1H , 2H , ...
|
Third to all other columns; atom counts of individual isotopes for an isotope peak. |
Too low values for threshold
may lead to unnecessary calculation of low probable isotope peaks - to the extent that not enough memory is available
for either of the two algorithms.
It is highly recommended to check argument chemforms
with check_chemform
prior to running
isopattern
; argument chemforms
must conform to chemical formulas as defined in check_chemform
.
Element names must be followed by numbers (atom counts of that element), i.e. C1H4 is a valid argument whereas CH4 is not.
Otherwise, numbers may only be used in square brackets to denote individual isotopes defined in the element name column of iso_list, such as [14]C or [18]O.
For example, [13]C2C35H67N1O13 is the molecular formula of erythromycin labeled at two C-positions with [13]C;
C37H67N1O13 is the molecular formula of the unlabeled compound.
For correct adduct isotope pattern calculations, please check adducts
.
Martin Loos, Christian Gerber
Loos, M., Gerber, C., Corona, F., Hollender, J., Singer, H. (2015). Accelerated isotope fine structure calculation using pruned transition trees, Analytical Chemistry 87(11), 5738-5744.
https://pubs.acs.org/doi/abs/10.1021/acs.analchem.5b00941
https://www.envipat.eawag.ch/index.php
isopattern
chemforms
check_chemform
getR
envelope
vdetect
check_several
############################ # batch of chemforms ####### data(isotopes) data(chemforms) pattern<-isopattern( isotopes, chemforms, threshold=0.1, plotit=TRUE, charge=FALSE, emass=0.00054858, algo=1 ) ############################ # Single chemical formula ## data(isotopes) pattern<-isopattern( isotopes, "C100H200S2Cl5", threshold=0.1, plotit=TRUE, charge=FALSE, emass=0.00054858, algo=1 ) ############################
############################ # batch of chemforms ####### data(isotopes) data(chemforms) pattern<-isopattern( isotopes, chemforms, threshold=0.1, plotit=TRUE, charge=FALSE, emass=0.00054858, algo=1 ) ############################ # Single chemical formula ## data(isotopes) pattern<-isopattern( isotopes, "C100H200S2Cl5", threshold=0.1, plotit=TRUE, charge=FALSE, emass=0.00054858, algo=1 ) ############################
Dataframe with stable isotopes.
data(isotopes)
data(isotopes)
A data frame with 302 observations on the following 4 variables.
element
Chemical element
isotope
Stable isotopes of an element
mass
Relative atomic mass
abundance
Isotopic composition of an element
ratioC
Maximum number of atoms of an element for one C-atom in a molecule, based on 99.99 % of case molecules.
The ratioC
-value stems from a database survey conducted by Kind&Fiehn (2007); to disable, set value to 0.
The list serves as input into several package nontarget-functions. The first column of the data frame also
contains names of specific isotopes used for labeled compounds.
https://physics.nist.gov/cgi-bin/Compositions/stand_alone.pl
Kind, T. and Fiehn, O., 2007. Seven golden rules for heuristic filtering of molecular formulas obtained by accurate mass spectrometry. BMC Bioinformatics, 8:105.
data(isotopes)
data(isotopes)
Wrapper combining the functions getR
, isopattern
, envelope
and
vdetect
.
Uses chemical formulas from check_chemform
as argument.
isowrap(isotopes, checked, resmass, resolution = FALSE, nknots = 6, spar = 0.2, threshold = 0.1, charge = 1, emass = 0.00054858, algo=2, ppm = FALSE, dmz = "get", frac = 1/4, env = "Gaussian", detect = "centroid", plotit = FALSE,verbose = TRUE )
isowrap(isotopes, checked, resmass, resolution = FALSE, nknots = 6, spar = 0.2, threshold = 0.1, charge = 1, emass = 0.00054858, algo=2, ppm = FALSE, dmz = "get", frac = 1/4, env = "Gaussian", detect = "centroid", plotit = FALSE,verbose = TRUE )
isotopes |
Dataframe listing all relevant isotopes, such as |
checked |
Output dataframe from |
resmass |
For resolution interpolation: dataframe with two columns, resolution and mass; see |
resolution |
Single resolution value. Only used if argument |
nknots |
Number of knots, see |
spar |
Smoothing parameter, see |
threshold |
Abundance below which isotope peaks are omitted, see |
charge |
z in m/z, see |
emass |
Electrone mass. Only relevant if |
algo |
Which algorithm to use? Type |
ppm |
Set stick discretization, see details section of |
dmz |
Set stick discretization, see details section of |
frac |
Set stick discretization, see details section of |
env |
Peak shape function, see |
detect |
Return either |
plotit |
Should results be plotted, |
verbose |
Verbose, |
List with length equal to length of list profiles
, with equal names of list entries.
Each entry in that list contains the centroids, intensoids or valley of the envelope in two columns:
m/z |
m/z |
abundance |
area(centroid) or abundance (intensoid, valley) |
Martin Loos, Christian Gerber
data(isotopes); data(resolution_list); data(chemforms); chemforms<-chemforms[1:10]; checked<-check_chemform( isotopes, chemforms ); resmass<-resolution_list[[1]] centro<-isowrap( isotopes, checked, resmass=resolution_list[[4]], resolution=FALSE, nknots=4, spar=0.2, threshold=0.1, charge=1, emass=0.00054858, algo=2, ppm=FALSE, dmz="get", # retrieve dm from R=m/dm frac=1/4, env="Gaussian", detect="centroid", plotit=TRUE )
data(isotopes); data(resolution_list); data(chemforms); chemforms<-chemforms[1:10]; checked<-check_chemform( isotopes, chemforms ); resmass<-resolution_list[[1]] centro<-isowrap( isotopes, checked, resmass=resolution_list[[4]], resolution=FALSE, nknots=4, spar=0.2, threshold=0.1, charge=1, emass=0.00054858, algo=2, ppm=FALSE, dmz="get", # retrieve dm from R=m/dm frac=1/4, env="Gaussian", detect="centroid", plotit=TRUE )
Combine chemical formulas
mergeform(formula1,formula2)
mergeform(formula1,formula2)
formula1 |
Vector of first chemical formula(s), character string(s) |
formula2 |
Second chemical formula, single character string |
Useful for adduct calculations, check adducts
.
Chemical formulas must conform to what is described in check_chemform
.
Merged chemical formula(s), character string
Martin Loos
formula1<-c("C10[13]C2H10Cl10") formula2<-c("C2H5Na1") mergeform(formula1,formula2)
formula1<-c("C10[13]C2H10Cl10") formula2<-c("C2H5Na1") mergeform(formula1,formula2)
Multiply all atom numbers in a chemical formula by a factor
multiform(formula_in,fact)
multiform(formula_in,fact)
formula_in |
Chemical formula to be multiplied, vector of character strings |
fact |
Factor to multiply with |
Useful for adduct calculations, check adducts
.
Chemical formulas must conform to what is described in check_chemform
.
Multiplied chemical formula, character string
Martin Loos
formula_in <- "C10[13]C2H10Cl10" multiform(formula_in, 3)
formula_in <- "C10[13]C2H10Cl10" multiform(formula_in, 3)
List of different resolutions R=f(m/z) for various (high-resolution) mass spectrometers. For each of the instruments, different resolution settings are available. Here, R is defined as R=(m/z)/(dm/z), with dm/z = peak width at half maximum. Serves as input to getR to interpolate R from given molecular masses.
data(resolution_list)
data(resolution_list)
The format is: List with 29 data sets: Instrument_(massRange_instrumentMode_slicerMode)_Resolution@m/z
Elite/R240000@400
Elite/R120000@400
Elite/R60000@400
Elite/R30000@400
OrbitrapXL,Velos,VelosPro/R120000@400
OrbitrapXL,Velos,VelosPro/R60000@400
OrbitrapXL,Velos,VelosPro/R30000@400
OrbitrapXL,Velos,VelosPro/R15000@400
OrbitrapXL,Velos,VelosPro/R7500@400
Q-Exactive,ExactivePlus/280K@200
Q-Exactive,ExactivePlus/R140000@200
Q-Exactive,ExactivePlus/R70000@200
Q-Exactive,ExactivePlus/R35000@200
Q-Exactive,ExactivePlus/R17500@200
Exactive/R100000@200
Exactive/R50000@200
Exactive/R25000@200
Exactive/R12500@200
OTFusion,QExactiveHF/480000@200
OTFusion,QExactiveHF/240000@200
OTFusion,QExactiveHF/120000@200
OTFusion,QExactiveHF/60000@200
OTFusion,QExactiveHF/30000@200
OTFusion,QExactiveHF/15000@200
QTOF_XevoG2-S/R25000@200
Sciex_TripleTOF5600_R25000@200
Sciex_TripleTOF6600_R25000@200
Sciex_QTOFX500R_R25000@200
Agilent_low_extended_highSens_QTOF6550_R25000@200
Data assembled from individual measurements.
data(resolution_list) resolution_list
data(resolution_list) resolution_list
Subtract one chemical formula from another
subform(formula1,formula2)
subform(formula1,formula2)
formula1 |
Chemical formula to subtract from |
formula2 |
Chemical formula to subtract |
Useful for adduct calculations, check adducts
.
Chemical formulas must conform to what is described in check_chemform
.
Prior check if formula2
is contained in formula2
at all? See check_ded
.
A unified and filtered peaklist
Martin Loos
formula1<-c("C10[13]C2H10Cl10") formula2<-c("C2H5[13]C1") subform(formula1,formula2)
formula1<-c("C10[13]C2H10Cl10") formula2<-c("C2H5[13]C1") subform(formula1,formula2)
Checks envelopes calculated by envelope
for valleys and extracts centroids or intensoids.
vdetect(profiles,detect="centroid",plotit=TRUE,verbose=TRUE)
vdetect(profiles,detect="centroid",plotit=TRUE,verbose=TRUE)
profiles |
List of stick profiles as generated by |
detect |
To return either |
plotit |
Should results be plotted, |
verbose |
Verbose, |
List with length equal to length of list profiles
, with equal names of list entries.
Each entry in that list contains the centroids, intensoids or valleys of the envelope in two columns:
m/z |
m/z |
abundance |
Area (centroid) or abundance (intensoid, valley) |
Valley: local profile minimum, i.e. any envelope stick flanked by two other sticks of higher abundance.
Stick: see envelope
.
Centroid mass: intensity-weighted sum of the m/z of sticks between two valleys.
Centroid intensity: profile area between two valleys (mean of upper and lower sum of stick intensities), normalized to the maximum centroid area of the envelope.
Intensoid mass: m/z of the most intense stick between two valleys.
Intensoid intensity: intensity of the most intensive stick between two valleys, normalized to the most intense intensoid.
Too low stick discretization leads to imprecision in valley, centroid and intensoid characteristics.
Martin Loos, Christian Gerber
############################ # batch of chemforms ####### data(isotopes) data(chemforms) chemforms<-chemforms[1:5] pattern<-isopattern( isotopes, chemforms, threshold=0.1, plotit=TRUE, charge=FALSE, emass=0.00054858, algo=2 ) profiles<-envelope( pattern, ppm=FALSE, dmz=0.0001, frac=1/4, env="Gaussian", resolution=1E6, plotit=TRUE ) centro<-vdetect( profiles, detect="centroid", plotit=TRUE ) ############################
############################ # batch of chemforms ####### data(isotopes) data(chemforms) chemforms<-chemforms[1:5] pattern<-isopattern( isotopes, chemforms, threshold=0.1, plotit=TRUE, charge=FALSE, emass=0.00054858, algo=2 ) profiles<-envelope( pattern, ppm=FALSE, dmz=0.0001, frac=1/4, env="Gaussian", resolution=1E6, plotit=TRUE ) centro<-vdetect( profiles, detect="centroid", plotit=TRUE ) ############################