ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
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"
35 
36 #include "Math/IFunction.h"
37 #include <string>
38 #include <cassert>
39 
40 
41 namespace ROOT {
42 
43  namespace Math {
44 
45 
46 /**
47  GSLMultiFit, internal class for implementing GSL non linear least square GSL fitting
48 
49  @ingroup MultiMin
50 */
51 class GSLMultiFit {
52 
53 public:
54 
55  /**
56  Default constructor
57  No need to specify the type so far since only one solver exists so far
58  */
59  GSLMultiFit (const gsl_multifit_fdfsolver_type * type = 0) :
60  fSolver(0),
61  fVec(0),
62  fCov(0),
63  fType(type)
64  {
65  if (fType == 0) fType = gsl_multifit_fdfsolver_lmsder; // default value
66  }
67 
68  /**
69  Destructor (no operations)
70  */
72  if (fSolver) gsl_multifit_fdfsolver_free(fSolver);
73  if (fVec != 0) gsl_vector_free(fVec);
74  if (fCov != 0) gsl_matrix_free(fCov);
75  }
76 
77 private:
78  // usually copying is non trivial, so we make this unaccessible
79 
80  /**
81  Copy constructor
82  */
84 
85  /**
86  Assignment operator
87  */
89  if (this == &rhs) return *this; // time saving self-test
90  return *this;
91  }
92 
93 
94 public:
95 
96  /// create the minimizer from the type and size of number of fitting points and number of parameters
97  void CreateSolver(unsigned int npoints, unsigned int npar) {
98  if (fSolver) gsl_multifit_fdfsolver_free(fSolver);
99  fSolver = gsl_multifit_fdfsolver_alloc(fType, npoints, npar);
100  }
101 
102  /// set the solver parameters
103  template<class Func>
104  int Set(const std::vector<Func> & funcVec, const double * x) {
105  // create a vector of the fit contributions
106  // create function wrapper from an iterator of functions
107  unsigned int npts = funcVec.size();
108  if (npts == 0) return -1;
109 
110  unsigned int npar = funcVec[0].NDim();
111  // Remove unused typedef to remove warning in GCC48
112  // http://gcc.gnu.org/gcc-4.8/porting_to.html
113  // typedef typename std::vector<Func> FuncVec;
114  //FuncIt funcIter = funcVec.begin();
115  fFunc.SetFunction(funcVec, npts, npar);
116  // create solver object
117  CreateSolver( npts, npar );
118  // set initial values
119  if (fVec != 0) gsl_vector_free(fVec);
120  fVec = gsl_vector_alloc( npar );
121  std::copy(x,x+npar, fVec->data);
122  assert(fSolver != 0);
123  return gsl_multifit_fdfsolver_set(fSolver, fFunc.GetFunc(), fVec);
124  }
125 
126  std::string Name() const {
127  if (fSolver == 0) return "undefined";
128  return std::string(gsl_multifit_fdfsolver_name(fSolver) );
129  }
130 
131  int Iterate() {
132  if (fSolver == 0) return -1;
133  return gsl_multifit_fdfsolver_iterate(fSolver);
134  }
135 
136  /// parameter values at the minimum
137  const double * X() const {
138  if (fSolver == 0) return 0;
139  gsl_vector * x = gsl_multifit_fdfsolver_position(fSolver);
140  return x->data;
141  }
142 
143  /// gradient value at the minimum
144  const double * Gradient() const {
145  if (fSolver == 0) return 0;
146  gsl_multifit_gradient(fSolver->J, fSolver->f,fVec);
147  return fVec->data;
148  }
149 
150  /// return covariance matrix of the parameters
151  const double * CovarMatrix() const {
152  if (fSolver == 0) return 0;
153  if (fCov != 0) gsl_matrix_free(fCov);
154  unsigned int npar = fSolver->fdf->p;
155  fCov = gsl_matrix_alloc( npar, npar );
156  static double kEpsrel = 0.0001;
157  int ret = gsl_multifit_covar(fSolver->J, kEpsrel, fCov);
158  if (ret != GSL_SUCCESS) return 0;
159  return fCov->data;
160  }
161 
162  /// test gradient (ask from solver gradient vector)
163  int TestGradient(double absTol) const {
164  if (fSolver == 0) return -1;
165  Gradient();
166  return gsl_multifit_test_gradient( fVec, absTol);
167  }
168 
169  /// test using abs and relative tolerance
170  /// |dx| < absTol + relTol*|x| for every component
171  int TestDelta(double absTol, double relTol) const {
172  if (fSolver == 0) return -1;
173  return gsl_multifit_test_delta(fSolver->dx, fSolver->x, absTol, relTol);
174  }
175 
176  // calculate edm 1/2 * ( grad^T * Cov * grad )
177  double Edm() const {
178  // product C * g
179  double edm = -1;
180  const double * g = Gradient();
181  if (g == 0) return edm;
182  const double * c = CovarMatrix();
183  if (c == 0) return edm;
184  gsl_vector * tmp = gsl_vector_alloc( fSolver->fdf->p );
185  int status = gsl_blas_dgemv(CblasNoTrans, 1.0, fCov, fVec, 0.,tmp);
186  if (status == 0) status |= gsl_blas_ddot(fVec, tmp, &edm);
187  gsl_vector_free(tmp);
188  if (status != 0) return -1;
189  // need to divide by 2 ??
190  return 0.5*edm;
191 
192  }
193 
194 
195 private:
196 
198  gsl_multifit_fdfsolver * fSolver;
199  // cached vector to avoid re-allocating every time a new one
200  mutable gsl_vector * fVec;
201  mutable gsl_matrix * fCov;
202  const gsl_multifit_fdfsolver_type * fType;
203 
204 };
205 
206  } // end namespace Math
207 
208 } // end namespace ROOT
209 
210 
211 #endif /* ROOT_Math_GSLMultiFit */
const gsl_multifit_fdfsolver_type * fType
Definition: GSLMultiFit.h:202
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:171
GSLMultiFit, internal class for implementing GSL non linear least square GSL fitting.
Definition: GSLMultiFit.h:51
return c
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:59
const double * Gradient() const
gradient value at the minimum
Definition: GSLMultiFit.h:144
#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:104
gsl_multifit_fdfsolver * fSolver
Definition: GSLMultiFit.h:198
std::string Name() const
Definition: GSLMultiFit.h:126
const double * CovarMatrix() const
return covariance matrix of the parameters
Definition: GSLMultiFit.h:151
GSLMultiFit(const GSLMultiFit &)
Copy constructor.
Definition: GSLMultiFit.h:83
Double_t x[n]
Definition: legend1.C:17
const double * X() const
parameter values at the minimum
Definition: GSLMultiFit.h:137
GSLMultiFit & operator=(const GSLMultiFit &rhs)
Assignment operator.
Definition: GSLMultiFit.h:88
double Edm() const
Definition: GSLMultiFit.h:177
int TestGradient(double absTol) const
test gradient (ask from solver gradient vector)
Definition: GSLMultiFit.h:163
GSLMultiFitFunctionWrapper fFunc
Definition: GSLMultiFit.h:197
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
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:97
~GSLMultiFit()
Destructor (no operations)
Definition: GSLMultiFit.h:71