OpenMS
IsotopeWaveletTransform< PeakType > Class Template Reference

A class implementing the isotope wavelet transform. If you just want to find features using the isotope wavelet, take a look at the FeatureFinderAlgorithmIsotopeWavelet class. Usually, you only have to consider the class at hand if you plan to change the basic implementation of the transform. More...

#include <OpenMS/TRANSFORMATIONS/FEATUREFINDER/IsotopeWaveletTransform.h>

Collaboration diagram for IsotopeWaveletTransform< PeakType >:
[legend]

Classes

struct  BoxElement
 Internally used data structure. More...
 
class  TransSpectrum
 Internally (only by GPUs) used data structure . It allows efficient data exchange between CPU and GPU and avoids unnecessary memory moves. The class is tailored on the isotope wavelet transform and is in general not applicable on similar - but different - situations. More...
 

Public Types

typedef std::multimap< UInt, BoxElementBox
 Key: RT index, value: BoxElement. More...
 

Public Member Functions

 IsotopeWaveletTransform (const double min_mz, const double max_mz, const UInt max_charge, const Size max_scan_size=0, const bool hr_data=false, const String intenstype="ref")
 Constructor. More...
 
virtual ~IsotopeWaveletTransform ()
 Destructor. More...
 
virtual void getTransform (MSSpectrum &c_trans, const MSSpectrum &c_ref, const UInt c)
 Computes the isotope wavelet transform of charge state c. More...
 
virtual void getTransformHighRes (MSSpectrum &c_trans, const MSSpectrum &c_ref, const UInt c)
 Computes the isotope wavelet transform of charge state c. More...
 
virtual void identifyCharge (const MSSpectrum &candidates, const MSSpectrum &ref, const UInt scan_index, const UInt c, const double ampl_cutoff, const bool check_PPMs)
 Given an isotope wavelet transformed spectrum candidates, this function assigns to every significant pattern its corresponding charge state and a score indicating the reliability of the prediction. The result of this process is stored internally. Important: Before calling this function, apply updateRanges() to the original map. More...
 
virtual void initializeScan (const MSSpectrum &c_ref, const UInt c=0)
 
void updateBoxStates (const PeakMap &map, const Size scan_index, const UInt RT_interleave, const UInt RT_votes_cutoff, const Int front_bound=-1, const Int end_bound=-1)
 A function keeping track of currently open and closed sweep line boxes. This function is used by the isotope wavelet feature finder and must be called for each processed scan. More...
 
FeatureMap mapSeeds2Features (const PeakMap &map, const UInt RT_votes_cutoff)
 Filters the candidates further more and maps the internally used data structures to the OpenMS framework. More...
 
virtual std::multimap< double, BoxgetClosedBoxes ()
 Returns the closed boxes. More...
 
double getLinearInterpolation (const typename MSSpectrum::const_iterator &left_iter, const double mz_pos, const typename MSSpectrum::const_iterator &right_iter)
 Computes a linear (intensity) interpolation. More...
 
double getLinearInterpolation (const double mz_a, const double intens_a, const double mz_pos, const double mz_b, const double intens_b)
 Computes a linear (intensity) interpolation. More...
 
double getSigma () const
 
void setSigma (const double sigma)
 
virtual void computeMinSpacing (const MSSpectrum &c_ref)
 
double getMinSpacing () const
 
Size getMaxScanSize () const
 

Protected Member Functions

 IsotopeWaveletTransform ()
 Default Constructor. More...
 
void sampleTheCMarrWavelet_ (const MSSpectrum &scan, const Int wavelet_length, const Int mz_index, const UInt charge)
 
virtual double scoreThis_ (const TransSpectrum &candidate, const UInt peak_cutoff, const double seed_mz, const UInt c, const double ampl_cutoff)
 Given a candidate for an isotopic pattern, this function computes the corresponding score. More...
 
virtual double scoreThis_ (const MSSpectrum &candidate, const UInt peak_cutoff, const double seed_mz, const UInt c, const double ampl_cutoff)
 Given a candidate for an isotopic pattern, this function computes the corresponding score. More...
 
virtual bool checkPositionForPlausibility_ (const TransSpectrum &candidate, const MSSpectrum &ref, const double seed_mz, const UInt c, const UInt scan_index, const bool check_PPMs, const double transintens, const double prev_score)
 A ugly but necessary function to handle "off-by-1-Dalton predictions" due to idiosyncrasies of the data set (in comparison to the averagine model) More...
 
virtual bool checkPositionForPlausibility_ (const MSSpectrum &candidate, const MSSpectrum &ref, const double seed_mz, const UInt c, const UInt scan_index, const bool check_PPMs, const double transintens, const double prev_score)
 A ugly but necessary function to handle "off-by-1-Dalton predictions" due to idiosyncrasies of the data set (in comparison to the averagine model) More...
 
virtual std::pair< double, doublecheckPPMTheoModel_ (const MSSpectrum &ref, const double c_mz, const UInt c)
 
double getAvIntens_ (const TransSpectrum &scan)
 Computes the average (transformed) intensity (neglecting negative values) of scan. More...
 
double getAvIntens_ (const MSSpectrum &scan)
 Computes the average intensity (neglecting negative values) of scan. More...
 
double getSdIntens_ (const TransSpectrum &scan, const double mean)
 Computes the standard deviation (neglecting negative values) of the (transformed) intensities of scan. More...
 
double getSdIntens_ (const MSSpectrum &scan, const double mean)
 Computes the standard deviation (neglecting negative values) of the intensities of scan. More...
 
virtual void push2Box_ (const double mz, const UInt scan, UInt c, const double score, const double intens, const double rt, const UInt MZ_begin, const UInt MZ_end, const double ref_intens)
 Inserts a potential isotopic pattern into an open box or - if no such box exists - creates a new one. More...
 
virtual void push2TmpBox_ (const double mz, const UInt scan, UInt charge, const double score, const double intens, const double rt, const UInt MZ_begin, const UInt MZ_end)
 Essentially the same function as. More...
 
double getAvMZSpacing_ (const MSSpectrum &scan)
 Computes the average MZ spacing of scan. More...
 
void clusterSeeds_ (const TransSpectrum &candidates, const MSSpectrum &ref, const UInt scan_index, const UInt c, const bool check_PPMs)
 Clusters the seeds stored by push2TmpBox_. More...
 
virtual void clusterSeeds_ (const MSSpectrum &candidates, const MSSpectrum &ref, const UInt scan_index, const UInt c, const bool check_PPMs)
 Clusters the seeds stored by push2TmpBox_. More...
 
void extendBox_ (const PeakMap &map, const Box &box)
 A currently still necessary function that extends the box box in order to capture also signals whose isotopic pattern is nearly diminishing. More...
 
double peptideMassRule_ (const double c_mass) const
 Returns the monoisotopic mass (with corresponding decimal values) we would expect at c_mass. More...
 
double getPPMs_ (const double mass_a, const double mass_b) const
 Returns the parts-per-million deviation of the masses. More...
 

Protected Attributes

std::multimap< double, Boxopen_boxes_
 
std::multimap< double, Boxclosed_boxes_
 
std::multimap< double, Boxend_boxes_
 
std::multimap< double, Boxfront_boxes_
 
std::vector< std::multimap< double, Box > > * tmp_boxes_
 
double av_MZ_spacing_
 
double sigma_
 
std::vector< doublec_mzs_
 
std::vector< doublec_spacings_
 
std::vector< doublepsi_
 
std::vector< doubleprod_
 
std::vector< doublexs_
 
std::vector< doubleinterpol_xs_
 
std::vector< doubleinterpol_ys_
 
Size max_scan_size_
 
UInt max_num_peaks_per_pattern_
 
UInt max_charge_
 
UInt data_length_
 
bool hr_data_
 
String intenstype_
 
Int from_max_to_left_
 
Int from_max_to_right_
 
std::vector< int > indices_
 
double min_spacing_
 
double max_mz_cutoff_
 
std::vector< float > scores_
 
std::vector< float > zeros_
 

Detailed Description

template<typename PeakType>
class OpenMS::IsotopeWaveletTransform< PeakType >

A class implementing the isotope wavelet transform. If you just want to find features using the isotope wavelet, take a look at the FeatureFinderAlgorithmIsotopeWavelet class. Usually, you only have to consider the class at hand if you plan to change the basic implementation of the transform.


Class Documentation

◆ OpenMS::IsotopeWaveletTransform::BoxElement

struct OpenMS::IsotopeWaveletTransform::BoxElement
template<typename PeakType>
struct OpenMS::IsotopeWaveletTransform< PeakType >::BoxElement

Internally used data structure.

Collaboration diagram for IsotopeWaveletTransform< PeakType >::BoxElement:
[legend]
Class Members
UInt c Note, this is not the charge (it is charge-1!!!)
double intens The transformed intensity at the monoisotopic mass.
double mz The monoisotopic position.
UInt MZ_begin Index.
UInt MZ_end Index.
double ref_intens
double RT The elution time (not the scan index)
UInt RT_index The elution time (map) index.
double score The associated score.

Member Typedef Documentation

◆ Box

typedef std::multimap<UInt, BoxElement> Box

Key: RT index, value: BoxElement.

Constructor & Destructor Documentation

◆ IsotopeWaveletTransform() [1/2]

IsotopeWaveletTransform ( const double  min_mz,
const double  max_mz,
const UInt  max_charge,
const Size  max_scan_size = 0,
const bool  hr_data = false,
const String  intenstype = "ref" 
)

Constructor.

Parameters
min_mzThe smallest m/z value occurring in your map.
max_mzThe largest m/z value occurring in your map.
max_chargeThe highest charge state you would like to consider.

References OpenMS::Constants::DEFAULT_NUM_OF_INTERPOLATION_POINTS, IsotopeWavelet::getMzPeakCutOffAtMonoPos(), IsotopeWavelet::getNumPeakCutOff(), IsotopeWavelet::init(), and OpenMS::Constants::IW_NEUTRON_MASS.

◆ ~IsotopeWaveletTransform()

Destructor.

◆ IsotopeWaveletTransform() [2/2]

Default Constructor.

Note
Provided just for inheritance reasons. You should always use the other constructor.

Member Function Documentation

◆ checkPositionForPlausibility_() [1/2]

bool checkPositionForPlausibility_ ( const MSSpectrum candidate,
const MSSpectrum ref,
const double  seed_mz,
const UInt  c,
const UInt  scan_index,
const bool  check_PPMs,
const double  transintens,
const double  prev_score 
)
protectedvirtual

A ugly but necessary function to handle "off-by-1-Dalton predictions" due to idiosyncrasies of the data set (in comparison to the averagine model)

Parameters
candidateThe wavelet transformed spectrum containing the candidate.
refThe original spectrum containing the candidate.
seed_mzThe m/z position of the candidate pattern.
cThe predicted charge state minus 1 (e.g. c=2 means charge state 3) of the candidate.
scan_indexThe index of the scan under consideration (w.r.t. the original map).

References OpenMS::Constants::c, IsotopeWavelet::getMzPeakCutOffAtMonoPos(), IsotopeWavelet::getNumPeakCutOff(), MSSpectrum::getRT(), OpenMS::Constants::IW_QUARTER_NEUTRON_MASS, and MSSpectrum::MZBegin().

◆ checkPositionForPlausibility_() [2/2]

bool checkPositionForPlausibility_ ( const TransSpectrum candidate,
const MSSpectrum ref,
const double  seed_mz,
const UInt  c,
const UInt  scan_index,
const bool  check_PPMs,
const double  transintens,
const double  prev_score 
)
protectedvirtual

A ugly but necessary function to handle "off-by-1-Dalton predictions" due to idiosyncrasies of the data set (in comparison to the averagine model)

Parameters
candidateThe wavelet transformed spectrum containing the candidate.
refThe original spectrum containing the candidate.
seed_mzThe m/z position of the candidate pattern.
cThe predicted charge state minus 1 (e.g. c=2 means charge state 3) of the candidate.
scan_indexThe index of the scan under consideration (w.r.t. the original map).

References IsotopeWaveletTransform< PeakType >::TransSpectrum::begin(), OpenMS::Constants::c, IsotopeWaveletTransform< PeakType >::TransSpectrum::end(), IsotopeWavelet::getMzPeakCutOffAtMonoPos(), IsotopeWavelet::getNumPeakCutOff(), MSSpectrum::getRT(), OpenMS::Constants::IW_QUARTER_NEUTRON_MASS, IsotopeWaveletTransform< PeakType >::TransSpectrum::MZBegin(), and MSSpectrum::MZBegin().

◆ checkPPMTheoModel_()

std::pair< double, double > checkPPMTheoModel_ ( const MSSpectrum ref,
const double  c_mz,
const UInt  c 
)
protectedvirtual

◆ clusterSeeds_() [1/2]

void clusterSeeds_ ( const MSSpectrum candidates,
const MSSpectrum ref,
const UInt  scan_index,
const UInt  c,
const bool  check_PPMs 
)
protectedvirtual

Clusters the seeds stored by push2TmpBox_.

Parameters
candidatesA isotope wavelet transformed spectrum.
refThe corresponding original spectrum (w.r.t. candidates).
scan_indexThe index of the scan under consideration (w.r.t. the original map).

References OpenMS::Constants::c, IsotopeWaveletTransform< PeakType >::BoxElement::c, IsotopeWaveletTransform< PeakType >::BoxElement::intens, IsotopeWaveletTransform< PeakType >::BoxElement::mz, IsotopeWaveletTransform< PeakType >::BoxElement::RT, and IsotopeWaveletTransform< PeakType >::BoxElement::score.

◆ clusterSeeds_() [2/2]

void clusterSeeds_ ( const TransSpectrum candidates,
const MSSpectrum ref,
const UInt  scan_index,
const UInt  c,
const bool  check_PPMs 
)
protected

Clusters the seeds stored by push2TmpBox_.

Parameters
candidatesA isotope wavelet transformed spectrum.
refThe corresponding original spectrum (w.r.t. candidates).
scan_indexThe index of the scan under consideration (w.r.t. the original map).

References OpenMS::Constants::c, IsotopeWaveletTransform< PeakType >::BoxElement::c, IsotopeWaveletTransform< PeakType >::BoxElement::intens, IsotopeWaveletTransform< PeakType >::BoxElement::mz, IsotopeWaveletTransform< PeakType >::BoxElement::RT, and IsotopeWaveletTransform< PeakType >::BoxElement::score.

◆ computeMinSpacing()

void computeMinSpacing ( const MSSpectrum c_ref)
virtual

◆ extendBox_()

void extendBox_ ( const PeakMap map,
const Box box 
)
protected

A currently still necessary function that extends the box box in order to capture also signals whose isotopic pattern is nearly diminishing.

Parameters
mapThe experimental map.
boxThe box to be extended.

References MSExperiment::begin().

◆ getAvIntens_() [1/2]

double getAvIntens_ ( const MSSpectrum scan)
inlineprotected

Computes the average intensity (neglecting negative values) of scan.

◆ getAvIntens_() [2/2]

double getAvIntens_ ( const TransSpectrum scan)
inlineprotected

Computes the average (transformed) intensity (neglecting negative values) of scan.

References IsotopeWaveletTransform< PeakType >::TransSpectrum::getTransIntensity(), and IsotopeWaveletTransform< PeakType >::TransSpectrum::size().

◆ getAvMZSpacing_()

double getAvMZSpacing_ ( const MSSpectrum scan)
inlineprotected

Computes the average MZ spacing of scan.

Parameters
scanThe scan we are interested in.

◆ getClosedBoxes()

virtual std::multimap< double, Box > getClosedBoxes ( )
inlinevirtual

Returns the closed boxes.

References IsotopeWaveletTransform< PeakType >::closed_boxes_.

◆ getLinearInterpolation() [1/2]

double getLinearInterpolation ( const double  mz_a,
const double  intens_a,
const double  mz_pos,
const double  mz_b,
const double  intens_b 
)
inline

Computes a linear (intensity) interpolation.

Parameters
mz_aThe m/z value of the point left to the query.
intens_aThe intensity value of the point left to the query.
mz_posThe query point.
mz_bThe m/z value of the point right to the query.
intens_bThe intensity value of the point left to the query.

◆ getLinearInterpolation() [2/2]

double getLinearInterpolation ( const typename MSSpectrum::const_iterator &  left_iter,
const double  mz_pos,
const typename MSSpectrum::const_iterator &  right_iter 
)
inline

Computes a linear (intensity) interpolation.

Parameters
left_iterThe point left to the query.
mz_posThe query point.
right_iterThe point right to the query.

◆ getMaxScanSize()

Size getMaxScanSize ( ) const
inline

◆ getMinSpacing()

double getMinSpacing ( ) const
inline

◆ getPPMs_()

double getPPMs_ ( const double  mass_a,
const double  mass_b 
) const
inlineprotected

Returns the parts-per-million deviation of the masses.

Parameters
mass_aThe first mass.
mass_bThe second mass.

◆ getSdIntens_() [1/2]

double getSdIntens_ ( const MSSpectrum scan,
const double  mean 
)
inlineprotected

Computes the standard deviation (neglecting negative values) of the intensities of scan.

References OpenMS::Math::mean().

◆ getSdIntens_() [2/2]

double getSdIntens_ ( const TransSpectrum scan,
const double  mean 
)
inlineprotected

Computes the standard deviation (neglecting negative values) of the (transformed) intensities of scan.

References IsotopeWaveletTransform< PeakType >::TransSpectrum::getTransIntensity(), OpenMS::Math::mean(), and IsotopeWaveletTransform< PeakType >::TransSpectrum::size().

◆ getSigma()

double getSigma ( ) const
inline

◆ getTransform()

void getTransform ( MSSpectrum c_trans,
const MSSpectrum c_ref,
const UInt  c 
)
virtual

Computes the isotope wavelet transform of charge state c.

Parameters
c_transThe transform.
c_refThe reference spectrum.
cThe charge state minus 1 (e.g. c=2 means charge state 3) at which you want to compute the transform.

References OpenMS::Constants::c, IsotopeWavelet::getLambdaL(), IsotopeWavelet::getMzPeakCutOffAtMonoPos(), IsotopeWavelet::getValueByLambda(), and OpenMS::Constants::IW_QUARTER_NEUTRON_MASS.

◆ getTransformHighRes()

void getTransformHighRes ( MSSpectrum c_trans,
const MSSpectrum c_ref,
const UInt  c 
)
virtual

Computes the isotope wavelet transform of charge state c.

Parameters
c_transThe transform.
c_refThe reference spectrum.
cThe charge state minus 1 (e.g. c=2 means charge state 3) at which you want to compute the transform.

References OpenMS::Constants::c, IsotopeWavelet::getLambdaL(), IsotopeWavelet::getMzPeakCutOffAtMonoPos(), IsotopeWavelet::getValueByLambda(), and OpenMS::Constants::IW_QUARTER_NEUTRON_MASS.

◆ identifyCharge()

void identifyCharge ( const MSSpectrum candidates,
const MSSpectrum ref,
const UInt  scan_index,
const UInt  c,
const double  ampl_cutoff,
const bool  check_PPMs 
)
virtual

Given an isotope wavelet transformed spectrum candidates, this function assigns to every significant pattern its corresponding charge state and a score indicating the reliability of the prediction. The result of this process is stored internally. Important: Before calling this function, apply updateRanges() to the original map.

Parameters
candidatesA isotope wavelet transformed spectrum. Entry "number i" in this vector must correspond to the charge-"(i-1)"-transform of its mass signal. (This is exactly the output of the function
See also
getTransforms.)
Parameters
refThe reference scan (the untransformed raw data) corresponding to candidates.
cThe corresponding charge state minus 1 (e.g. c=2 means charge state 3)
scan_indexThe index of the scan (w.r.t. to some map) currently under consideration.
ampl_cutoffThe thresholding parameter. This parameter is the only (and hence a really important) parameter of the isotope wavelet transform. On the basis of ampl_cutoff the program tries to distinguish between noise and signal. Please note that it is not a "simple" hard thresholding parameter in the sense of drawing a virtual line in the spectrum, which is then used as a guillotine cut. Maybe you should play around a bit with this parameter to get a feeling about its range. For peptide mass fingerprints on small data sets (like single MALDI-scans e.g.), it makes sense to start ampl_cutoff=0 or even ampl_cutoff=-1, indicating no thresholding at all. Note that also ampl_cutoff=0 triggers (a moderate) thresholding based on the average intensity in the wavelet transform.
check_PPMsIf enabled, the algorithm will check each monoisotopic mass candidate for its plausibility by computing the ppm difference between this mass and the averagine model.

References ConstRefVector< ContainerT >::begin(), OpenMS::Constants::c, ConstRefVector< ContainerT >::end(), IsotopeWavelet::getMzPeakCutOffAtMonoPos(), IsotopeWavelet::getNumPeakCutOff(), MSSpectrum::getRT(), OpenMS::Constants::IW_NEUTRON_MASS, OpenMS::Constants::IW_QUARTER_NEUTRON_MASS, MSSpectrum::MZBegin(), MSSpectrum::MZEnd(), and ConstRefVector< ContainerT >::sortByIntensity().

◆ initializeScan()

◆ mapSeeds2Features()

FeatureMap mapSeeds2Features ( const PeakMap map,
const UInt  RT_votes_cutoff 
)

Filters the candidates further more and maps the internally used data structures to the OpenMS framework.

Parameters
mapThe original map containing the data set to be analyzed.
max_chargeThe maximal charge state under consideration.
RT_votes_cutoffSee the IsotopeWaveletFF class.

References ConvexHull2D::addPoints(), IsotopeWavelet::getLambdaL(), IsotopeWavelet::getMzPeakCutOffAtMonoPos(), OpenMS::Constants::IW_NEUTRON_MASS, OpenMS::Constants::IW_QUARTER_NEUTRON_MASS, MSSpectrum::MZBegin(), ExposedVector< VectorElement >::push_back(), BaseFeature::setCharge(), Feature::setConvexHulls(), Peak2D::setIntensity(), Peak2D::setMZ(), Feature::setOverallQuality(), Peak2D::setRT(), and MSExperiment::size().

◆ peptideMassRule_()

double peptideMassRule_ ( const double  c_mass) const
inlineprotected

Returns the monoisotopic mass (with corresponding decimal values) we would expect at c_mass.

Parameters
c_massThe mass for which we would like to know the averagine decimal places.

References OpenMS::Constants::PEPTIDE_MASS_RULE_BOUND, and OpenMS::Constants::PEPTIDE_MASS_RULE_FACTOR.

◆ push2Box_()

void push2Box_ ( const double  mz,
const UInt  scan,
UInt  c,
const double  score,
const double  intens,
const double  rt,
const UInt  MZ_begin,
const UInt  MZ_end,
const double  ref_intens 
)
protectedvirtual

Inserts a potential isotopic pattern into an open box or - if no such box exists - creates a new one.

Parameters
mzThe position of the pattern.
scanThe index of the scan, we are currently analyzing (w.r.t. the data map). This information is necessary for the post-processing (sweep lining).
chargeThe estimated charge state minus 1 (e.g. c=2 means charge state 3) of the pattern.
scoreThe pattern's score.
intensThe intensity at the monoisotopic peak.
rtThe retention time of the scan (similar to scan, but here: no index, but the real value).
MZ_beginThe starting index of the pattern (m/z) w.r.t. the current scan.
MZ_endThe end index (w.r.t. the monoisotopic position!) of the pattern (m/z) w.r.t. the current scan.

References OpenMS::Constants::c, IsotopeWaveletTransform< PeakType >::BoxElement::c, IsotopeWaveletTransform< PeakType >::BoxElement::intens, OpenMS::Constants::IW_HALF_NEUTRON_MASS, IsotopeWaveletTransform< PeakType >::BoxElement::mz, IsotopeWaveletTransform< PeakType >::BoxElement::MZ_begin, IsotopeWaveletTransform< PeakType >::BoxElement::MZ_end, IsotopeWaveletTransform< PeakType >::BoxElement::ref_intens, IsotopeWaveletTransform< PeakType >::BoxElement::RT, IsotopeWaveletTransform< PeakType >::BoxElement::RT_index, and IsotopeWaveletTransform< PeakType >::BoxElement::score.

◆ push2TmpBox_()

void push2TmpBox_ ( const double  mz,
const UInt  scan,
UInt  charge,
const double  score,
const double  intens,
const double  rt,
const UInt  MZ_begin,
const UInt  MZ_end 
)
protectedvirtual

Essentially the same function as.

See also
push2Box_. In contrast to
push2Box this function stores its candidates only temporarily. In particular, this function is only used within a single scan transform. After the wavelet transform is computed on that scan, all candidates are pushed by this function and finally clustered together by
clusterSeeds_. Afterwards, a final push by
push2Box_ is performed storing the clustered candidates.
Parameters
mzThe position of the pattern.
scanThe index of the scan, we are currently analyzing (w.r.t. the data map). This information is necessary for the post-processing (sweep lining).
chargeThe estimated charge state minus 1 (e.g. c=2 means charge state 3) of the pattern.
scoreThe pattern's score.
intensThe intensity at the monoisotopic peak.
rtThe retention time of the scan (similar to scan, but here: no index, but the real value).
MZ_beginThe starting index of the pattern (m/z) w.r.t. the current scan.
MZ_endThe end index (w.r.t. the monoisotopic position!) of the pattern (m/z) w.r.t. the current scan.

References OpenMS::Constants::c, IsotopeWaveletTransform< PeakType >::BoxElement::c, IsotopeWaveletTransform< PeakType >::BoxElement::intens, OpenMS::Constants::IW_HALF_NEUTRON_MASS, IsotopeWaveletTransform< PeakType >::BoxElement::mz, IsotopeWaveletTransform< PeakType >::BoxElement::MZ_begin, IsotopeWaveletTransform< PeakType >::BoxElement::MZ_end, IsotopeWaveletTransform< PeakType >::BoxElement::ref_intens, IsotopeWaveletTransform< PeakType >::BoxElement::RT, IsotopeWaveletTransform< PeakType >::BoxElement::RT_index, and IsotopeWaveletTransform< PeakType >::BoxElement::score.

◆ sampleTheCMarrWavelet_()

void sampleTheCMarrWavelet_ ( const MSSpectrum scan,
const Int  wavelet_length,
const Int  mz_index,
const UInt  charge 
)
inlineprotected

◆ scoreThis_() [1/2]

double scoreThis_ ( const MSSpectrum candidate,
const UInt  peak_cutoff,
const double  seed_mz,
const UInt  c,
const double  ampl_cutoff 
)
protectedvirtual

Given a candidate for an isotopic pattern, this function computes the corresponding score.

Parameters
candidateA isotope wavelet transformed spectrum.
peak_cutoffThe number of peaks we will consider for the isotopic pattern.
seed_mzThe predicted position of the monoisotopic peak.
cThe charge state minus 1 (e.g. c=2 means charge state 3) for which the score should be determined.
ampl_cutoffThe threshold.

References OpenMS::Constants::c, OpenMS::Constants::IW_HALF_NEUTRON_MASS, OpenMS::Constants::IW_NEUTRON_MASS, and MSSpectrum::MZBegin().

◆ scoreThis_() [2/2]

double scoreThis_ ( const TransSpectrum candidate,
const UInt  peak_cutoff,
const double  seed_mz,
const UInt  c,
const double  ampl_cutoff 
)
protectedvirtual

Given a candidate for an isotopic pattern, this function computes the corresponding score.

Parameters
candidateA isotope wavelet transformed spectrum.
peak_cutoffThe number of peaks we will consider for the isotopic pattern.
seed_mzThe predicted position of the monoisotopic peak.
cThe charge state minus 1 (e.g. c=2 means charge state 3) for which the score should be determined.
ampl_cutoffThe threshold.

References IsotopeWaveletTransform< PeakType >::TransSpectrum::begin(), OpenMS::Constants::c, IsotopeWaveletTransform< PeakType >::TransSpectrum::getMZ(), IsotopeWaveletTransform< PeakType >::TransSpectrum::getTransIntensity(), OpenMS::Constants::IW_HALF_NEUTRON_MASS, OpenMS::Constants::IW_NEUTRON_MASS, IsotopeWaveletTransform< PeakType >::TransSpectrum::MZBegin(), and IsotopeWaveletTransform< PeakType >::TransSpectrum::size().

◆ setSigma()

void setSigma ( const double  sigma)
inline

◆ updateBoxStates()

void updateBoxStates ( const PeakMap map,
const Size  scan_index,
const UInt  RT_interleave,
const UInt  RT_votes_cutoff,
const Int  front_bound = -1,
const Int  end_bound = -1 
)

A function keeping track of currently open and closed sweep line boxes. This function is used by the isotope wavelet feature finder and must be called for each processed scan.

Parameters
mapThe original map containing the data set to be analyzed.
scan_indexThe index of the scan currently under consideration w.r.t. its MS map. This information is necessary to sweep across the map after each scan has been evaluated.
RT_votes_cutoffSee the IsotopeWaveletFF class.

References MSExperiment::clear(), and MSExperiment::size().

Member Data Documentation

◆ av_MZ_spacing_

double av_MZ_spacing_
protected

◆ c_mzs_

std::vector<double> c_mzs_
protected

◆ c_spacings_

std::vector<double> c_spacings_
protected

◆ closed_boxes_

std::multimap<double, Box> closed_boxes_
protected

◆ data_length_

UInt data_length_
protected

◆ end_boxes_

std::multimap<double, Box> end_boxes_
protected

◆ from_max_to_left_

Int from_max_to_left_
protected

◆ from_max_to_right_

Int from_max_to_right_
protected

◆ front_boxes_

std::multimap<double, Box> front_boxes_
protected

◆ hr_data_

bool hr_data_
protected

◆ indices_

std::vector<int> indices_
protected

◆ intenstype_

String intenstype_
protected

◆ interpol_xs_

std::vector<double> interpol_xs_
protected

◆ interpol_ys_

std::vector<double> interpol_ys_
protected

◆ max_charge_

UInt max_charge_
protected

◆ max_mz_cutoff_

double max_mz_cutoff_
protected

◆ max_num_peaks_per_pattern_

UInt max_num_peaks_per_pattern_
protected

◆ max_scan_size_

Size max_scan_size_
protected

◆ min_spacing_

double min_spacing_
protected

◆ open_boxes_

std::multimap<double, Box> open_boxes_
protected

◆ prod_

std::vector<double> prod_
protected

◆ psi_

std::vector<double> psi_
protected

◆ scores_

std::vector<float> scores_
protected

◆ sigma_

◆ tmp_boxes_

std::vector<std::multimap<double, Box> >* tmp_boxes_
protected

◆ xs_

std::vector<double> xs_
protected

◆ zeros_

std::vector<float> zeros_
protected