All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
fit_base.h
Go to the documentation of this file.
1 /*
2  -------------------------------------------------------------------
3 
4  Copyright (C) 2006-2014, Andrew W. Steiner
5 
6  This file is part of O2scl.
7 
8  O2scl is free software; you can redistribute it and/or modify
9  it under the terms of the GNU General Public License as published by
10  the Free Software Foundation; either version 3 of the License, or
11  (at your option) any later version.
12 
13  O2scl is distributed in the hope that it will be useful,
14  but WITHOUT ANY WARRANTY; without even the implied warranty of
15  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16  GNU General Public License for more details.
17 
18  You should have received a copy of the GNU General Public License
19  along with O2scl. If not, see <http://www.gnu.org/licenses/>.
20 
21  -------------------------------------------------------------------
22 */
23 #ifndef O2SCL_FIT_BASE_H
24 #define O2SCL_FIT_BASE_H
25 
26 /** \file fit_base.h
27  \brief File defining \ref o2scl::fit_base and fitting functions
28 */
29 
30 #include <string>
31 
32 #include <o2scl/jacobian.h>
33 #include <o2scl/fparser.h>
34 #include <o2scl/mm_funct.h>
35 
36 #ifndef DOXYGEN_NO_O2NS
37 namespace o2scl {
38 #endif
39 
40  /// Array of multi-dimensional functions typedef (C++11 version)
41  typedef std::function<
42  double(size_t,const boost::numeric::ublas::vector<double> &,
43  double)> fit_funct11;
44 
45  /** \brief String fitting function
46 
47  Default template arguments
48  - \c vec_t - \ref boost::numeric::ublas::vector < double >
49  */
50  template<class vec_t=boost::numeric::ublas::vector<double> >
52 
53  public:
54 
55  /** \brief Specify a fitting function through a string
56  */
57  fit_funct11_strings(std::string formula, std::string parms,
58  std::string var, int nauxp=0,
59  std::string auxp="") {
60  if(nauxp<1) {
61  std::string both=parms+","+var;
62  fpw.Parse(formula,both);
63  naux=0;
64  st_auxp="";
65  } else {
66  std::string all=parms+","+var+","+auxp;
67  fpw.Parse(formula,all);
68  naux=nauxp;
69  aux=new double[naux];
70  st_auxp=auxp;
71  }
72  st_form=formula;
73  st_parms=parms;
74  st_var=var;
75  }
76 
77  virtual ~fit_funct11_strings() {
78  if (naux>0) {
79  delete[] aux;
80  }
81  }
82 
83  /** \brief Set the values of the auxilliary parameters that were
84  specified in \c auxp in the constructor
85  */
86  int set_aux_parms(vec_t &ap) {
87  for(int i=0;i<naux;i++) {
88  aux[i]=ap[i];
89  }
90  return 0;
91  }
92 
93  /** \brief Using parameters in \c p, predict \c y given \c x
94  */
95  virtual double operator()(size_t np, const vec_t &p, double x) {
96 
97  size_t i;
98  double y;
99  if(naux<1) {
100  double *both=new double[np+1];
101  for(i=0;i<np;i++) both[i]=p[i];
102  both[np]=x;
103  y=fpw.Eval(both);
104  delete[] both;
105  } else {
106  double *all=new double[np+1+naux];
107  for(i=0;i<np;i++) all[i]=p[i];
108  all[np]=x;
109  for(i=np+1;i<np+1+naux;i++) all[i]=aux[i-np-1];
110  y=fpw.Eval(all);
111  delete[] all;
112  }
113  return y;
114  }
115 
116 #ifndef DOXYGEN_INTERNAL
117 
118  protected:
119 
120  /// The function evaluation object
122 
123  /// \name The auxillary parameters
124  //@{
125  int naux;
126  double *aux;
127  std::string st_auxp;
128  //@}
129 
130  /// The formula
131  std::string st_form;
132  /// The parameters
133  std::string st_parms;
134  /// The variables
135  std::string st_var;
136 
137  fit_funct11_strings() {};
138 
139  /// Specify the strings which define the fitting function
140  int set_function(std::string formula, std::string parms,
141  std::string var, int nauxp=0, std::string auxp="");
142 
143  private:
144 
146  fit_funct11_strings& operator=(const fit_funct11_strings&);
147 
148 #endif
149 
150  };
151 
152  /** \brief Generalized fitting function [abstract base]
153 
154  Default template arguments
155  - \c vec_t - \ref boost::numeric::ublas::vector < double >
156  - \c mat_t - \ref boost::numeric::ublas::matrix < double >
157  */
158  template<class vec_t=boost::numeric::ublas::vector<double>,
159  class mat_t=boost::numeric::ublas::matrix<double> >
161 
162  public:
163 
164  gen_fit_funct() {}
165 
166  virtual ~gen_fit_funct() {}
167 
168  /** \brief Using parameters in \c p, compute the
169  relative deviations in \c f
170  */
171  virtual void operator()(size_t np, const vec_t &p, size_t nd,
172  vec_t &f)=0;
173 
174  /** \brief Using parameters in \c p, compute the Jacobian
175  in \c J
176  */
177  virtual void jac(size_t np, vec_t &p, size_t nd, vec_t &f,
178  mat_t &J)=0;
179 
180  /// Return the number of data points
181  virtual size_t get_ndata()=0;
182 
183 #ifndef DOXYGEN_INTERNAL
184 
185  private:
186 
187  gen_fit_funct(const gen_fit_funct &);
188  gen_fit_funct& operator=(const gen_fit_funct&);
189 
190 #endif
191 
192  };
193 
194  /** \brief Standard fitting function based on one-dimensional
195  data with a numerical Jacobian
196 
197  This class specifies the deviations (in <tt>operator()</tt> ) and
198  Jacobian (in \ref jac()) for a fitting class like \ref fit_nonlin.
199  It assumes a one-dimensional data set with no uncertainty in the
200  abcissae and a fitting function specified in a form similar to
201  \ref fit_funct11.
202  \comment
203  For some reason the reference to operator() above doesn't work
204  in doxygen.
205  \endcomment
206 
207  The default method for numerically computing the Jacobian is
208  from \ref jacobian_gsl. This default is identical to the GSL
209  approach, except that the default value of \ref
210  jacobian_gsl::epsmin is non-zero. See \ref jacobian_gsl for more
211  details.
212 
213  Default template arguments
214  - \c vec_t - \ref boost::numeric::ublas::vector < double >
215  - \c mat_t - \ref boost::numeric::ublas::matrix < double >
216  - \c func_t - \ref fit_funct11
217 
218  \future Allow a user-specified Jacobian or make that into
219  a separate class?
220  \future Default constructor?
221  */
222  template<class vec_t=boost::numeric::ublas::vector<double>,
223  class mat_t=boost::numeric::ublas::matrix<double>,
224  class fit_func_t=fit_funct11> class chi_fit_funct :
225  public gen_fit_funct<vec_t,mat_t> {
226 
227  public:
228 
229  /** \brief Create an object with specified data and specified
230  fitting function
231  */
232  chi_fit_funct(size_t ndat, const vec_t &xdat, const vec_t &ydat,
233  const vec_t &yerr, fit_func_t &fun) {
234  ndat_=ndat;
235  xdat_=&xdat;
236  ydat_=&ydat;
237  yerr_=&yerr;
238  fun_=&fun;
239 
240  mfm=std::bind(std::mem_fn<int(size_t,const vec_t &,vec_t &)>
242  std::placeholders::_1,std::placeholders::_2,
243  std::placeholders::_3);
244 
245  //mfm.set_function(this,&chi_fit_funct::jac_mm_funct);
247  double sqrt_dbl_eps=sqrt(std::numeric_limits<double>::epsilon());
248  auto_jac.epsrel=sqrt_dbl_eps;
249  }
250 
251  /** \brief Set the data to be fit
252  */
253  void set_data(size_t ndat, const vec_t &xdat, const vec_t &ydat,
254  const vec_t &yerr) {
255  ndat_=ndat;
256  xdat_=&xdat;
257  ydat_=&ydat;
258  yerr_=&yerr;
259  return;
260  }
261 
262  /** \brief Set the fitting function
263  */
264  void set_func(fit_func_t &fun) {
265  fun_=&fun;
266  return;
267  }
268 
269  /** \brief Return \f$ \chi^2 \f$
270  */
271  virtual double chi2(size_t np, const vec_t &p) {
272  double ret=0.0;
273  for(size_t i=0;i<ndat_;i++) {
274  double yi=((*fun_)(np,p,(*xdat_)[i])-(*ydat_)[i])/((*yerr_)[i]);
275  ret+=yi*yi;
276  }
277  return ret;
278  }
279 
280  /** \brief Using parameters in \c p, compute the
281  relative deviations in \c f
282  */
283  virtual void operator()(size_t np, const vec_t &p, size_t nd, vec_t &f) {
284 
285  for(size_t i=0;i<nd;i++) {
286  double yi=(*fun_)(np,p,(*xdat_)[i]);
287  f[i]=(yi-(*ydat_)[i])/((*yerr_)[i]);
288  }
289  return;
290  }
291 
292  /** \brief Using parameters in \c p, compute the Jacobian
293  in \c J
294  */
295  virtual void jac(size_t np, vec_t &p, size_t nd, vec_t &f,
296  mat_t &J) {
297 
298  auto_jac(np,p,nd,f,J);
299 
300  return;
301  }
302 
303  /// Return the number of data points
304  virtual size_t get_ndata() {
305  return ndat_;
306  }
307 
308  /// Automatic Jacobian object
310  vec_t,mat_t> auto_jac;
311 
312 #ifndef DOXYGEN_INTERNAL
313 
314  protected:
315 
316  /// Reformulate <tt>operator()</tt> into a \ref mm_funct11 object
317  int jac_mm_funct(size_t np, const vec_t &p, vec_t &f) {
318  operator()(np,p,ndat_,f);
319  return 0;
320  }
321 
322  /// Function object for Jacobian object
323  std::function<int(size_t,const vec_t &,vec_t &)> mfm;
324 
325  /// \name Data and uncertainties
326  //@{
327  size_t ndat_;
328  const vec_t *xdat_;
329  const vec_t *ydat_;
330  const vec_t *yerr_;
331  //@}
332 
333  /// Fitting function
334  fit_func_t *fun_;
335 
336  private:
337 
338  chi_fit_funct(const chi_fit_funct &);
339  chi_fit_funct& operator=(const chi_fit_funct&);
340 
341 #endif
342 
343  };
344 
345  // ----------------------------------------------------------------
346  // Fitter base
347  // ----------------------------------------------------------------
348 
349  /** \brief Non-linear least-squares fitting [abstract base]
350  */
351  template<class func_t=fit_funct11,
354 
355  public:
356 
357  fit_base() {
358  ntrial=500;
359  tol_abs=1.0e-4;
360  tol_rel=1.0e-4;
361  verbose=0;
362  }
363 
364  virtual ~fit_base() {}
365 
366  /// Maximum number of iterations (default 500)
367  size_t ntrial;
368 
369  /// Absolute tolerance (default 1.0e-4)
370  double tol_abs;
371 
372  /// (default 1.0e-4)
373  double tol_rel;
374 
375  /** \brief Print out iteration information.
376 
377  Depending on the value of the variable verbose, this prints out
378  the iteration information. If verbose=0, then no information is
379  printed, while if verbose>1, then after each iteration, the
380  present values of x and y are output to std::cout along with the
381  iteration number. If verbose>=2 then each iteration waits for a
382  character.
383  */
384  virtual int print_iter(size_t nv, vec_t &x, double y, int iter,
385  double value=0.0, double limit=0.0) {
386  if (verbose<=0) return 0;
387 
388  size_t i;
389  char ch;
390 
391  std::cout << "Iteration: " << iter << std::endl;
392  std::cout << "x: ";
393  for(i=0;i<nv;i++) std::cout << x[i] << " ";
394  std::cout << std::endl;
395  std::cout << "y: " << y << " Val: " << value << " Lim: " << limit
396  << std::endl;
397  if (verbose>1) {
398  std::cout << "Press a key and type enter to continue. ";
399  std::cin >> ch;
400  }
401 
402  return 0;
403  }
404 
405  /** \brief Fit function \c fitfun using parameters in \c parms as
406  initial guesses
407 
408  The covariance matrix for the parameters is returned in \c covar
409  and the value of \f$ \chi^2 \f$ is returned in \c chi2.
410 
411  */
412  virtual int fit(size_t npar, vec_t &parms, mat_t &covar,
413  double &chi2, func_t &fitfun)=0;
414 
415  /** \brief An integer describing the verbosity of the output
416  */
417  int verbose;
418 
419  /// Return string denoting type ("fit_base")
420  virtual const char *type() { return "fit_base"; }
421 
422  /// The number of data points
423  size_t n_dat;
424 
425  /// The number of parameters
426  size_t n_par;
427 
428 
429  };
430 
431 #ifndef DOXYGEN_NO_O2NS
432 }
433 #endif
434 
435 #endif
Parse a mathematical function specified in a string.
Definition: fparser.h:29
Generalized fitting function [abstract base].
Definition: fit_base.h:160
virtual int set_function(func_t &f)
Set the function to compute the Jacobian of.
Definition: jacobian.h:74
virtual void jac(size_t np, vec_t &p, size_t nd, vec_t &f, mat_t &J)
Using parameters in p, compute the Jacobian in J.
Definition: fit_base.h:295
String fitting function.
Definition: fit_base.h:51
void set_func(fit_func_t &fun)
Set the fitting function.
Definition: fit_base.h:264
Standard fitting function based on one-dimensional data with a numerical Jacobian.
Definition: fit_base.h:224
int set_aux_parms(vec_t &ap)
Set the values of the auxilliary parameters that were specified in auxp in the constructor.
Definition: fit_base.h:86
virtual void jac(size_t np, vec_t &p, size_t nd, vec_t &f, mat_t &J)=0
Using parameters in p, compute the Jacobian in J.
virtual const char * type()
Return string denoting type ("fit_base")
Definition: fit_base.h:420
size_t n_dat
The number of data points.
Definition: fit_base.h:423
void set_data(size_t ndat, const vec_t &xdat, const vec_t &ydat, const vec_t &yerr)
Set the data to be fit.
Definition: fit_base.h:253
std::function< double(size_t, const boost::numeric::ublas::vector< double > &, double)> fit_funct11
Array of multi-dimensional functions typedef (C++11 version)
Definition: fit_base.h:43
double epsrel
The relative stepsize for finite-differencing (default )
Definition: jacobian.h:162
size_t ntrial
Maximum number of iterations (default 500)
Definition: fit_base.h:367
virtual void operator()(size_t np, const vec_t &p, size_t nd, vec_t &f)=0
Using parameters in p, compute the relative deviations in f.
int set_function(std::string formula, std::string parms, std::string var, int nauxp=0, std::string auxp="")
Specify the strings which define the fitting function.
fit_func_t * fun_
Fitting function.
Definition: fit_base.h:334
std::function< int(size_t, const vec_t &, vec_t &)> mfm
Function object for Jacobian object.
Definition: fit_base.h:323
double tol_rel
(default 1.0e-4)
Definition: fit_base.h:373
virtual int print_iter(size_t nv, vec_t &x, double y, int iter, double value=0.0, double limit=0.0)
Print out iteration information.
Definition: fit_base.h:384
int verbose
An integer describing the verbosity of the output.
Definition: fit_base.h:417
std::string st_form
The formula.
Definition: fit_base.h:131
std::string st_parms
The parameters.
Definition: fit_base.h:133
size_t n_par
The number of parameters.
Definition: fit_base.h:426
std::string st_var
The variables.
Definition: fit_base.h:135
FunctionParser fpw
The function evaluation object.
Definition: fit_base.h:121
virtual size_t get_ndata()
Return the number of data points.
Definition: fit_base.h:304
virtual size_t get_ndata()=0
Return the number of data points.
virtual void operator()(size_t np, const vec_t &p, size_t nd, vec_t &f)
Using parameters in p, compute the relative deviations in f.
Definition: fit_base.h:283
Non-linear least-squares fitting [abstract base].
Definition: fit_base.h:353
virtual double chi2(size_t np, const vec_t &p)
Return .
Definition: fit_base.h:271
chi_fit_funct(size_t ndat, const vec_t &xdat, const vec_t &ydat, const vec_t &yerr, fit_func_t &fun)
Create an object with specified data and specified fitting function.
Definition: fit_base.h:232
jacobian_gsl< std::function< int(size_t, const vec_t &, vec_t &)>, vec_t, mat_t > auto_jac
Automatic Jacobian object.
Definition: fit_base.h:310
double tol_abs
Absolute tolerance (default 1.0e-4)
Definition: fit_base.h:370
int jac_mm_funct(size_t np, const vec_t &p, vec_t &f)
Reformulate operator() into a mm_funct11 object.
Definition: fit_base.h:317
fit_funct11_strings(std::string formula, std::string parms, std::string var, int nauxp=0, std::string auxp="")
Specify a fitting function through a string.
Definition: fit_base.h:57
Simple automatic Jacobian.
Definition: jacobian.h:126
virtual int fit(size_t npar, vec_t &parms, mat_t &covar, double &chi2, func_t &fitfun)=0
Fit function fitfun using parameters in parms as initial guesses.
virtual double operator()(size_t np, const vec_t &p, double x)
Using parameters in p, predict y given x.
Definition: fit_base.h:95

Documentation generated with Doxygen. Provided under the GNU Free Documentation License (see License Information).
Hosted at Get Object-oriented Scientific Computing
Lib at SourceForge.net. Fast, secure and Free Open Source software
downloads..