OpenMS
TraceFitter.h
Go to the documentation of this file.
1 // Copyright (c) 2002-2023, The OpenMS Team -- EKU Tuebingen, ETH Zurich, and FU Berlin
2 // SPDX-License-Identifier: BSD-3-Clause
3 //
4 // --------------------------------------------------------------------------
5 // $Maintainer: Timo Sachsenberg$
6 // $Authors: Stephan Aiche, Marc Sturm $
7 // --------------------------------------------------------------------------
8 
9 #pragma once
10 
13 #include <OpenMS/KERNEL/Peak1D.h>
14 
15 // forward decl
16 namespace Eigen
17 {
18  template<typename _Scalar, int _Rows, int _Cols, int _Options, int _MaxRows, int _MaxCols>
19  class Matrix;
20  using MatrixXd = Matrix<double, -1, -1, 0, -1, -1>;
21  using VectorXd = Matrix<double, -1, 1, 0, -1, 1>;
22 }
23 
24 namespace OpenMS
25 {
26 
36  class OPENMS_DLLAPI TraceFitter :
37  public DefaultParamHandler
38  {
39 
40 public:
42  //TODO: This is copy and paste from LevMarqFitter1d.h. Make a generic wrapper for LM optimization
44  {
45 public:
46  int inputs() const;
47  int values() const;
48 
49  GenericFunctor(int dimensions, int num_data_points);
50 
51  virtual ~GenericFunctor();
52 
53  virtual int operator()(const Eigen::VectorXd& x, Eigen::VectorXd& fvec) = 0;
54  // compute Jacobian matrix for the different parameters
55  virtual int df(const Eigen::VectorXd& x, Eigen::MatrixXd& J) = 0;
56 
57 protected:
58  const int m_inputs, m_values;
59  };
60 
63 
65  TraceFitter(const TraceFitter& source);
66 
69 
71  ~TraceFitter() override;
72 
77 
81  virtual double getLowerRTBound() const = 0;
82 
86  virtual double getUpperRTBound() const = 0;
87 
91  virtual double getHeight() const = 0;
92 
96  virtual double getCenter() const = 0;
97 
101  virtual double getFWHM() const = 0;
102 
106  virtual double getValue(double rt) const = 0;
107 
115 
122  virtual bool checkMinimalRTSpan(const std::pair<double, double>& rt_bounds, const double min_rt_span) = 0;
123 
129  virtual bool checkMaximalRTSpan(const double max_rt_span) = 0;
130 
134  virtual double getArea() = 0;
135 
145  virtual String getGnuplotFormula(const FeatureFinderAlgorithmPickedHelperStructs::MassTrace& trace, const char function_name, const double baseline, const double rt_shift) = 0;
146 
147 protected:
148  struct ModelData
149  {
151  bool weighted;
152  };
153 
154  void updateMembers_() override;
155 
161  virtual void getOptimizedParameters_(const Eigen::VectorXd& s) = 0;
165  void optimize_(Eigen::VectorXd& x_init, GenericFunctor& functor);
166 
170  bool weighted_;
171 
172  };
173 
174 }
175 
A base class for all classes handling default parameters.
Definition: DefaultParamHandler.h:66
A more convenient string class.
Definition: String.h:34
Definition: TraceFitter.h:44
virtual int df(const Eigen::VectorXd &x, Eigen::MatrixXd &J)=0
GenericFunctor(int dimensions, int num_data_points)
virtual int operator()(const Eigen::VectorXd &x, Eigen::VectorXd &fvec)=0
const int m_inputs
Definition: TraceFitter.h:58
Abstract fitter for RT profile fitting.
Definition: TraceFitter.h:38
TraceFitter()
default constructor
virtual void getOptimizedParameters_(const Eigen::VectorXd &s)=0
virtual double getLowerRTBound() const =0
virtual void fit(FeatureFinderAlgorithmPickedHelperStructs::MassTraces &traces)=0
FeatureFinderAlgorithmPickedHelperStructs::MassTraces * traces_ptr
Definition: TraceFitter.h:150
bool weighted_
Whether to weight mass traces by theoretical intensity during the optimization.
Definition: TraceFitter.h:170
virtual double getHeight() const =0
void optimize_(Eigen::VectorXd &x_init, GenericFunctor &functor)
virtual String getGnuplotFormula(const FeatureFinderAlgorithmPickedHelperStructs::MassTrace &trace, const char function_name, const double baseline, const double rt_shift)=0
bool weighted
Definition: TraceFitter.h:151
virtual double getValue(double rt) const =0
virtual double getCenter() const =0
TraceFitter & operator=(const TraceFitter &source)
assignment operator
double computeTheoretical(const FeatureFinderAlgorithmPickedHelperStructs::MassTrace &trace, Size k) const
virtual double getUpperRTBound() const =0
virtual double getFWHM() const =0
TraceFitter(const TraceFitter &source)
copy constructor
void updateMembers_() override
This method is used to update extra member variables at the end of the setParameters() method.
virtual bool checkMaximalRTSpan(const double max_rt_span)=0
virtual double getArea()=0
SignedSize max_iterations_
Maximum number of iterations.
Definition: TraceFitter.h:168
~TraceFitter() override
destructor
virtual bool checkMinimalRTSpan(const std::pair< double, double > &rt_bounds, const double min_rt_span)=0
Definition: TraceFitter.h:149
ptrdiff_t SignedSize
Signed Size type e.g. used as pointer difference.
Definition: Types.h:108
size_t Size
Size type e.g. used as variable which can hold result of size()
Definition: Types.h:101
Definition: IsobaricIsotopeCorrector.h:15
Matrix< double, -1, -1, 0, -1, -1 > MatrixXd
Definition: IsobaricIsotopeCorrector.h:18
Matrix< double, -1, 1, 0, -1, 1 > VectorXd
Definition: IsobaricIsotopeCorrector.h:19
Definition: IsobaricIsotopeCorrector.h:17
const double k
Definition: Constants.h:132
Main OpenMS namespace.
Definition: FeatureDeconvolution.h:22
Helper struct for mass traces used in FeatureFinderAlgorithmPicked.
Definition: FeatureFinderAlgorithmPickedHelperStructs.h:54
Helper struct for a collection of mass traces used in FeatureFinderAlgorithmPicked.
Definition: FeatureFinderAlgorithmPickedHelperStructs.h:85