ROOT  6.06/09
Reference Guide
GSLMultiFit.h
Go to the documentation of this file.
1 // @(#)root/mathmore:$Id$
2 // Author: L. Moneta Wed Dec 20 17:26:06 2006
3 
4 /**********************************************************************
5  * *
6  * Copyright (c) 2006 LCG ROOT Math Team, CERN/PH-SFT *
7  * *
8  * This library is free software; you can redistribute it and/or *
9  * modify it under the terms of the GNU General Public License *
10  * as published by the Free Software Foundation; either version 2 *
11  * of the License, or (at your option) any later version. *
12  * *
13  * This library 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 GNU *
16  * General Public License for more details. *
17  * *
18  * You should have received a copy of the GNU General Public License *
19  * along with this library (see file COPYING); if not, write *
20  * to the Free Software Foundation, Inc., 59 Temple Place, Suite *
21  * 330, Boston, MA 02111-1307 USA, or contact the author. *
22  * *
23  **********************************************************************/
24 
25 // Header file for class GSLMultiFit
26 
27 #ifndef ROOT_Math_GSLMultiFit
28 #define ROOT_Math_GSLMultiFit
29 
30 #include "gsl/gsl_vector.h"
31 #include "gsl/gsl_matrix.h"
32 #include "gsl/gsl_multifit_nlin.h"
33 #include "gsl/gsl_blas.h"
34 #include "gsl/gsl_version.h"
36 
37 #include "Math/IFunction.h"
38 #include <string>
39 #include <cassert>
40 
41 
42 namespace ROOT {
43 
44  namespace Math {
45 
46 
47 /**
48  GSLMultiFit, internal class for implementing GSL non linear least square GSL fitting
49 
50  @ingroup MultiMin
51 */
52 class GSLMultiFit {
53 
54 public:
55 
56  /**
57  Default constructor
58  No need to specify the type so far since only one solver exists so far
59  */
60  GSLMultiFit (const gsl_multifit_fdfsolver_type * type = 0) :
61  fSolver(0),
62  fVec(0),
63  fCov(0),
64  fType(type)
65  {
66  if (fType == 0) fType = gsl_multifit_fdfsolver_lmsder; // default value
67  }
68 
69  /**
70  Destructor (no operations)
71  */
73  if (fSolver) gsl_multifit_fdfsolver_free(fSolver);
74  if (fVec != 0) gsl_vector_free(fVec);
75  if (fCov != 0) gsl_matrix_free(fCov);
76  }
77 
78 private:
79  // usually copying is non trivial, so we make this unaccessible
80 
81  /**
82  Copy constructor
83  */
85 
86  /**
87  Assignment operator
88  */
90  if (this == &rhs) return *this; // time saving self-test
91  return *this;
92  }
93 
94 
95 public:
96 
97  /// create the minimizer from the type and size of number of fitting points and number of parameters
98  void CreateSolver(unsigned int npoints, unsigned int npar) {
99  if (fSolver) gsl_multifit_fdfsolver_free(fSolver);
100  fSolver = gsl_multifit_fdfsolver_alloc(fType, npoints, npar);
101  }
102 
103  /// set the solver parameters
104  template<class Func>
105  int Set(const std::vector<Func> & funcVec, const double * x) {
106  // create a vector of the fit contributions
107  // create function wrapper from an iterator of functions
108  unsigned int npts = funcVec.size();
109  if (npts == 0) return -1;
110 
111  unsigned int npar = funcVec[0].NDim();
112  // Remove unused typedef to remove warning in GCC48
113  // http://gcc.gnu.org/gcc-4.8/porting_to.html
114  // typedef typename std::vector<Func> FuncVec;
115  //FuncIt funcIter = funcVec.begin();
116  fFunc.SetFunction(funcVec, npts, npar);
117  // create solver object
118  CreateSolver( npts, npar );
119  // set initial values
120  if (fVec != 0) gsl_vector_free(fVec);
121  fVec = gsl_vector_alloc( npar );
122  std::copy(x,x+npar, fVec->data);
123  assert(fSolver != 0);
124  return gsl_multifit_fdfsolver_set(fSolver, fFunc.GetFunc(), fVec);
125  }
126 
127  std::string Name() const {
128  if (fSolver == 0) return "undefined";
129  return std::string(gsl_multifit_fdfsolver_name(fSolver) );
130  }
131 
132  int Iterate() {
133  if (fSolver == 0) return -1;
134  return gsl_multifit_fdfsolver_iterate(fSolver);
135  }
136 
137  /// parameter values at the minimum
138  const double * X() const {
139  if (fSolver == 0) return 0;
140  gsl_vector * x = gsl_multifit_fdfsolver_position(fSolver);
141  return x->data;
142  }
143 
144  /// gradient value at the minimum
145  const double * Gradient() const {
146  if (fSolver == 0) return 0;
147 #if GSL_MAJOR_VERSION > 1
148  fType->gradient(fSolver->state, fVec);
149 #else
150  gsl_multifit_gradient(fSolver->J, fSolver->f,fVec);
151 #endif
152  return fVec->data;
153  }
154 
155  /// return covariance matrix of the parameters
156  const double * CovarMatrix() const {
157  if (fSolver == 0) return 0;
158  if (fCov != 0) gsl_matrix_free(fCov);
159  unsigned int npar = fSolver->fdf->p;
160  fCov = gsl_matrix_alloc( npar, npar );
161  static double kEpsrel = 0.0001;
162 #if GSL_MAJOR_VERSION > 1
163  gsl_matrix* J = gsl_matrix_alloc(npar,npar);
164  gsl_multifit_fdfsolver_jac (fSolver, J);
165  int ret = gsl_multifit_covar(J, kEpsrel, fCov);
166  gsl_matrix_free(J);
167 #else
168  int ret = gsl_multifit_covar(fSolver->J, kEpsrel, fCov);
169 #endif
170  if (ret != GSL_SUCCESS) return 0;
171  return fCov->data;
172  }
173 
174  /// test gradient (ask from solver gradient vector)
175  int TestGradient(double absTol) const {
176  if (fSolver == 0) return -1;
177  Gradient();
178  return gsl_multifit_test_gradient( fVec, absTol);
179  }
180 
181  /// test using abs and relative tolerance
182  /// |dx| < absTol + relTol*|x| for every component
183  int TestDelta(double absTol, double relTol) const {
184  if (fSolver == 0) return -1;
185  return gsl_multifit_test_delta(fSolver->dx, fSolver->x, absTol, relTol);
186  }
187 
188  // calculate edm 1/2 * ( grad^T * Cov * grad )
189  double Edm() const {
190  // product C * g
191  double edm = -1;
192  const double * g = Gradient();
193  if (g == 0) return edm;
194  const double * c = CovarMatrix();
195  if (c == 0) return edm;
196  gsl_vector * tmp = gsl_vector_alloc( fSolver->fdf->p );
197  int status = gsl_blas_dgemv(CblasNoTrans, 1.0, fCov, fVec, 0.,tmp);
198  if (status == 0) status |= gsl_blas_ddot(fVec, tmp, &edm);
199  gsl_vector_free(tmp);
200  if (status != 0) return -1;
201  // need to divide by 2 ??
202  return 0.5*edm;
203 
204  }
205 
206 
207 private:
208 
210  gsl_multifit_fdfsolver * fSolver;
211  // cached vector to avoid re-allocating every time a new one
212  mutable gsl_vector * fVec;
213  mutable gsl_matrix * fCov;
214  const gsl_multifit_fdfsolver_type * fType;
215 
216 };
217 
218  } // end namespace Math
219 
220 } // end namespace ROOT
221 
222 
223 #endif /* ROOT_Math_GSLMultiFit */
const gsl_multifit_fdfsolver_type * fType
Definition: GSLMultiFit.h:214
const double absTol
int TestDelta(double absTol, double relTol) const
test using abs and relative tolerance |dx| < absTol + relTol*|x| for every component ...
Definition: GSLMultiFit.h:183
Namespace for new ROOT classes and functions.
Definition: ROOT.py:1
GSLMultiFit, internal class for implementing GSL non linear least square GSL fitting.
Definition: GSLMultiFit.h:52
GSLMultiFit(const gsl_multifit_fdfsolver_type *type=0)
Default constructor No need to specify the type so far since only one solver exists so far...
Definition: GSLMultiFit.h:60
const double * Gradient() const
gradient value at the minimum
Definition: GSLMultiFit.h:145
#define assert(cond)
Definition: unittest.h:542
TAlienJobStatus * status
Definition: TAlienJob.cxx:51
int Set(const std::vector< Func > &funcVec, const double *x)
set the solver parameters
Definition: GSLMultiFit.h:105
gsl_multifit_fdfsolver * fSolver
Definition: GSLMultiFit.h:210
std::string Name() const
Definition: GSLMultiFit.h:127
const double * CovarMatrix() const
return covariance matrix of the parameters
Definition: GSLMultiFit.h:156
GSLMultiFit(const GSLMultiFit &)
Copy constructor.
Definition: GSLMultiFit.h:84
Double_t x[n]
Definition: legend1.C:17
const double * X() const
parameter values at the minimum
Definition: GSLMultiFit.h:138
GSLMultiFit & operator=(const GSLMultiFit &rhs)
Assignment operator.
Definition: GSLMultiFit.h:89
double Edm() const
Definition: GSLMultiFit.h:189
int TestGradient(double absTol) const
test gradient (ask from solver gradient vector)
Definition: GSLMultiFit.h:175
GSLMultiFitFunctionWrapper fFunc
Definition: GSLMultiFit.h:209
void SetFunction(const FuncVector &f, unsigned int nres, unsigned int npar)
Fill gsl function structure from a C++ function iterator and size and number of residuals.
int type
Definition: TGX11.cxx:120
Namespace for new Math classes and functions.
wrapper to a multi-dim function withtout derivatives for multi-dimensional minimization algorithm ...
void CreateSolver(unsigned int npoints, unsigned int npar)
create the minimizer from the type and size of number of fitting points and number of parameters ...
Definition: GSLMultiFit.h:98
~GSLMultiFit()
Destructor (no operations)
Definition: GSLMultiFit.h:72