OpenMS  2.8.0
TwoDOptimization.h
Go to the documentation of this file.
1 // --------------------------------------------------------------------------
2 // OpenMS -- Open-Source Mass Spectrometry
3 // --------------------------------------------------------------------------
4 // Copyright The OpenMS Team -- Eberhard Karls University Tuebingen,
5 // ETH Zurich, and Freie Universitaet Berlin 2002-2021.
6 //
7 // This software is released under a three-clause BSD license:
8 // * Redistributions of source code must retain the above copyright
9 // notice, this list of conditions and the following disclaimer.
10 // * Redistributions in binary form must reproduce the above copyright
11 // notice, this list of conditions and the following disclaimer in the
12 // documentation and/or other materials provided with the distribution.
13 // * Neither the name of any author or any participating institution
14 // may be used to endorse or promote products derived from this software
15 // without specific prior written permission.
16 // For a full list of authors, refer to the file AUTHORS.
17 // --------------------------------------------------------------------------
18 // THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
19 // AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
20 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
21 // ARE DISCLAIMED. IN NO EVENT SHALL ANY OF THE AUTHORS OR THE CONTRIBUTING
22 // INSTITUTIONS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
23 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
24 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS;
25 // OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY,
26 // WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR
27 // OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF
28 // ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
29 //
30 // --------------------------------------------------------------------------
31 // $Maintainer: Timo Sachsenberg $
32 // $Authors: $
33 // --------------------------------------------------------------------------
34 
35 #pragma once
36 
37 //#define DEBUG_2D
38 #undef DEBUG_2D
39 
40 #ifdef DEBUG_2D
41 #include <iostream>
42 #include <fstream>
43 #endif
44 
45 #include <vector>
46 #include <utility>
47 #include <cmath>
48 #include <set>
49 
60 
61 
62 #ifndef OPENMS_SYSTEM_STOPWATCH_H
63 #endif
64 
67 
68 namespace OpenMS
69 {
70 
85  class OPENMS_DLLAPI TwoDOptimization :
86  public DefaultParamHandler
87  {
88 public:
89 
91 
94 
97 
99  ~TwoDOptimization() override{}
100 
103 
104 
106  inline double getMZTolerance() const {return tolerance_mz_; }
108  inline void setMZTolerance(double tolerance_mz)
109  {
110  tolerance_mz_ = tolerance_mz;
111  param_.setValue("2d:tolerance_mz", tolerance_mz);
112  }
113 
115  inline double getMaxPeakDistance() const {return max_peak_distance_; }
117  inline void setMaxPeakDistance(double max_peak_distance)
118  {
119  max_peak_distance_ = max_peak_distance;
120  param_.setValue("2d:max_peak_distance", max_peak_distance);
121  }
122 
124  inline UInt getMaxIterations() const {return max_iteration_; }
126  inline void setMaxIterations(UInt max_iteration)
127  {
128  max_iteration_ = max_iteration;
129  param_.setValue("iterations", max_iteration);
130  }
131 
133  inline const OptimizationFunctions::PenaltyFactorsIntensity& getPenalties() const {return penalties_; }
136  {
137  penalties_ = penalties;
138  param_.setValue("penalties:position", penalties.pos);
139  param_.setValue("penalties:height", penalties.height);
140  param_.setValue("penalties:left_width", penalties.lWidth);
141  param_.setValue("penalties:right_width", penalties.rWidth);
142  }
143 
162  PeakMap& ms_exp, bool real2D = true);
163 
164 
165 protected:
167  struct Data
168  {
169  std::vector<std::pair<SignedSize, SignedSize> > signal2D;
170  std::multimap<double, IsotopeCluster>::iterator iso_map_iter;
172  std::map<Int, std::vector<PeakIndex> > matching_peaks;
176  std::vector<double> positions;
177  std::vector<double> signal;
178  };
179 
180  class OPENMS_DLLAPI TwoDOptFunctor
181  {
182  public:
183  int inputs() const { return m_inputs; }
184  int values() const { return m_values; }
185 
186  TwoDOptFunctor(unsigned dimensions, unsigned num_data_points, const TwoDOptimization::Data* data)
187  : m_inputs(dimensions), m_values(num_data_points), m_data(data) {}
188 
190  // compute Jacobian matrix for the different parameters
191  int df(const Eigen::VectorXd &x, Eigen::MatrixXd &J);
192 
193  private:
194  const int m_inputs, m_values;
196  };
197 
198 
200  std::multimap<double, IsotopeCluster> iso_map_;
201 
203  std::multimap<double, IsotopeCluster>::const_iterator curr_region_;
204 
207 
210 
212  // std::map<Int, std::vector<PeakMap::SpectrumType::Iterator > > matching_peaks_;
213  std::map<Int, std::vector<PeakIndex> > matching_peaks_;
214 
215 
218 
220  bool real_2D_;
221 
222 
225 
226 
231  std::vector<double>::iterator searchInScan_(std::vector<double>::iterator scan_begin,
232  std::vector<double>::iterator scan_end,
233  double current_mz);
234 
237  InputSpectrumIterator& last,
238  PeakMap& ms_exp);
239 
242  InputSpectrumIterator& last,
243  PeakMap& ms_exp);
244 
245 
248  InputSpectrumIterator& first,
249  InputSpectrumIterator& last,
250  Size iso_map_idx,
251  double noise_level,
253 
255  void findMatchingPeaks_(std::multimap<double, IsotopeCluster>::iterator& it,
256  PeakMap& ms_exp);
257 
259 
261  void updateMembers_() override;
262  };
263 
264 }
265 
A base class for all classes handling default parameters.
Definition: DefaultParamHandler.h:93
In-Memory representation of a mass spectrometry run.
Definition: MSExperiment.h:73
Base::const_iterator const_iterator
Definition: MSExperiment.h:118
std::vector< SpectrumType >::const_iterator ConstIterator
Non-mutable iterator.
Definition: MSExperiment.h:106
Definition: TwoDOptimization.h:181
int values() const
Definition: TwoDOptimization.h:184
int df(const Eigen::VectorXd &x, Eigen::MatrixXd &J)
const int m_inputs
Definition: TwoDOptimization.h:194
int operator()(const Eigen::VectorXd &x, Eigen::VectorXd &fvec)
const TwoDOptimization::Data * m_data
Definition: TwoDOptimization.h:195
TwoDOptFunctor(unsigned dimensions, unsigned num_data_points, const TwoDOptimization::Data *data)
Definition: TwoDOptimization.h:186
int inputs() const
Definition: TwoDOptimization.h:183
This class provides the two-dimensional optimization of the picked peak parameters.
Definition: TwoDOptimization.h:87
double getMZTolerance() const
Non-mutable access to the matching epsilon.
Definition: TwoDOptimization.h:106
std::multimap< double, IsotopeCluster >::iterator iso_map_iter
Definition: TwoDOptimization.h:170
std::vector< std::pair< SignedSize, SignedSize > > signal2D
Definition: TwoDOptimization.h:169
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.
OptimizationFunctions::PenaltyFactorsIntensity penalties_
Penalty factors for some parameters in the optimization.
Definition: TwoDOptimization.h:224
void setPenalties(OptimizationFunctions::PenaltyFactorsIntensity &penalties)
Mutable access to the minimal number of adjacent scans.
Definition: TwoDOptimization.h:135
PeakMap::ConstIterator raw_data_first
Definition: TwoDOptimization.h:174
std::map< Int, std::vector< PeakIndex > > matching_peaks
Definition: TwoDOptimization.h:172
PeakMap picked_peaks
Definition: TwoDOptimization.h:173
OptimizationFunctions::PenaltyFactorsIntensity penalties
Definition: TwoDOptimization.h:175
void findMatchingPeaks_(std::multimap< double, IsotopeCluster >::iterator &it, PeakMap &ms_exp)
Identify matching peak in a peak cluster.
std::map< Int, std::vector< PeakIndex > > matching_peaks_
Indices of peaks in the adjacent scans matching peaks in the scan with no. ref_scan.
Definition: TwoDOptimization.h:213
void optimizeRegions_(InputSpectrumIterator &first, InputSpectrumIterator &last, PeakMap &ms_exp)
TwoDOptimization()
Constructor.
void optimizeRegionsScanwise_(InputSpectrumIterator &first, InputSpectrumIterator &last, PeakMap &ms_exp)
bool real_2D_
Optimization considering all scans of a cluster or optimization of each scan separately.
Definition: TwoDOptimization.h:220
std::vector< double > positions
Definition: TwoDOptimization.h:176
MSExperiment::const_iterator InputSpectrumIterator
Definition: TwoDOptimization.h:90
double max_peak_distance_
upper bound for distance between two peaks belonging to the same region
Definition: TwoDOptimization.h:206
double tolerance_mz_
threshold for the difference in the peak position of two matching peaks
Definition: TwoDOptimization.h:209
UInt max_iteration_
Convergence Parameter: Maximal number of iterations.
Definition: TwoDOptimization.h:217
double getMaxPeakDistance() const
Non-mutable access to the maximal peak distance in a cluster.
Definition: TwoDOptimization.h:115
UInt getMaxIterations() const
Non-mutable access to the maximal number of iterations.
Definition: TwoDOptimization.h:124
Size total_nr_peaks
Definition: TwoDOptimization.h:171
std::vector< double >::iterator searchInScan_(std::vector< double >::iterator scan_begin, std::vector< double >::iterator scan_end, double current_mz)
std::multimap< double, IsotopeCluster > iso_map_
stores the retention time of each isotopic cluster
Definition: TwoDOptimization.h:200
void setMaxPeakDistance(double max_peak_distance)
Mutable access to the maximal peak distance in a cluster.
Definition: TwoDOptimization.h:117
TwoDOptimization(const TwoDOptimization &opt)
Copy constructor.
TwoDOptimization & operator=(const TwoDOptimization &opt)
Assignment operator.
void updateMembers_() override
update members method from DefaultParamHandler to update the members
std::multimap< double, IsotopeCluster >::const_iterator curr_region_
Pointer to the current region.
Definition: TwoDOptimization.h:203
void optimize(InputSpectrumIterator first, InputSpectrumIterator last, PeakMap &ms_exp, bool real2D=true)
Find two dimensional peak clusters and optimize their peak parameters.
void setMZTolerance(double tolerance_mz)
Mutable access to the matching epsilon.
Definition: TwoDOptimization.h:108
const OptimizationFunctions::PenaltyFactorsIntensity & getPenalties() const
Non-mutable access to the minimal number of adjacent scans.
Definition: TwoDOptimization.h:133
void setMaxIterations(UInt max_iteration)
Mutable access to the maximal number of iterations.
Definition: TwoDOptimization.h:126
~TwoDOptimization() override
Destructor.
Definition: TwoDOptimization.h:99
std::vector< double > signal
Definition: TwoDOptimization.h:177
Helper struct (contains the size of an area and a raw data container)
Definition: TwoDOptimization.h:168
unsigned int UInt
Unsigned integer type.
Definition: Types.h:94
size_t Size
Size type e.g. used as variable which can hold result of size()
Definition: Types.h:127
Definition: IsobaricIsotopeCorrector.h:43
Main OpenMS namespace.
Definition: FeatureDeconvolution.h:47
Class for the penalty factors used during the optimization.
Definition: OptimizePeakDeconvolution.h:59
double height
Definition: OptimizePeakDeconvolution.h:76
double rWidth
Penalty factor for the peak shape's right width parameter.
Definition: OptimizePick.h:90
double pos
Penalty factor for the peak shape's position.
Definition: OptimizePick.h:86
double lWidth
Penalty factor for the peak shape's left width parameter.
Definition: OptimizePick.h:88