OpenMS  2.7.0
OptimizePick.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: Eva Lange $
33 // --------------------------------------------------------------------------
34 
35 #pragma once
36 
38 #include <OpenMS/KERNEL/Peak1D.h>
39 
40 #include <Eigen/Core>
41 #include <unsupported/Eigen/NonLinearOptimization>
42 
43 #include <vector>
44 
45 namespace OpenMS
46 {
47  namespace OptimizationFunctions
48  {
50  typedef std::vector<Peak1D> RawDataVector;
52  typedef RawDataVector::iterator PeakIterator;
53 
62  struct OPENMS_DLLAPI PenaltyFactors
63  {
65  pos(0), lWidth(0), rWidth(0) {}
67  pos(p.pos), lWidth(p.lWidth), rWidth(p.rWidth) {}
69  {
70  pos = p.pos;
71  lWidth = p.lWidth;
72  rWidth = p.rWidth;
73 
74  return *this;
75  }
76 
78 
80  double pos;
82  double lWidth;
84  double rWidth;
85  };
86  }
87 
88 
95  class OPENMS_DLLAPI OptimizePick
96  {
97 public:
98 
99  struct Data
100  {
102  std::vector<double> positions;
103  std::vector<double> signal;
105  std::vector<PeakShape> peaks;
106 
108 
109  };
110 
112  {
113  public:
114  int inputs() const { return m_inputs; }
115  int values() const { return m_values; }
116 
117  OptPeakFunctor(unsigned dimensions, unsigned num_data_points, const OptimizePick::Data * data)
118  : m_inputs(dimensions), m_values(num_data_points), m_data(data) {}
119 
120  int operator()(const Eigen::VectorXd &x, Eigen::VectorXd &fvec);
121  // compute Jacobian matrix for the different parameters
122  int df(const Eigen::VectorXd &x, Eigen::MatrixXd &J);
123 
124  private:
125  const int m_inputs, m_values;
126  const Data * m_data;
127  };
128 
130  typedef std::vector<Peak1D> RawDataVector;
132  typedef RawDataVector::iterator PeakIterator;
133 
134 
137  max_iteration_(400)
138  {}
139 
142  const int max_iteration_);
143 
146 
148  inline const struct OptimizationFunctions::PenaltyFactors & getPenalties() const { return penalties_; }
150  inline struct OptimizationFunctions::PenaltyFactors & getPenalties() { return penalties_; }
152  inline void setPenalties(const struct OptimizationFunctions::PenaltyFactors & penalties) { penalties_ = penalties; }
153 
155  inline UInt getNumberIterations() const { return max_iteration_; }
157  inline unsigned int & getNumberIterations() { return max_iteration_; }
159  inline void setNumberIterations(const int max_iteration) { max_iteration_ = max_iteration; }
160 
162  void optimize(std::vector<PeakShape> & peaks, Data & data);
163 
164 
165 protected:
167  struct OptimizationFunctions::PenaltyFactors penalties_;
168 
170  unsigned int max_iteration_;
171  };
172 }
173 
Definition: OptimizePick.h:112
OptPeakFunctor(unsigned dimensions, unsigned num_data_points, const OptimizePick::Data *data)
Definition: OptimizePick.h:117
int values() const
Definition: OptimizePick.h:115
int df(const Eigen::VectorXd &x, Eigen::MatrixXd &J)
const Data * m_data
Definition: OptimizePick.h:126
const int m_inputs
Definition: OptimizePick.h:125
int operator()(const Eigen::VectorXd &x, Eigen::VectorXd &fvec)
int inputs() const
Definition: OptimizePick.h:114
This class provides the non-linear optimization of the peak parameters.
Definition: OptimizePick.h:96
void setNumberIterations(const int max_iteration)
Mutable access to the number of iterations.
Definition: OptimizePick.h:159
void setPenalties(const struct OptimizationFunctions::PenaltyFactors &penalties)
Mutable access to the penalty factors.
Definition: OptimizePick.h:152
std::vector< Peak1D > RawDataVector
Profile data vector type.
Definition: OptimizePick.h:130
OptimizationFunctions::PenaltyFactors penalties
Definition: OptimizePick.h:107
OptimizePick()
Constructor.
Definition: OptimizePick.h:136
RawDataVector::iterator PeakIterator
Profile data iterator type.
Definition: OptimizePick.h:132
unsigned int max_iteration_
Maximum number of iterations during optimization.
Definition: OptimizePick.h:170
std::vector< double > positions
Positions and intensity values of the profile data.
Definition: OptimizePick.h:102
OptimizePick(const struct OptimizationFunctions::PenaltyFactors &penalties_, const int max_iteration_)
Constructor to set the penalty factors, the number of optimization iterations as well as the threshol...
unsigned int & getNumberIterations()
Mutable access to the number of iterations.
Definition: OptimizePick.h:157
struct OptimizationFunctions::PenaltyFactors & getPenalties()
Mutable access to the penalty factors.
Definition: OptimizePick.h:150
std::vector< PeakShape > peaks
This container contains the peak shapes to be optimized.
Definition: OptimizePick.h:105
const struct OptimizationFunctions::PenaltyFactors & getPenalties() const
Non-mutable access to the penalty factors.
Definition: OptimizePick.h:148
~OptimizePick()
Destructor.
void optimize(std::vector< PeakShape > &peaks, Data &data)
Start the optimization of the peak shapes peaks. The original peak shapes will be substituted by the ...
UInt getNumberIterations() const
Non-mutable access to the number of iterations.
Definition: OptimizePick.h:155
std::vector< double > signal
Definition: OptimizePick.h:103
Definition: OptimizePick.h:100
unsigned int UInt
Unsigned integer type.
Definition: Types.h:94
std::vector< Peak1D > RawDataVector
Profile data vector type.
Definition: OptimizePick.h:50
RawDataVector::iterator PeakIterator
Profile data iterator type.
Definition: OptimizePick.h:52
Main OpenMS namespace.
Definition: FeatureDeconvolution.h:47
Class for the penalty factors used during the optimization.
Definition: OptimizePick.h:63
double rWidth
Penalty factor for the peak shape's right width parameter.
Definition: OptimizePick.h:84
double pos
Penalty factor for the peak shape's position.
Definition: OptimizePick.h:80
PenaltyFactors()
Definition: OptimizePick.h:64
double lWidth
Penalty factor for the peak shape's left width parameter.
Definition: OptimizePick.h:82
PenaltyFactors & operator=(const PenaltyFactors &p)
Definition: OptimizePick.h:68
PenaltyFactors(const PenaltyFactors &p)
Definition: OptimizePick.h:66
~PenaltyFactors()
Definition: OptimizePick.h:77