OpenMS
EmgFitter1D.h
Go to the documentation of this file.
1 // Copyright (c) 2002-present, 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: $
7 // --------------------------------------------------------------------------
8 
9 #pragma once
10 
13 
14 namespace OpenMS
15 {
21  class OPENMS_DLLAPI EmgFitter1D :
22  public LevMarqFitter1D
23  {
24 public:
25 
28 
30  EmgFitter1D(const EmgFitter1D& source);
31 
33  ~EmgFitter1D() override;
34 
36  virtual EmgFitter1D& operator=(const EmgFitter1D& source);
37 
39  QualityType fit1d(const RawDataArrayType& range, std::unique_ptr<InterpolationModel>& model) override;
40 
41 protected:
43  struct Data
44  {
45  typedef Peak1D PeakType;
46  typedef std::vector<PeakType> RawDataArrayType;
47 
50  };
51 
54  {
55 public:
56  EgmFitterFunctor(int dimensions, const EmgFitter1D::Data* data) :
57  LevMarqFitter1D::GenericFunctor(dimensions,
58  static_cast<int>(data->n)),
59  m_data(data)
60  {}
61 
62  int operator()(const Eigen::VectorXd& x, Eigen::VectorXd& fvec) const override;
63  // compute Jacobian matrix for the different parameters
64  int df(const Eigen::VectorXd& x, Eigen::MatrixXd& J) const override;
65 
66 protected:
72  };
73 
75  virtual void setInitialParameters_(const RawDataArrayType& set);
79 
88 
89  void updateMembers_() override;
90  };
91 
92 }
93 
Definition: IsobaricIsotopeCorrector.h:17
Definition: EmgFitter1D.h:54
static const EmgFitter1D::CoordinateType sqrt_2
Definition: EmgFitter1D.h:71
int df(const Eigen::VectorXd &x, Eigen::MatrixXd &J) const override
const EmgFitter1D::Data * m_data
Definition: EmgFitter1D.h:67
static const EmgFitter1D::CoordinateType emg_const
Definition: EmgFitter1D.h:70
int operator()(const Eigen::VectorXd &x, Eigen::VectorXd &fvec) const override
EgmFitterFunctor(int dimensions, const EmgFitter1D::Data *data)
Definition: EmgFitter1D.h:56
static const EmgFitter1D::CoordinateType c
Definition: EmgFitter1D.h:68
static const EmgFitter1D::CoordinateType sqrt2pi
Definition: EmgFitter1D.h:69
Exponentially modified gaussian distribution fitter (1-dim.) using Levenberg-Marquardt algorithm (Eig...
Definition: EmgFitter1D.h:23
~EmgFitter1D() override
destructor
std::vector< PeakType > RawDataArrayType
Definition: EmgFitter1D.h:46
virtual void setInitialParameters_(const RawDataArrayType &set)
Compute start parameter.
EmgFitter1D()
Default constructor.
RawDataArrayType set
Definition: EmgFitter1D.h:49
Size n
Definition: EmgFitter1D.h:48
CoordinateType height_
Parameter of emg - peak height.
Definition: EmgFitter1D.h:81
QualityType fit1d(const RawDataArrayType &range, std::unique_ptr< InterpolationModel > &model) override
return interpolation model
CoordinateType retention_
Parameter of emg - peak retention time.
Definition: EmgFitter1D.h:87
virtual EmgFitter1D & operator=(const EmgFitter1D &source)
assignment operator
void updateMembers_() override
This method is used to update extra member variables at the end of the setParameters() method.
Peak1D PeakType
Definition: EmgFitter1D.h:45
void setInitialParametersMOM_(const RawDataArrayType &set)
EmgFitter1D(const EmgFitter1D &source)
copy constructor
CoordinateType width_
Parameter of emg - peak width.
Definition: EmgFitter1D.h:83
CoordinateType symmetry_
Parameter of emg - peak symmetry.
Definition: EmgFitter1D.h:85
Helper struct (contains the size of an area and a raw data container)
Definition: EmgFitter1D.h:44
std::vector< PeakType > RawDataArrayType
Peak type data container type using for the temporary storage of the input data.
Definition: Fitter1D.h:43
Feature::QualityType QualityType
Quality of a feature.
Definition: Fitter1D.h:39
Feature::CoordinateType CoordinateType
Single coordinate.
Definition: Fitter1D.h:37
Definition: LevMarqFitter1D.h:38
Abstract class for 1D-model fitter using Levenberg-Marquardt algorithm for parameter optimization.
Definition: LevMarqFitter1D.h:31
A 1-dimensional raw data point or peak.
Definition: Peak1D.h:28
size_t Size
Size type e.g. used as variable which can hold result of size()
Definition: Types.h:97
Main OpenMS namespace.
Definition: openswathalgo/include/OpenMS/OPENSWATHALGO/DATAACCESS/ISpectrumAccess.h:19