OpenMS
CoarseIsotopePatternGenerator Class Reference

Isotope pattern generator for coarse isotope distributions. More...

#include <OpenMS/CHEMISTRY/ISOTOPEDISTRIBUTION/CoarseIsotopePatternGenerator.h>

Inheritance diagram for CoarseIsotopePatternGenerator:
[legend]
Collaboration diagram for CoarseIsotopePatternGenerator:
[legend]

Public Member Functions

 CoarseIsotopePatternGenerator (const Size max_isotope=0, const bool round_masses=false)
 
 ~CoarseIsotopePatternGenerator () override
 
IsotopeDistribution run (const EmpiricalFormula &) const override
 Creates an isotope distribution from an empirical sum formula. More...
 
IsotopeDistribution estimateFromPeptideWeight (double average_weight)
 Estimate Peptide Isotopedistribution from weight and number of isotopes that should be reported. More...
 
IsotopeDistribution estimateFromPeptideMonoWeight (double mono_weight)
 Estimate Peptide Isotopedistribution from monoisotopic weight and number of isotopes that should be reported. More...
 
IsotopeDistribution estimateFromPeptideWeightAndS (double average_weight, UInt S)
 Estimate peptide IsotopeDistribution from average weight and exact number of sulfurs. More...
 
IsotopeDistribution estimateFromRNAWeight (double average_weight)
 Estimate Nucleotide Isotopedistribution from weight and number of isotopes that should be reported. More...
 
IsotopeDistribution estimateFromRNAMonoWeight (double mono_weight)
 Estimate Nucleotide Isotopedistribution from monoisotopic weight and number of isotopes that should be reported. More...
 
IsotopeDistribution estimateFromDNAWeight (double average_weight)
 Estimate Nucleotide Isotopedistribution from weight and number of isotopes that should be reported averagine model from Zubarev, R. A.; Demirev, P. A. in "Isotope depletion of large biomolecules: Implications for molecular mass measurements.". More...
 
IsotopeDistribution estimateFromWeightAndComp (double average_weight, double C, double H, double N, double O, double S, double P)
 Estimate Isotopedistribution from weight, average composition, and number of isotopes that should be reported. More...
 
IsotopeDistribution estimateFromMonoWeightAndComp (double mono_weight, double C, double H, double N, double O, double S, double P)
 Estimate Isotopedistribution from monoisotopic weight, average composition, and number of isotopes that should be reported. More...
 
IsotopeDistribution estimateFromWeightAndCompAndS (double average_weight, UInt S, double C, double H, double N, double O, double P)
 Estimate IsotopeDistribution from weight, exact number of sulfurs, and average remaining composition. More...
 
IsotopeDistribution estimateForFragmentFromPeptideWeight (double average_weight_precursor, double average_weight_fragment, const std::set< UInt > &precursor_isotopes)
 Estimate peptide fragment IsotopeDistribution from the precursor's average weight, fragment's average weight, and a list of isolated precursor isotopes. More...
 
IsotopeDistribution estimateForFragmentFromPeptideWeightAndS (double average_weight_precursor, UInt S_precursor, double average_weight_fragment, UInt S_fragment, const std::set< UInt > &precursor_isotopes) const
 Estimate peptide fragment IsotopeDistribution from the precursor's average weight, number of sulfurs in the precursor, fragment's average weight, number of sulfurs in the fragment, and a list of isolated precursor isotopes. More...
 
IsotopeDistribution estimateForFragmentFromRNAWeight (double average_weight_precursor, double average_weight_fragment, const std::set< UInt > &precursor_isotopes)
 Estimate RNA fragment IsotopeDistribution from the precursor's average weight, fragment's average weight, and a list of isolated precursor isotopes. More...
 
IsotopeDistribution estimateForFragmentFromDNAWeight (double average_weight_precursor, double average_weight_fragment, const std::set< UInt > &precursor_isotopes)
 Estimate DNA fragment IsotopeDistribution from the precursor's average weight, fragment's average weight, and a list of isolated precursor isotopes. More...
 
IsotopeDistribution estimateForFragmentFromWeightAndComp (double average_weight_precursor, double average_weight_fragment, const std::set< UInt > &precursor_isotopes, double C, double H, double N, double O, double S, double P) const
 Estimate fragment IsotopeDistribution from the precursor's average weight, fragment's average weight, a list of isolated precursor isotopes, and average composition. More...
 
IsotopeDistribution calcFragmentIsotopeDist (const IsotopeDistribution &fragment_isotope_dist, const IsotopeDistribution &comp_fragment_isotope_dist, const std::set< UInt > &precursor_isotopes, const double fragment_mono_mass) const
 Calculate isotopic distribution for a fragment molecule. More...
 
CoarseIsotopePatternGeneratoroperator= (const CoarseIsotopePatternGenerator &iso)
 
IsotopeDistribution::ContainerType convolve (const IsotopeDistribution::ContainerType &left, const IsotopeDistribution::ContainerType &right) const
 convolves the distributions left and right and stores the result in result More...
 
Accessors
void setMaxIsotope (const Size &max_isotope)
 sets the maximal isotope with max_isotope More...
 
void setRoundMasses (const bool round_masses)
 sets the round_masses_ flag to round masses to integer values (true) or return accurate masses (false) More...
 
Size getMaxIsotope () const
 returns the currently set maximum isotope More...
 
bool getRoundMasses () const
 returns the current value of the flag to return expected masses (true) or atomic numbers (false). More...
 
- Public Member Functions inherited from IsotopePatternGenerator
 IsotopePatternGenerator ()
 
 IsotopePatternGenerator (double probability_cutoff)
 
virtual ~IsotopePatternGenerator ()
 

Static Public Member Functions

static IsotopeDistribution approximateFromPeptideWeight (double mass, UInt num_peaks=20, UInt charge=1)
 roughly approximate peptide IsotopeDistribution from monoisotopic weight using Poisson distribution. m/z values approximated by adding one neutron mass (divided by charge) for every peak, starting at the given monoisotopic weight. Foundation from: Bellew et al, https://dx.doi.org/10.1093/bioinformatics/btl276 More...
 
static std::vector< double > approximateIntensities (double mass, UInt num_peaks=20)
 roughly approximate intensity distribution of peptidic isotope patterns from monoisotopic weight using Poisson distribution. Foundation from: Bellew et al, https://dx.doi.org/10.1093/bioinformatics/btl276 More...
 

Protected Member Functions

IsotopeDistribution::ContainerType convolvePow_ (const IsotopeDistribution::ContainerType &input, Size factor) const
 convolves the distribution input factor times and stores the result in result More...
 
IsotopeDistribution::ContainerType convolveSquare_ (const IsotopeDistribution::ContainerType &input) const
 convolves the distribution input with itself and stores the result in result More...
 
IsotopeDistribution::ContainerType correctMass_ (const IsotopeDistribution::ContainerType &input, const double mono_weight) const
 converts the masses of distribution input from atomic numbers to accurate masses More...
 
IsotopeDistribution calcFragmentIsotopeDist_ (const IsotopeDistribution::ContainerType &fragment_isotope_dist, const IsotopeDistribution::ContainerType &comp_fragment_isotope_dist, const std::set< UInt > &precursor_isotopes) const
 calculates the fragment distribution for a fragment molecule and stores it in result. More...
 
IsotopeDistribution::ContainerType fillGaps_ (const IsotopeDistribution::ContainerType &id) const
 fill a gapped isotope pattern (i.e. certain masses are missing), with zero probability masses More...
 

Protected Attributes

Size max_isotope_
 maximal isotopes which is used to calculate the distribution More...
 
bool round_masses_
 flag to determine whether masses should be rounded or not More...
 
- Protected Attributes inherited from IsotopePatternGenerator
double min_prob_
 

Detailed Description

Isotope pattern generator for coarse isotope distributions.

This algorithm implements IsotopePatternGenerator and generates theoretical pattern distributions for empirical formulas with resolution of 1Da. It convolves the empirical abundances of each element in a molecular formula, thus producing accurate intensities (probabilities) for each isotopic peak. However, it will assume that every isotope has an atomic mass that is rounded to the closest integer in Daltons, therefore it produces coarse distributions (it does not discriminate between 13C, N15 and O18 peaks). For example, for a molecule that contains both Carbon and Nitrogen, it will add up the probabilities for 13C and 15N, ignoring the fact that their masses are (slightly) different. Therefore, the probability distributions generated by the CoarseIsotopePatternGenerator are accurate, but the masses (m/z) are only approximately accurate. In case you need fine resolution, please consider using the FineIsotopePatternGenerator.

The output is a list of pairs containing nominal isotope probabilities paired with a number that is either an approximately accurate or rounded (integer) mass. The accurate masses assume the nominal isotopes are mostly due to (13)Carbon. To return accurate vs rounded masses, use setRoundMasses accordingly. The default is to return accurate masses (note that setting this option will not influence the probabilities and still produce a coarse distributions spaced at ca 1Da). For example, using rounded mass, for a C100 molecule, you will get:

1200 : 0.341036528
1201 : 0.368855864
1202 : 0.197477505
1203 : 0.0697715357

while accurate mass will produce:

1200 : 0.341036528
1201.00335 : 0.368855864
1202.00671 : 0.197477505
1203.01006 : 0.0697715357

The other important value which needs to be set is the max isotope value. This value can be set using the setMaxIsotope method. It is an upper bound for the number of isotopes which are calculated If e.g., set to 3, only the first three isotopes, Monoisotopic mass, +1 and +2 are calculated.

Note
By default all possible isotopes are calculated, which leads to a large number of values, if the mass value is large!
If you need fine isotope distributions, consider using the FineIsotopePatternGenerator.

See also method run()

Constructor & Destructor Documentation

◆ CoarseIsotopePatternGenerator()

CoarseIsotopePatternGenerator ( const Size  max_isotope = 0,
const bool  round_masses = false 
)

◆ ~CoarseIsotopePatternGenerator()

Member Function Documentation

◆ approximateFromPeptideWeight()

static IsotopeDistribution approximateFromPeptideWeight ( double  mass,
UInt  num_peaks = 20,
UInt  charge = 1 
)
static

roughly approximate peptide IsotopeDistribution from monoisotopic weight using Poisson distribution. m/z values approximated by adding one neutron mass (divided by charge) for every peak, starting at the given monoisotopic weight. Foundation from: Bellew et al, https://dx.doi.org/10.1093/bioinformatics/btl276

This method is around 50 times faster than estimateFromPeptideWeight, but only an approximation. The following are the intensities of the first 6 peaks generated for a monoisotopic mass of 1000:

estimateFromPeptideWeight: 0.571133000;0.306181000;0.095811100;0.022036900;0.004092170;0.000644568 approximateFromPeptideWeight: 0.573753000;0.318752000;0.088542200;0.016396700;0.002277320;0.000253036

KL divergences of the first 20 intensities of estimateFromPeptideWeight and this approximation range from 4.97E-5 for a monoisotopic mass of 20 to 0.0144 for a mass of 2500. For comparison, when comparing an observed pattern with a theoretical ground truth, the observed pattern is said to be an isotopic pattern if the KL between the two is below 0.05 for 2 peaks and below 0.6 for >=6 peaks by Guo Ci Teo et al.

Parameters
average_weightm/z of monoisotopic peak (with charge = 1) to approximate the distribution of intensities for
num_peaksHow many peaks should be generated (independent of this->max_isotope)

◆ approximateIntensities()

static std::vector<double> approximateIntensities ( double  mass,
UInt  num_peaks = 20 
)
static

roughly approximate intensity distribution of peptidic isotope patterns from monoisotopic weight using Poisson distribution. Foundation from: Bellew et al, https://dx.doi.org/10.1093/bioinformatics/btl276

This method is around 100 times faster than estimateFromPeptideWeight, but only an approximation of the intensities. It does not return IsotopeDistribution but a vector of intensities. For an assessment of accuracy, see approximateFromPeptideWeight.

Parameters
average_weightm/z of monoisotopic peak (with charge = 1) to approximate the distribution of intensities for
num_peaksHow many peaks should be generated (independent of this->max_isotope)

◆ calcFragmentIsotopeDist()

IsotopeDistribution calcFragmentIsotopeDist ( const IsotopeDistribution fragment_isotope_dist,
const IsotopeDistribution comp_fragment_isotope_dist,
const std::set< UInt > &  precursor_isotopes,
const double  fragment_mono_mass 
) const

Calculate isotopic distribution for a fragment molecule.

This calculates the isotopic distribution for a fragment molecule given the isotopic distribution of the fragment and complementary fragment (as if they were precursors), and which precursor isotopes were isolated.

Note
Do consider normalising the distribution afterwards to get conditional probabilities.

Equations come from Rockwood, AL; Kushnir, MA; Nelson, GJ. in "Dissociation of Individual Isotopic Peaks: Predicting Isotopic Distributions of Product Ions in MSn"

Parameters
fragment_isotope_distthe isotopic distribution of the fragment (as if it was a precursor).
comp_fragment_isotope_distthe isotopic distribution of the complementary fragment (as if it was a precursor).
precursor_isotopesa list of which precursor isotopes were isolated. 0 corresponds to the mono-isotopic molecule (M0), 1->M1, etc.
fragment_mono_massthe monoisotopic mass of the fragment.
Precondition
fragment_isotope_dist and comp_fragment_isotope_dist are gapless (no missing isotopes between the min/max isotopes of the dist)

◆ calcFragmentIsotopeDist_()

IsotopeDistribution calcFragmentIsotopeDist_ ( const IsotopeDistribution::ContainerType fragment_isotope_dist,
const IsotopeDistribution::ContainerType comp_fragment_isotope_dist,
const std::set< UInt > &  precursor_isotopes 
) const
protected

calculates the fragment distribution for a fragment molecule and stores it in result.

Parameters
fragment_isotope_distthe isotopic distribution of the fragment (as if it was a precursor).
comp_fragment_isotope_distthe isotopic distribution of the complementary fragment (as if it was a precursor).
precursor_isotopeswhich precursor isotopes were isolated. 0 corresponds to the mono-isotopic molecule (M0), 1->M1, etc.

◆ convolve()

convolves the distributions left and right and stores the result in result

◆ convolvePow_()

IsotopeDistribution::ContainerType convolvePow_ ( const IsotopeDistribution::ContainerType input,
Size  factor 
) const
protected

convolves the distribution input factor times and stores the result in result

◆ convolveSquare_()

IsotopeDistribution::ContainerType convolveSquare_ ( const IsotopeDistribution::ContainerType input) const
protected

convolves the distribution input with itself and stores the result in result

◆ correctMass_()

IsotopeDistribution::ContainerType correctMass_ ( const IsotopeDistribution::ContainerType input,
const double  mono_weight 
) const
protected

converts the masses of distribution input from atomic numbers to accurate masses

◆ estimateForFragmentFromDNAWeight()

IsotopeDistribution estimateForFragmentFromDNAWeight ( double  average_weight_precursor,
double  average_weight_fragment,
const std::set< UInt > &  precursor_isotopes 
)

Estimate DNA fragment IsotopeDistribution from the precursor's average weight, fragment's average weight, and a list of isolated precursor isotopes.

The max_depth of the isotopic distribution is set to max(precursor_isotopes)+1.

Parameters
average_weight_precursoraverage weight of the precursor nucleotide
average_weight_fragmentaverage weight of the fragment
precursor_isotopesthe precursor isotopes that were isolated. 0 corresponds to the mono-isotopic molecule (M0), 1->M1, etc.
Precondition
average_weight_precursor >= average_weight_fragment
average_weight_precursor > 0
average_weight_fragment > 0
precursor_isotopes.size() > 0

◆ estimateForFragmentFromPeptideWeight()

IsotopeDistribution estimateForFragmentFromPeptideWeight ( double  average_weight_precursor,
double  average_weight_fragment,
const std::set< UInt > &  precursor_isotopes 
)

Estimate peptide fragment IsotopeDistribution from the precursor's average weight, fragment's average weight, and a list of isolated precursor isotopes.

The max_depth of the isotopic distribution is set to max(precursor_isotopes)+1.

Parameters
average_weight_precursoraverage weight of the precursor peptide
average_weight_fragmentaverage weight of the fragment
precursor_isotopesthe precursor isotopes that were isolated. 0 corresponds to the mono-isotopic molecule (M0), 1->M1, etc.
Precondition
average_weight_precursor >= average_weight_fragment
average_weight_fragment > 0
average_weight_precursor > 0
precursor_isotopes.size() > 0

◆ estimateForFragmentFromPeptideWeightAndS()

IsotopeDistribution estimateForFragmentFromPeptideWeightAndS ( double  average_weight_precursor,
UInt  S_precursor,
double  average_weight_fragment,
UInt  S_fragment,
const std::set< UInt > &  precursor_isotopes 
) const

Estimate peptide fragment IsotopeDistribution from the precursor's average weight, number of sulfurs in the precursor, fragment's average weight, number of sulfurs in the fragment, and a list of isolated precursor isotopes.

The max_depth of the isotopic distribution is set to max(precursor_isotopes)+1.

Parameters
average_weight_precursoraverage weight of the precursor peptide
S_precursorThe exact number of Sulfurs in the precursor peptide
average_weight_fragmentaverage weight of the fragment
S_fragmentThe exact number of Sulfurs in the fragment
precursor_isotopesthe precursor isotopes that were isolated
Precondition
S_fragment <= average_weight_fragment / average_weight(sulfur)
S_precursor - S_fragment <= (average_weight_precursor - average_weight_fragment) / average_weight(sulfur)
average_weight_precursor >= average_weight_fragment
average_weight_precursor > 0
average_weight_fragment > 0
precursor_isotopes.size() > 0

◆ estimateForFragmentFromRNAWeight()

IsotopeDistribution estimateForFragmentFromRNAWeight ( double  average_weight_precursor,
double  average_weight_fragment,
const std::set< UInt > &  precursor_isotopes 
)

Estimate RNA fragment IsotopeDistribution from the precursor's average weight, fragment's average weight, and a list of isolated precursor isotopes.

The max_depth of the isotopic distribution is set to max(precursor_isotopes)+1.

Parameters
average_weight_precursoraverage weight of the precursor nucleotide
average_weight_fragmentaverage weight of the fragment
precursor_isotopesthe precursor isotopes that were isolated. 0 corresponds to the mono-isotopic molecule (M0), 1->M1, etc.
Precondition
average_weight_precursor >= average_weight_fragment
average_weight_precursor > 0
average_weight_fragment > 0
precursor_isotopes.size() > 0

◆ estimateForFragmentFromWeightAndComp()

IsotopeDistribution estimateForFragmentFromWeightAndComp ( double  average_weight_precursor,
double  average_weight_fragment,
const std::set< UInt > &  precursor_isotopes,
double  C,
double  H,
double  N,
double  O,
double  S,
double  P 
) const

Estimate fragment IsotopeDistribution from the precursor's average weight, fragment's average weight, a list of isolated precursor isotopes, and average composition.

The max_depth of the isotopic distribution is set to max(precursor_isotopes)+1.

Parameters
average_weight_precursoraverage weight of the precursor molecule
average_weight_fragmentaverage weight of the fragment molecule
precursor_isotopesthe precursor isotopes that were isolated. 0 corresponds to the mono-isotopic molecule (M0), 1->M1, etc.
CThe approximate relative stoichiometry of Carbons to other elements in this molecule
HThe approximate relative stoichiometry of Hydrogens to other elements in this molecule
NThe approximate relative stoichiometry of Nitrogens to other elements in this molecule
OThe approximate relative stoichiometry of Oxygens to other elements in this molecule
SThe approximate relative stoichiometry of Sulfurs to other elements in this molecule
PThe approximate relative stoichiometry of Phosphoruses to other elements in this molecule
Precondition
S, C, H, N, O, P >= 0
average_weight_precursor >= average_weight_fragment
average_weight_precursor > 0
average_weight_fragment > 0
precursor_isotopes.size() > 0

◆ estimateFromDNAWeight()

IsotopeDistribution estimateFromDNAWeight ( double  average_weight)

Estimate Nucleotide Isotopedistribution from weight and number of isotopes that should be reported averagine model from Zubarev, R. A.; Demirev, P. A. in "Isotope depletion of large biomolecules: Implications for molecular mass measurements.".

◆ estimateFromMonoWeightAndComp()

IsotopeDistribution estimateFromMonoWeightAndComp ( double  mono_weight,
double  C,
double  H,
double  N,
double  O,
double  S,
double  P 
)

Estimate Isotopedistribution from monoisotopic weight, average composition, and number of isotopes that should be reported.

◆ estimateFromPeptideMonoWeight()

IsotopeDistribution estimateFromPeptideMonoWeight ( double  mono_weight)

Estimate Peptide Isotopedistribution from monoisotopic weight and number of isotopes that should be reported.

Implementation using the averagine model proposed by Senko et al. in "Determination of Monoisotopic Masses and Ion Populations for Large Biomolecules from Resolved Isotopic Distributions" But this function takes monoisotopic mass. Thus determination of monoisotopic mass is not performed.

◆ estimateFromPeptideWeight()

IsotopeDistribution estimateFromPeptideWeight ( double  average_weight)

Estimate Peptide Isotopedistribution from weight and number of isotopes that should be reported.

Implementation using the averagine model proposed by Senko et al. in "Determination of Monoisotopic Masses and Ion Populations for Large Biomolecules from Resolved Isotopic Distributions"

Referenced by IsotopeMarker::apply().

◆ estimateFromPeptideWeightAndS()

IsotopeDistribution estimateFromPeptideWeightAndS ( double  average_weight,
UInt  S 
)

Estimate peptide IsotopeDistribution from average weight and exact number of sulfurs.

Parameters
average_weightAverage weight to estimate an EmpiricalFormula for
SThe exact number of Sulfurs in this molecule
Precondition
S <= average_weight / average_weight(sulfur)
average_weight >= 0

◆ estimateFromRNAMonoWeight()

IsotopeDistribution estimateFromRNAMonoWeight ( double  mono_weight)

Estimate Nucleotide Isotopedistribution from monoisotopic weight and number of isotopes that should be reported.

averagine model from Zubarev, R. A.; Demirev, P. A. in "Isotope depletion of large biomolecules: Implications for molecular mass measurements."

◆ estimateFromRNAWeight()

IsotopeDistribution estimateFromRNAWeight ( double  average_weight)

Estimate Nucleotide Isotopedistribution from weight and number of isotopes that should be reported.

averagine model from Zubarev, R. A.; Demirev, P. A. in "Isotope depletion of large biomolecules: Implications for molecular mass measurements."

◆ estimateFromWeightAndComp()

IsotopeDistribution estimateFromWeightAndComp ( double  average_weight,
double  C,
double  H,
double  N,
double  O,
double  S,
double  P 
)

Estimate Isotopedistribution from weight, average composition, and number of isotopes that should be reported.

◆ estimateFromWeightAndCompAndS()

IsotopeDistribution estimateFromWeightAndCompAndS ( double  average_weight,
UInt  S,
double  C,
double  H,
double  N,
double  O,
double  P 
)

Estimate IsotopeDistribution from weight, exact number of sulfurs, and average remaining composition.

Parameters
average_weightAverage weight to estimate an IsotopeDistribution for
SThe exact numbers of Sulfurs in this molecule
CThe approximate relative stoichiometry of Carbons to other elements (excluding Sulfur) in this molecule
HThe approximate relative stoichiometry of Hydrogens to other elements (excluding Sulfur) in this molecule
NThe approximate relative stoichiometry of Nitrogens to other elements (excluding Sulfur) in this molecule
OThe approximate relative stoichiometry of Oxygens to other elements (excluding Sulfur) in this molecule
PThe approximate relative stoichiometry of Phosphoruses to other elements (excluding Sulfur) in this molecule
Precondition
S, C, H, N, O, P >= 0
average_weight >= 0

◆ fillGaps_()

fill a gapped isotope pattern (i.e. certain masses are missing), with zero probability masses

◆ getMaxIsotope()

Size getMaxIsotope ( ) const

returns the currently set maximum isotope

◆ getRoundMasses()

bool getRoundMasses ( ) const

returns the current value of the flag to return expected masses (true) or atomic numbers (false).

◆ operator=()

◆ run()

IsotopeDistribution run ( const EmpiricalFormula ) const
overridevirtual

Creates an isotope distribution from an empirical sum formula.

Iterates through all elements, convolves them according to the number of atoms from that element and sums up the result.

Implements IsotopePatternGenerator.

◆ setMaxIsotope()

void setMaxIsotope ( const Size max_isotope)

sets the maximal isotope with max_isotope

sets the maximal isotope which is included in the distribution and used to limit the calculations. This is useful as distributions with numerous isotopes tend to have a lot of numerical zeros at the end

◆ setRoundMasses()

void setRoundMasses ( const bool  round_masses)

sets the round_masses_ flag to round masses to integer values (true) or return accurate masses (false)

Member Data Documentation

◆ max_isotope_

Size max_isotope_
protected

maximal isotopes which is used to calculate the distribution

◆ round_masses_

bool round_masses_
protected

flag to determine whether masses should be rounded or not