Home  · Classes  · Annotated Classes  · Modules  · Members  · Namespaces  · Related Pages
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-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: Clemens Groepl $
32 // $Authors: $
33 // --------------------------------------------------------------------------
34 
35 #ifndef OPENMS_TRANSFORMATIONS_FEATUREFINDER_LEVMARQFITTER1D_H
36 #define OPENMS_TRANSFORMATIONS_FEATUREFINDER_LEVMARQFITTER1D_H
37 
39 
40 #include <gsl/gsl_rng.h> // gsl random number generators
41 #include <gsl/gsl_randist.h> // gsl random number distributions
42 #include <gsl/gsl_vector.h> // gsl vector and matrix definitions
43 #include <gsl/gsl_multifit_nlin.h> // gsl multidimensional fitting
44 #include <gsl/gsl_blas.h> // gsl linear algebra stuff
45 
46 namespace OpenMS
47 {
48 
52  class OPENMS_DLLAPI LevMarqFitter1D :
53  public Fitter1D
54  {
55 
56 public:
57 
58  typedef std::vector<double> ContainerType;
59 
62  Fitter1D()
63  {
64  this->defaults_.setValue("max_iteration", 500, "Maximum number of iterations using by Levenberg-Marquardt algorithm.", StringList::create("advanced"));
65  this->defaults_.setValue("deltaAbsError", 0.0001, "Absolute error used by the Levenberg-Marquardt algorithm.", StringList::create("advanced"));
66  this->defaults_.setValue("deltaRelError", 0.0001, "Relative error used by the Levenberg-Marquardt algorithm.", StringList::create("advanced"));
67  }
68 
70  LevMarqFitter1D(const LevMarqFitter1D & source) :
71  Fitter1D(source),
72  max_iteration_(source.max_iteration_),
73  abs_error_(source.abs_error_),
74  rel_error_(source.rel_error_)
75  {
76  }
77 
79  virtual ~LevMarqFitter1D()
80  {
81  }
82 
84  virtual LevMarqFitter1D & operator=(const LevMarqFitter1D & source)
85  {
86  if (&source == this) return *this;
87 
88  Fitter1D::operator=(source);
89  max_iteration_ = source.max_iteration_;
90  abs_error_ = source.abs_error_;
91  rel_error_ = source.rel_error_;
92 
93  return *this;
94  }
95 
96 protected:
97 
105  CoordinateType abs_error_;
109 
113  virtual void printState_(Int iter, gsl_multifit_fdfsolver * s) = 0;
114 
117  {
118  return gsl_strerror(gsl_status_);
119  }
120 
126  void optimize_(const RawDataArrayType & set, Int num_params, CoordinateType x_init[],
127  Int (* residual)(const gsl_vector * x, void * params, gsl_vector * f),
128  Int (* jacobian)(const gsl_vector * x, void * params, gsl_matrix * J),
129  Int (* evaluate)(const gsl_vector * x, void * params, gsl_vector * f, gsl_matrix * J),
130  void * advanced_params
131  )
132  {
133 
134  const gsl_multifit_fdfsolver_type * T;
135  gsl_multifit_fdfsolver * s;
136 
137  Int status;
138  Int iter = 0;
139  const UInt n = (UInt)set.size();
140 
141  // number of parameters to be optimized
142  UInt p = num_params;
143 
144  // gsl always expects N>=p or default gsl error handler invoked,
145  // cause Jacobian be rectangular M x N with M>=N
146  if (n < p) throw Exception::UnableToFit(__FILE__, __LINE__, __PRETTY_FUNCTION__, "UnableToFit-FinalSet", "Skipping feature, gsl always expects N>=p");
147 
148  // allocate space for a covariance matrix of size p by p
149  gsl_matrix * covar = gsl_matrix_alloc(p, p);
150  gsl_multifit_function_fdf f;
151 
152  gsl_vector_view x = gsl_vector_view_array(x_init, p);
153 
154  gsl_rng_env_setup();
155 
156  // set up the function to be fit
157  f.f = (residual); // the function of residuals
158  f.df = (jacobian); // the gradient of this function
159  f.fdf = (evaluate); // combined function and gradient
160  f.n = set.size(); // number of points in the data set
161  f.p = p; // number of parameters in the fit function
162  f.params = advanced_params; // // structure with the data and error bars
163 
164  T = gsl_multifit_fdfsolver_lmsder;
165  s = gsl_multifit_fdfsolver_alloc(T, n, p);
166  gsl_multifit_fdfsolver_set(s, &f, &x.vector);
167 
168 #ifdef DEBUG_FEATUREFINDER
169  printState_(iter, s);
170 #endif
171 
172  // this is the loop for fitting
173  do
174  {
175  iter++;
176 
177  // perform a single iteration of the fitting routine
178  status = gsl_multifit_fdfsolver_iterate(s);
179 
180 #ifdef DEBUG_FEATUREFINDER
181  // customized routine to print out current parameters
182  printState_(iter, s);
183 #endif
184 
185  /* check if solver is stuck */
186  if (status) break;
187 
188  // test for convergence with an absolute and relative error
189  status = gsl_multifit_test_delta(s->dx, s->x, abs_error_, rel_error_);
190  }
191  while (status == GSL_CONTINUE && iter < max_iteration_);
192 
193  // This function uses Jacobian matrix J to compute the covariance matrix of the best-fit parameters, covar.
194  // The parameter epsrel (0.0) is used to remove linear-dependent columns when J is rank deficient.
195  gsl_multifit_covar(s->J, 0.0, covar);
196 
197 #ifdef DEBUG_FEATUREFINDER
198  gsl_matrix_fprintf(stdout, covar, "covar %g");
199 #endif
200 
201 #define FIT(i) gsl_vector_get(s->x, i)
202 #define ERR(i) sqrt(gsl_matrix_get(covar, i, i))
203 
204  // Set GSl status
205  gsl_status_ = status;
206 
207 #ifdef DEBUG_FEATUREFINDER
208  {
209  // chi-squared value
210  DoubleReal chi = gsl_blas_dnrm2(s->f);
211  DoubleReal dof = n - p;
212  DoubleReal c = GSL_MAX_DBL(1, chi / sqrt(dof));
213 
214  printf("chisq/dof = %g\n", pow(chi, 2.0) / dof);
215 
216  for (Size i = 0; i < p; ++i)
217  {
218  std::cout << i;
219  printf(".Parameter = %.5f +/- %.5f\n", FIT(i), c * ERR(i));
220  }
221  }
222 #endif
223 
224  // set optimized parameters
225  for (Size i = 0; i < p; ++i)
226  {
227  x_init[i] = FIT(i);
228  }
229 
230  gsl_multifit_fdfsolver_free(s);
231  gsl_matrix_free(covar);
232 
233  }
234 
236  {
238  max_iteration_ = this->param_.getValue("max_iteration");
239  abs_error_ = this->param_.getValue("deltaAbsError");
240  rel_error_ = this->param_.getValue("deltaRelError");
241  }
242 
243  };
244 }
245 
246 #endif // OPENMS_TRANSFORMATIONS_FEATUREFINDER_LEVMARQFITTER1D_H
Abstract class for 1D-model fitter using Levenberg-Marquardt algorithm for parameter optimization...
Definition: LevMarqFitter1D.h:52
A more convenient string class.
Definition: String.h:56
virtual ~LevMarqFitter1D()
destructor
Definition: LevMarqFitter1D.h:79
Int max_iteration_
Maximum number of iterations.
Definition: LevMarqFitter1D.h:103
#define FIT(i)
virtual LevMarqFitter1D & operator=(const LevMarqFitter1D &source)
assignment operator
Definition: LevMarqFitter1D.h:84
LevMarqFitter1D(const LevMarqFitter1D &source)
copy constructor
Definition: LevMarqFitter1D.h:70
unsigned int UInt
Unsigned integer type.
Definition: Types.h:92
const double c
#define ERR(i)
virtual Fitter1D & operator=(const Fitter1D &source)
assignment operator
CoordinateType rel_error_
Relative error.
Definition: LevMarqFitter1D.h:108
int evaluate(const gsl_vector *x, void *params, gsl_vector *f, gsl_matrix *J)
Driver function for the evaluation of function and jacobian.
void optimize_(const RawDataArrayType &set, Int num_params, CoordinateType 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), void *advanced_params)
Optimize start parameter.
Definition: LevMarqFitter1D.h:126
LevMarqFitter1D()
Default constructor.
Definition: LevMarqFitter1D.h:61
bool symmetric_
Parameter indicates symmetric peaks.
Definition: LevMarqFitter1D.h:101
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...
const String getGslStatus_()
Return GSL status as string.
Definition: LevMarqFitter1D.h:116
int residual(const gsl_vector *x, void *params, gsl_vector *f)
Evaluation of the target function for nonlinear optimization.
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
Abstract base class for all 1D-dimensional model fitter.
Definition: Fitter1D.h:59
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...
std::vector< double > ContainerType
Definition: LevMarqFitter1D.h:58
virtual void updateMembers_()
This method is used to update extra member variables at the end of the setParameters() method...
int Int
Signed integer type.
Definition: Types.h:100
CoordinateType abs_error_
Absolute error.
Definition: LevMarqFitter1D.h:106
Int gsl_status_
GSL status.
Definition: LevMarqFitter1D.h:99
std::vector< PeakType > RawDataArrayType
Raw data container type using for the temporary storage of the input data.
Definition: Fitter1D.h:76
void updateMembers_()
This method is used to update extra member variables at the end of the setParameters() method...
Definition: LevMarqFitter1D.h:235

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