Home  · Classes  · Annotated Classes  · Modules  · Members  · Namespaces  · Related Pages
TraceFitter.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-2013.
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: Stephan Aiche$
32 // $Authors: Stephan Aiche, Marc Sturm $
33 // --------------------------------------------------------------------------
34 
35 #ifndef OPENMS_TRANSFORMATIONS_FEATUREFINDER_TRACEFITTER_H
36 #define OPENMS_TRANSFORMATIONS_FEATUREFINDER_TRACEFITTER_H
37 
39 
41 
43 
44 #include <gsl/gsl_rng.h>
45 #include <gsl/gsl_vector.h>
46 #include <gsl/gsl_multifit_nlin.h>
47 #include <gsl/gsl_blas.h>
48 
49 namespace OpenMS
50 {
51 
61  template <class PeakType>
62  class TraceFitter :
63  public DefaultParamHandler
64  {
65 
66 public:
69  DefaultParamHandler("TraceFitter")
70  {
71  this->defaults_.setValue("max_iteration", 500, "Maximum number of iterations using by Levenberg-Marquardt algorithm.", StringList::create("advanced"));
72  this->defaults_.setValue("epsilon_abs", 0.0001, "Absolute error used by the Levenberg-Marquardt algorithm.", StringList::create("advanced"));
73  this->defaults_.setValue("epsilon_rel", 0.0001, "Relative error used by the Levenberg-Marquardt algorithm.", StringList::create("advanced"));
74  }
75 
77  TraceFitter(const TraceFitter & source) :
78  DefaultParamHandler(source),
79  epsilon_abs_(source.epsilon_abs_),
80  epsilon_rel_(source.epsilon_rel_),
82  {
83  }
84 
86  virtual TraceFitter & operator=(const TraceFitter & source)
87  {
90  epsilon_abs_ = source.epsilon_abs_;
91  epsilon_rel_ = source.epsilon_rel_;
93 
94  return *this;
95  }
96 
98  virtual ~TraceFitter()
99  {
100  }
101 
106 
110  virtual DoubleReal getLowerRTBound() const = 0;
111 
115  virtual DoubleReal getUpperRTBound() const = 0;
116 
120  virtual DoubleReal getHeight() const = 0;
121 
125  virtual DoubleReal getCenter() const = 0;
126 
130  virtual DoubleReal getFWHM() const = 0;
131 
139 
146  virtual bool checkMinimalRTSpan(const std::pair<DoubleReal, DoubleReal> & rt_bounds, const DoubleReal min_rt_span) = 0;
147 
153  virtual bool checkMaximalRTSpan(const DoubleReal max_rt_span) = 0;
154 
160 
170  virtual String getGnuplotFormula(FeatureFinderAlgorithmPickedHelperStructs::MassTrace<PeakType> const & trace, const char function_name, const DoubleReal baseline, const DoubleReal rt_shift) = 0;
171 
172 protected:
173 
180  virtual void printState_(SignedSize iter, gsl_multifit_fdfsolver * s) = 0;
181 
182  virtual void updateMembers_()
183  {
184  max_iterations_ = this->param_.getValue("max_iteration");
185  epsilon_abs_ = this->param_.getValue("epsilon_abs");
186  epsilon_rel_ = this->param_.getValue("epsilon_rel");
187  }
188 
194  virtual void getOptimizedParameters_(gsl_multifit_fdfsolver * s) = 0;
195 
199  void optimize_(FeatureFinderAlgorithmPickedHelperStructs::MassTraces<PeakType> & traces, const Size num_params, double x_init[],
200  Int (* residual)(const gsl_vector * x, void * params, gsl_vector * f),
201  Int (* jacobian)(const gsl_vector * x, void * params, gsl_matrix * J),
202  Int (* evaluate)(const gsl_vector * x, void * params, gsl_vector * f, gsl_matrix * J))
203  {
204  const gsl_multifit_fdfsolver_type * T;
205  gsl_multifit_fdfsolver * s;
206 
207  const size_t data_count = traces.getPeakCount();
208 
209  // gsl always expects N>=p or default gsl error handler invoked,
210  // cause Jacobian be rectangular M x N with M>=N
211  if (data_count < num_params) throw Exception::UnableToFit(__FILE__, __LINE__, __PRETTY_FUNCTION__, "UnableToFit-FinalSet", "Skipping feature, gsl always expects N>=p");
212 
213  gsl_multifit_function_fdf func;
214  gsl_vector_view x = gsl_vector_view_array(x_init, num_params);
215  gsl_rng_env_setup();
216  func.f = (residual);
217  func.df = (jacobian);
218  func.fdf = (evaluate);
219  func.n = data_count;
220  func.p = num_params;
221  func.params = &traces;
222  T = gsl_multifit_fdfsolver_lmsder;
223  s = gsl_multifit_fdfsolver_alloc(T, data_count, num_params);
224  gsl_multifit_fdfsolver_set(s, &func, &x.vector);
225  SignedSize iter = 0;
226  Int gsl_status_;
227  do
228  {
229  iter++;
230  gsl_status_ = gsl_multifit_fdfsolver_iterate(s);
231  printState_(iter, s);
232  if (gsl_status_) break;
233  gsl_status_ = gsl_multifit_test_delta(s->dx, s->x, epsilon_abs_, epsilon_rel_);
234  }
235  while (gsl_status_ == GSL_CONTINUE && iter < max_iterations_);
236 
237  // get the parameters out of the fdfsolver
239 
240  gsl_multifit_fdfsolver_free(s);
241  }
242 
250 
251  };
252 
253 }
254 
255 #endif // #ifndef OPENMS_TRANSFORMATIONS_FEATUREFINDER_FEATUREFINDERALGORITHMPICKED_RTFITTING_H
const double k
virtual DoubleReal getFeatureIntensityContribution()=0
virtual String getGnuplotFormula(FeatureFinderAlgorithmPickedHelperStructs::MassTrace< PeakType > const &trace, const char function_name, const DoubleReal baseline, const DoubleReal rt_shift)=0
Param defaults_
Container for default parameters. This member should be filled in the constructor of derived classes!...
Definition: DefaultParamHandler.h:155
void setValue(const String &key, const DataValue &value, const String &description="", const StringList &tags=StringList())
Sets a value.
A more convenient string class.
Definition: String.h:56
DoubleReal epsilon_rel_
Relative error.
Definition: TraceFitter.h:247
void optimize_(FeatureFinderAlgorithmPickedHelperStructs::MassTraces< PeakType > &traces, const Size num_params, double x_init[], Int(*residual)(const gsl_vector *x, void *params, gsl_vector *f), Int(*jacobian)(const gsl_vector *x, void *params, gsl_matrix *J), Int(*evaluate)(const gsl_vector *x, void *params, gsl_vector *f, gsl_matrix *J))
Definition: TraceFitter.h:199
virtual void updateMembers_()
This method is used to update extra member variables at the end of the setParameters() method...
Definition: TraceFitter.h:182
Param param_
Container for current parameters.
Definition: DefaultParamHandler.h:148
Helper struct for a collection of mass traces used in FeatureFinderAlgorithmPicked.
Definition: FeatureFinderAlgorithmPickedHelperStructs.h:153
Abstract fitter for RT profile fitting.
Definition: TraceFitter.h:62
virtual bool checkMinimalRTSpan(const std::pair< DoubleReal, DoubleReal > &rt_bounds, const DoubleReal min_rt_span)=0
virtual DoubleReal getLowerRTBound() const =0
virtual void getOptimizedParameters_(gsl_multifit_fdfsolver *s)=0
ptrdiff_t SignedSize
Signed Size type e.g. used as pointer difference.
Definition: Types.h:151
Size getPeakCount() const
Returns the peak count of all traces.
Definition: FeatureFinderAlgorithmPickedHelperStructs.h:163
virtual DoubleReal computeTheoretical(const FeatureFinderAlgorithmPickedHelperStructs::MassTrace< PeakType > &trace, Size k)=0
int evaluate(const gsl_vector *x, void *params, gsl_vector *f, gsl_matrix *J)
Driver function for the evaluation of function and jacobian.
virtual TraceFitter & operator=(const TraceFitter &source)
assignment operator
Definition: TraceFitter.h:86
const DataValue & getValue(const String &key) const
Returns a value of a parameter.
TraceFitter(const TraceFitter &source)
copy constructor
Definition: TraceFitter.h:77
virtual bool checkMaximalRTSpan(const DoubleReal max_rt_span)=0
virtual DoubleReal getHeight() const =0
SignedSize max_iterations_
Maximum number of iterations.
Definition: TraceFitter.h:249
virtual DoubleReal getFWHM() const =0
virtual DoubleReal getUpperRTBound() const =0
static StringList create(const String &list, const char splitter= ',')
Returns a list that is created by splitting the given (comma-separated) string (String are not trimme...
TraceFitter()
default constructor.
Definition: TraceFitter.h:68
int residual(const gsl_vector *x, void *params, gsl_vector *f)
Evaluation of the target function for nonlinear optimization.
virtual void fit(FeatureFinderAlgorithmPickedHelperStructs::MassTraces< PeakType > &traces)=0
virtual ~TraceFitter()
destructor
Definition: TraceFitter.h:98
virtual void printState_(SignedSize iter, gsl_multifit_fdfsolver *s)=0
Exception used if an error occurred while fitting a model to a given dataset.
Definition: Exception.h:662
size_t Size
Size type e.g. used as variable which can hold result of size()
Definition: Types.h:144
int jacobian(const gsl_vector *x, void *params, gsl_matrix *J)
Compute the Jacobian of the residual, where each row of the matrix corresponds to a point in the data...
DoubleReal epsilon_abs_
Absolute error.
Definition: TraceFitter.h:245
virtual DefaultParamHandler & operator=(const DefaultParamHandler &rhs)
Assignment operator.
A base class for all classes handling default parameters.
Definition: DefaultParamHandler.h:90
virtual DoubleReal getCenter() const =0
int Int
Signed integer type.
Definition: Types.h:100
Helper struct for mass traces used in FeatureFinderAlgorithmPicked.
Definition: FeatureFinderAlgorithmPickedHelperStructs.h:83

OpenMS / TOPP release 1.11.1 Documentation generated on Thu Nov 14 2013 11:19:22 using doxygen 1.8.5