OpenMS
TwoDOptimization Class Reference

This class provides the two-dimensional optimization of the picked peak parameters. More...

#include <OpenMS/TRANSFORMATIONS/RAW2PEAK/TwoDOptimization.h>

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

Classes

struct  Data
 Helper struct (contains the size of an area and a raw data container) More...
 
class  TwoDOptFunctor
 

Public Types

typedef MSExperiment::const_iterator InputSpectrumIterator
 

Public Member Functions

 TwoDOptimization ()
 Constructor. More...
 
 TwoDOptimization (const TwoDOptimization &opt)
 Copy constructor. More...
 
 ~TwoDOptimization () override
 Destructor. More...
 
TwoDOptimizationoperator= (const TwoDOptimization &opt)
 Assignment operator. More...
 
double getMZTolerance () const
 Non-mutable access to the matching epsilon. More...
 
void setMZTolerance (double tolerance_mz)
 Mutable access to the matching epsilon. More...
 
double getMaxPeakDistance () const
 Non-mutable access to the maximal peak distance in a cluster. More...
 
void setMaxPeakDistance (double max_peak_distance)
 Mutable access to the maximal peak distance in a cluster. More...
 
UInt getMaxIterations () const
 Non-mutable access to the maximal number of iterations. More...
 
void setMaxIterations (UInt max_iteration)
 Mutable access to the maximal number of iterations. More...
 
const OptimizationFunctions::PenaltyFactorsIntensitygetPenalties () const
 Non-mutable access to the minimal number of adjacent scans. More...
 
void setPenalties (OptimizationFunctions::PenaltyFactorsIntensity &penalties)
 Mutable access to the minimal number of adjacent scans. More...
 
void optimize (InputSpectrumIterator first, InputSpectrumIterator last, PeakMap &ms_exp, bool real2D=true)
 Find two dimensional peak clusters and optimize their peak parameters. More...
 
- Public Member Functions inherited from DefaultParamHandler
 DefaultParamHandler (const String &name)
 Constructor with name that is displayed in error messages. More...
 
 DefaultParamHandler (const DefaultParamHandler &rhs)
 Copy constructor. More...
 
virtual ~DefaultParamHandler ()
 Destructor. More...
 
DefaultParamHandleroperator= (const DefaultParamHandler &rhs)
 Assignment operator. More...
 
virtual bool operator== (const DefaultParamHandler &rhs) const
 Equality operator. More...
 
void setParameters (const Param &param)
 Sets the parameters. More...
 
const ParamgetParameters () const
 Non-mutable access to the parameters. More...
 
const ParamgetDefaults () const
 Non-mutable access to the default parameters. More...
 
const StringgetName () const
 Non-mutable access to the name. More...
 
void setName (const String &name)
 Mutable access to the name. More...
 
const std::vector< String > & getSubsections () const
 Non-mutable access to the registered subsections. More...
 

Protected Member Functions

Auxiliary Functions for the search of matching regions
std::vector< double >::iterator searchInScan_ (std::vector< double >::iterator scan_begin, std::vector< double >::iterator scan_end, double current_mz)
 
void optimizeRegions_ (InputSpectrumIterator &first, InputSpectrumIterator &last, PeakMap &ms_exp)
 
void optimizeRegionsScanwise_ (InputSpectrumIterator &first, InputSpectrumIterator &last, PeakMap &ms_exp)
 
void getRegionEndpoints_ (PeakMap &exp, InputSpectrumIterator &first, InputSpectrumIterator &last, Size iso_map_idx, double noise_level, TwoDOptimization::Data &d)
 Get the indices of the first and last raw data point of this region. More...
 
void findMatchingPeaks_ (std::multimap< double, IsotopeCluster >::iterator &it, PeakMap &ms_exp)
 Identify matching peak in a peak cluster. More...
 
void updateMembers_ () override
 update members method from DefaultParamHandler to update the members More...
 
- Protected Member Functions inherited from DefaultParamHandler
void defaultsToParam_ ()
 Updates the parameters after the defaults have been set in the constructor. More...
 

Protected Attributes

std::multimap< double, IsotopeClusteriso_map_
 stores the retention time of each isotopic cluster More...
 
std::multimap< double, IsotopeCluster >::const_iterator curr_region_
 Pointer to the current region. More...
 
double max_peak_distance_
 upper bound for distance between two peaks belonging to the same region More...
 
double tolerance_mz_
 threshold for the difference in the peak position of two matching peaks More...
 
std::map< Int, std::vector< PeakIndex > > matching_peaks_
 Indices of peaks in the adjacent scans matching peaks in the scan with no. ref_scan. More...
 
UInt max_iteration_
 Convergence Parameter: Maximal number of iterations. More...
 
bool real_2D_
 Optimization considering all scans of a cluster or optimization of each scan separately. More...
 
OptimizationFunctions::PenaltyFactorsIntensity penalties_
 Penalty factors for some parameters in the optimization. More...
 
- Protected Attributes inherited from DefaultParamHandler
Param param_
 Container for current parameters. More...
 
Param defaults_
 Container for default parameters. This member should be filled in the constructor of derived classes! More...
 
std::vector< Stringsubsections_
 Container for registered subsections. This member should be filled in the constructor of derived classes! More...
 
String error_name_
 Name that is displayed in error messages during the parameter checking. More...
 
bool check_defaults_
 If this member is set to false no checking if parameters in done;. More...
 
bool warn_empty_defaults_
 If this member is set to false no warning is emitted when defaults are empty;. More...
 

Additional Inherited Members

- Static Public Member Functions inherited from DefaultParamHandler
static void writeParametersToMetaValues (const Param &write_this, MetaInfoInterface &write_here, const String &key_prefix="")
 Writes all parameters to meta values. More...
 

Detailed Description

This class provides the two-dimensional optimization of the picked peak parameters.

Given the picked peaks, this class optimizes the peak parameters of each isotope pattern using a non-linear optimization. The peaks of adjacent scans are adjusted to achieve that a peak occurring in several scans has always the same m/z position. For the optimization the Levenberg-Marquardt algorithm provided from the Eigen is used. The optimized parameters are the m/z values, the left and right width, which shall be equal for a peak in all scans, and the peaks' heights.

Todo:
Works only with defined types due to pointers to the data in the optimization namespace! Change that or remove templates (Alexandra)
Parameters of this class are:

NameTypeDefaultRestrictionsDescription
iterations int10  maximal number of iterations for the fitting step
penalties:position float0.0  If the position changes more than 0.2Da during the fitting it can be penalized
penalties:height float1.0  penalty term for the fitting of the intensity:If it gets negative during the fitting it can be penalized.
penalties:left_width float0.0  penalty term for the fitting of the left width:If the left width gets too broad or negative during the fitting it can be penalized.
penalties:right_width float0.0  penalty term for the fitting of the right width:If the right width gets too broad or negative during the fitting it can be penalized.
2d:tolerance_mz float2.2  mz tolerance for cluster construction
2d:max_peak_distance float1.2  maximal peak distance in mz in a cluster

Note:
  • If a section name is documented, the documentation is displayed as tooltip.
  • Advanced parameter names are italic.

Class Documentation

◆ OpenMS::TwoDOptimization::Data

struct OpenMS::TwoDOptimization::Data

Helper struct (contains the size of an area and a raw data container)

Collaboration diagram for TwoDOptimization::Data:
[legend]
Class Members
multimap< double, IsotopeCluster >::iterator iso_map_iter
map< Int, vector< PeakIndex > > matching_peaks
PenaltyFactorsIntensity penalties
PeakMap picked_peaks
vector< double > positions
ConstIterator raw_data_first
vector< double > signal
vector< pair< SignedSize, SignedSize > > signal2D
Size total_nr_peaks

Member Typedef Documentation

◆ InputSpectrumIterator

Constructor & Destructor Documentation

◆ TwoDOptimization() [1/2]

Constructor.

◆ TwoDOptimization() [2/2]

Copy constructor.

◆ ~TwoDOptimization()

~TwoDOptimization ( )
inlineoverride

Destructor.

Member Function Documentation

◆ findMatchingPeaks_()

void findMatchingPeaks_ ( std::multimap< double, IsotopeCluster >::iterator &  it,
PeakMap ms_exp 
)
protected

Identify matching peak in a peak cluster.

◆ getMaxIterations()

UInt getMaxIterations ( ) const
inline

Non-mutable access to the maximal number of iterations.

◆ getMaxPeakDistance()

double getMaxPeakDistance ( ) const
inline

Non-mutable access to the maximal peak distance in a cluster.

◆ getMZTolerance()

double getMZTolerance ( ) const
inline

Non-mutable access to the matching epsilon.

◆ getPenalties()

const OptimizationFunctions::PenaltyFactorsIntensity& getPenalties ( ) const
inline

Non-mutable access to the minimal number of adjacent scans.

◆ getRegionEndpoints_()

void getRegionEndpoints_ ( PeakMap exp,
InputSpectrumIterator first,
InputSpectrumIterator last,
Size  iso_map_idx,
double  noise_level,
TwoDOptimization::Data d 
)
protected

Get the indices of the first and last raw data point of this region.

◆ operator=()

TwoDOptimization& operator= ( const TwoDOptimization opt)

Assignment operator.

◆ optimize()

void optimize ( InputSpectrumIterator  first,
InputSpectrumIterator  last,
PeakMap ms_exp,
bool  real2D = true 
)

Find two dimensional peak clusters and optimize their peak parameters.

Note
For the peak spectra, the following meta data arrays (see MSSpectrum) have to be present and have to be named just as listed here:
  • intensity (index:1)
  • leftWidth (index:3)
  • rightWidth (index:4)
  • peakShape (index:5)
Parameters
firstbegin of the raw data spectra iterator range
lastend of the raw data spectra iterator range
ms_exppeak map corresponding to the raw data in the range from first to last
real2Dflag if the optimization should be two dimensional or on each scan separately
Exceptions
Exception::IllegalArgumentis thrown if required meta information from peak picking is missing (area, shape, left width, right width) or if the input data is invalid in some other way

◆ optimizeRegions_()

void optimizeRegions_ ( InputSpectrumIterator first,
InputSpectrumIterator last,
PeakMap ms_exp 
)
protected

Performs 2D optimization of all regions

◆ optimizeRegionsScanwise_()

void optimizeRegionsScanwise_ ( InputSpectrumIterator first,
InputSpectrumIterator last,
PeakMap ms_exp 
)
protected

Performs an optimization of all regions by calling OptimizePick

◆ searchInScan_()

std::vector<double>::iterator searchInScan_ ( std::vector< double >::iterator  scan_begin,
std::vector< double >::iterator  scan_end,
double  current_mz 
)
protected

◆ setMaxIterations()

void setMaxIterations ( UInt  max_iteration)
inline

Mutable access to the maximal number of iterations.

◆ setMaxPeakDistance()

void setMaxPeakDistance ( double  max_peak_distance)
inline

Mutable access to the maximal peak distance in a cluster.

◆ setMZTolerance()

void setMZTolerance ( double  tolerance_mz)
inline

Mutable access to the matching epsilon.

◆ setPenalties()

void setPenalties ( OptimizationFunctions::PenaltyFactorsIntensity penalties)
inline

Mutable access to the minimal number of adjacent scans.

References PenaltyFactorsIntensity::height, PenaltyFactors::lWidth, PenaltyFactors::pos, and PenaltyFactors::rWidth.

◆ updateMembers_()

void updateMembers_ ( )
overrideprotectedvirtual

update members method from DefaultParamHandler to update the members

Reimplemented from DefaultParamHandler.

Member Data Documentation

◆ curr_region_

std::multimap<double, IsotopeCluster>::const_iterator curr_region_
protected

Pointer to the current region.

◆ iso_map_

std::multimap<double, IsotopeCluster> iso_map_
protected

stores the retention time of each isotopic cluster

◆ matching_peaks_

std::map<Int, std::vector<PeakIndex> > matching_peaks_
protected

Indices of peaks in the adjacent scans matching peaks in the scan with no. ref_scan.

◆ max_iteration_

UInt max_iteration_
protected

Convergence Parameter: Maximal number of iterations.

◆ max_peak_distance_

double max_peak_distance_
protected

upper bound for distance between two peaks belonging to the same region

◆ penalties_

Penalty factors for some parameters in the optimization.

◆ real_2D_

bool real_2D_
protected

Optimization considering all scans of a cluster or optimization of each scan separately.

◆ tolerance_mz_

double tolerance_mz_
protected

threshold for the difference in the peak position of two matching peaks