OpenMS
Loading...
Searching...
No Matches
TraceFitter.h
Go to the documentation of this file.
1// Copyright (c) 2002-present, OpenMS Inc. -- 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
14#include <vector>
15
16namespace OpenMS
17{
18
28 class OPENMS_DLLAPI TraceFitter :
30 {
31
32public:
37 //TODO: This is copy and paste from LevMarqFitter1d.h. Make a generic wrapper for LM optimization
39 {
40public:
41 int inputs() const;
42 int values() const;
43
44 GenericFunctor(int dimensions, int num_data_points);
45
46 virtual ~GenericFunctor();
47
49 virtual int operator()(const double* x, double* fvec) = 0;
51 virtual int df(const double* x, double* J) = 0;
52
53protected:
54 const int m_inputs, m_values;
55 };
56
59
61 TraceFitter(const TraceFitter& source);
62
65
67 ~TraceFitter() override;
68
73
77 virtual double getLowerRTBound() const = 0;
78
82 virtual double getUpperRTBound() const = 0;
83
87 virtual double getHeight() const = 0;
88
92 virtual double getCenter() const = 0;
93
97 virtual double getFWHM() const = 0;
98
102 virtual double getValue(double rt) const = 0;
103
111
118 virtual bool checkMinimalRTSpan(const std::pair<double, double>& rt_bounds, const double min_rt_span) = 0;
119
125 virtual bool checkMaximalRTSpan(const double max_rt_span) = 0;
126
130 virtual double getArea() = 0;
131
141 virtual String getGnuplotFormula(const FeatureFinderAlgorithmPickedHelperStructs::MassTrace& trace, const char function_name, const double baseline, const double rt_shift) = 0;
142
143protected:
149
150 void updateMembers_() override;
151
157 virtual void getOptimizedParameters_(const std::vector<double>& s) = 0;
161 void optimize_(std::vector<double>& x_init, GenericFunctor& functor);
162
167
168 };
169
170}
171
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:39
virtual int df(const double *x, double *J)=0
Compute Jacobian matrix. x has size inputs(), J is values() x inputs() (column-major)
GenericFunctor(int dimensions, int num_data_points)
const int m_inputs
Definition TraceFitter.h:54
virtual int operator()(const double *x, double *fvec)=0
Compute residuals. x has size inputs(), fvec has size values()
Abstract fitter for RT profile fitting.
Definition TraceFitter.h:30
TraceFitter()
default constructor
virtual double getLowerRTBound() const =0
virtual void fit(FeatureFinderAlgorithmPickedHelperStructs::MassTraces &traces)=0
FeatureFinderAlgorithmPickedHelperStructs::MassTraces * traces_ptr
Definition TraceFitter.h:146
bool weighted_
Whether to weight mass traces by theoretical intensity during the optimization.
Definition TraceFitter.h:166
virtual double getHeight() const =0
virtual String getGnuplotFormula(const FeatureFinderAlgorithmPickedHelperStructs::MassTrace &trace, const char function_name, const double baseline, const double rt_shift)=0
bool weighted
Definition TraceFitter.h:147
virtual double getValue(double rt) const =0
virtual double getCenter() const =0
virtual void getOptimizedParameters_(const std::vector< double > &s)=0
double computeTheoretical(const FeatureFinderAlgorithmPickedHelperStructs::MassTrace &trace, Size k) const
void optimize_(std::vector< double > &x_init, GenericFunctor &functor)
virtual double getUpperRTBound() const =0
TraceFitter & operator=(const TraceFitter &source)
assignment operator
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:164
~TraceFitter() override
destructor
virtual bool checkMinimalRTSpan(const std::pair< double, double > &rt_bounds, const double min_rt_span)=0
Definition TraceFitter.h:145
ptrdiff_t SignedSize
Signed Size type e.g. used as pointer difference.
Definition Types.h:104
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
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