OpenMS  2.8.0
LevMarqFitter1D.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 
38 
39 
40 #include <algorithm>
41 
43 
44 // forward decl
45 namespace Eigen
46 {
47  template<typename _Scalar, int _Rows, int _Cols, int _Options, int _MaxRows, int _MaxCols>
48  class Matrix;
49  using MatrixXd = Matrix<double, -1, -1, 0, -1, -1>;
50  using VectorXd = Matrix<double, -1, 1, 0, -1, 1>;
51 }
52 
53 namespace OpenMS
54 {
55 
59  class OPENMS_DLLAPI LevMarqFitter1D :
60  public Fitter1D
61  {
62 
63 public:
64 
65  typedef std::vector<double> ContainerType;
66 
68  //TODO: This is copy and paste from TraceFitter.h. Make a generic wrapper for LM optimization
70  {
71  public:
72  int inputs() const { return m_inputs; }
73  int values() const { return m_values; }
74 
75  GenericFunctor(int dimensions, int num_data_points)
76  : m_inputs(dimensions), m_values(num_data_points) {}
77 
78  virtual ~GenericFunctor() {}
79 
80  virtual int operator()(const Eigen::VectorXd &x, Eigen::VectorXd &fvec) const = 0;
81 
82  // compute Jacobian matrix for the different parameters
83  virtual int df(const Eigen::VectorXd &x, Eigen::MatrixXd &J) const = 0;
84 
85  protected:
86  const int m_inputs, m_values;
87  };
88 
91  Fitter1D()
92  {
93  this->defaults_.setValue("max_iteration", 500, "Maximum number of iterations using by Levenberg-Marquardt algorithm.", {"advanced"});
94  }
95 
97  LevMarqFitter1D(const LevMarqFitter1D & source) :
98  Fitter1D(source),
99  max_iteration_(source.max_iteration_)
100  {
101  }
102 
104  ~LevMarqFitter1D() override
105  {
106  }
107 
109  virtual LevMarqFitter1D & operator=(const LevMarqFitter1D & source)
110  {
111  if (&source == this) return *this;
112 
113  Fitter1D::operator=(source);
114  max_iteration_ = source.max_iteration_;
115 
116  return *this;
117  }
118 
119 protected:
120 
125 
131  void optimize_(Eigen::VectorXd& x_init, GenericFunctor& functor) const;
132 
133  void updateMembers_() override;
134 
135  };
136 }
137 
Abstract base class for all 1D-dimensional model fitter.
Definition: Fitter1D.h:62
virtual Fitter1D & operator=(const Fitter1D &source)
assignment operator
Definition: LevMarqFitter1D.h:70
virtual ~GenericFunctor()
Definition: LevMarqFitter1D.h:78
virtual int df(const Eigen::VectorXd &x, Eigen::MatrixXd &J) const =0
int values() const
Definition: LevMarqFitter1D.h:73
GenericFunctor(int dimensions, int num_data_points)
Definition: LevMarqFitter1D.h:75
virtual int operator()(const Eigen::VectorXd &x, Eigen::VectorXd &fvec) const =0
const int m_inputs
Definition: LevMarqFitter1D.h:86
int inputs() const
Definition: LevMarqFitter1D.h:72
Abstract class for 1D-model fitter using Levenberg-Marquardt algorithm for parameter optimization.
Definition: LevMarqFitter1D.h:61
LevMarqFitter1D()
Default constructor.
Definition: LevMarqFitter1D.h:90
Int max_iteration_
Maximum number of iterations.
Definition: LevMarqFitter1D.h:124
std::vector< double > ContainerType
Definition: LevMarqFitter1D.h:65
bool symmetric_
Parameter indicates symmetric peaks.
Definition: LevMarqFitter1D.h:122
void updateMembers_() override
This method is used to update extra member variables at the end of the setParameters() method.
~LevMarqFitter1D() override
destructor
Definition: LevMarqFitter1D.h:104
virtual LevMarqFitter1D & operator=(const LevMarqFitter1D &source)
assignment operator
Definition: LevMarqFitter1D.h:109
void optimize_(Eigen::VectorXd &x_init, GenericFunctor &functor) const
Optimize start parameter.
LevMarqFitter1D(const LevMarqFitter1D &source)
copy constructor
Definition: LevMarqFitter1D.h:97
int Int
Signed integer type.
Definition: Types.h:102
Definition: IsobaricIsotopeCorrector.h:41
Matrix< double, -1, -1, 0, -1, -1 > MatrixXd
Definition: IsobaricIsotopeCorrector.h:44
Matrix< double, -1, 1, 0, -1, 1 > VectorXd
Definition: IsobaricIsotopeCorrector.h:45
Definition: IsobaricIsotopeCorrector.h:43
Main OpenMS namespace.
Definition: FeatureDeconvolution.h:47