Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
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 <vector>
40#include <cassert>
41
42
43namespace ROOT {
44
45 namespace Math {
46
47
48/**
49 GSLMultiFit, internal class for implementing GSL non linear least square GSL fitting
50
51 @ingroup MultiMin
52*/
54
55public:
56
57 /**
58 Default constructor
59 No need to specify the type so far since only one solver exists so far
60 */
61 GSLMultiFit (const gsl_multifit_fdfsolver_type * type = nullptr) :
62 fSolver(nullptr),
63 fVec(nullptr),
64 fTmp(nullptr),
65 fCov(nullptr),
66#if GSL_MAJOR_VERSION > 1
67 fJac(nullptr),
68#endif
69 fType(type)
70 {
71 if (fType == nullptr) fType = gsl_multifit_fdfsolver_lmsder; // default value
72 }
73
74 /**
75 Destructor (no operations)
76 */
78 if (fSolver) gsl_multifit_fdfsolver_free(fSolver);
79 if (fVec != nullptr) gsl_vector_free(fVec);
80 if (fTmp != nullptr) gsl_vector_free(fTmp);
81 if (fCov != nullptr) gsl_matrix_free(fCov);
82#if GSL_MAJOR_VERSION > 1
83 if (fJac != nullptr) gsl_matrix_free(fJac);
84#endif
85 }
86
87 // usually copying is non trivial, so we delete this
88 GSLMultiFit(const GSLMultiFit &) = delete;
89 GSLMultiFit & operator = (const GSLMultiFit & rhs) = delete;
92
93 /// create the minimizer from the type and size of number of fitting points and number of parameters
94 void CreateSolver(unsigned int npoints, unsigned int npar) {
95 if (fSolver) gsl_multifit_fdfsolver_free(fSolver);
96 fSolver = gsl_multifit_fdfsolver_alloc(fType, npoints, npar);
97 if (fVec != nullptr) gsl_vector_free(fVec);
98 fVec = gsl_vector_alloc( npar );
99 if (fTmp != nullptr) gsl_vector_free(fTmp);
100 fTmp = gsl_vector_alloc( npar );
101 if (fCov != nullptr) gsl_matrix_free(fCov);
102 fCov = gsl_matrix_alloc( npar, npar );
103#if GSL_MAJOR_VERSION > 1
104 if (fJac != nullptr) gsl_matrix_free(fJac);
105 fJac = gsl_matrix_alloc( npoints, npar );
106#endif
107 }
108
109 /// set the solver parameters
110 template<class Func>
111 int Set(const std::vector<Func> & funcVec, const double * x) {
112 // create a vector of the fit contributions
113 // create function wrapper from an iterator of functions
114 unsigned int npts = funcVec.size();
115 if (npts == 0) return -1;
116
117 unsigned int npar = funcVec[0].NDim();
118 // Remove unused typedef to remove warning in GCC48
119 // http://gcc.gnu.org/gcc-4.8/porting_to.html
120 // typedef typename std::vector<Func> FuncVec;
121 //FuncIt funcIter = funcVec.begin();
122 fFunc.SetFunction(funcVec, npts, npar);
123 // create solver object
124 CreateSolver( npts, npar );
125 assert(fSolver != nullptr);
126 // set initial values
127 assert(fVec != nullptr);
128 std::copy(x,x+npar, fVec->data);
129 assert(fTmp != nullptr);
130 gsl_vector_set_zero(fTmp);
131 assert(fCov != nullptr);
132 gsl_matrix_set_zero(fCov);
133#if GSL_MAJOR_VERSION > 1
134 assert(fJac != nullptr);
135 gsl_matrix_set_zero(fJac);
136#endif
137 return gsl_multifit_fdfsolver_set(fSolver, fFunc.GetFunc(), fVec);
138 }
139
140 std::string Name() const {
141 if (fSolver == nullptr) return "undefined";
142 return std::string(gsl_multifit_fdfsolver_name(fSolver) );
143 }
144
145 int Iterate() {
146 if (fSolver == nullptr) return -1;
147 return gsl_multifit_fdfsolver_iterate(fSolver);
148 }
149
150 /// parameter values at the minimum
151 const double * X() const {
152 if (fSolver == nullptr) return nullptr;
153 gsl_vector * x = gsl_multifit_fdfsolver_position(fSolver);
154 return x->data;
155 }
156
157 /// gradient value at the minimum
158 const double * Gradient() const {
159 if (fSolver == nullptr) return nullptr;
160#if GSL_MAJOR_VERSION > 1
161 fType->gradient(fSolver->state, fVec);
162#else
163 gsl_multifit_gradient(fSolver->J, fSolver->f,fVec);
164#endif
165 return fVec->data;
166 }
167
168 /// return covariance matrix of the parameters
169 const double * CovarMatrix() const {
170 if (fSolver == nullptr) return nullptr;
171 static double kEpsrel = 0.0001;
172#if GSL_MAJOR_VERSION > 1
173 gsl_multifit_fdfsolver_jac (fSolver, fJac);
174 int ret = gsl_multifit_covar(fJac, kEpsrel, fCov);
175#else
176 int ret = gsl_multifit_covar(fSolver->J, kEpsrel, fCov);
177#endif
178 if (ret != GSL_SUCCESS) return nullptr;
179 return fCov->data;
180 }
181
182 /// test gradient (ask from solver gradient vector)
183 int TestGradient(double absTol) const {
184 if (fSolver == nullptr) return -1;
185 Gradient();
186 return gsl_multifit_test_gradient( fVec, absTol);
187 }
188
189 /// test using abs and relative tolerance
190 /// |dx| < absTol + relTol*|x| for every component
191 int TestDelta(double absTol, double relTol) const {
192 if (fSolver == nullptr) return -1;
193 return gsl_multifit_test_delta(fSolver->dx, fSolver->x, absTol, relTol);
194 }
195
196 // calculate edm 1/2 * ( grad^T * Cov * grad )
197 double Edm() const {
198 // product C * g
199 double edm = -1;
200 const double * g = Gradient();
201 if (g == nullptr) return edm;
202 const double * c = CovarMatrix();
203 if (c == nullptr) return edm;
204 if (fTmp == nullptr) return edm;
205 int status = gsl_blas_dgemv(CblasNoTrans, 1.0, fCov, fVec, 0.,fTmp);
206 if (status == 0) status |= gsl_blas_ddot(fVec, fTmp, &edm);
207 if (status != 0) return -1;
208 // need to divide by 2 ??
209 return 0.5*edm;
210
211 }
212
213
214private:
215
217 gsl_multifit_fdfsolver * fSolver;
218 // cached vector to avoid re-allocating every time a new one
219 mutable gsl_vector * fVec;
220 mutable gsl_vector * fTmp;
221 mutable gsl_matrix * fCov;
222#if GSL_MAJOR_VERSION > 1
223 mutable gsl_matrix * fJac;
224#endif
225 const gsl_multifit_fdfsolver_type * fType;
226
227};
228
229 } // end namespace Math
230
231} // end namespace ROOT
232
233
234#endif /* ROOT_Math_GSLMultiFit */
#define c(i)
Definition RSha256.hxx:101
#define g(i)
Definition RSha256.hxx:105
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h Atom_t Int_t ULong_t ULong_t unsigned char prop_list Atom_t Atom_t Atom_t Time_t type
wrapper to a multi-dim function withtout derivatives for multi-dimensional minimization algorithm
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.
GSLMultiFit, internal class for implementing GSL non linear least square GSL fitting.
Definition GSLMultiFit.h:53
GSLMultiFit(GSLMultiFit &&)=delete
int TestGradient(double absTol) const
test gradient (ask from solver gradient vector)
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:94
GSLMultiFit(const GSLMultiFit &)=delete
~GSLMultiFit()
Destructor (no operations)
Definition GSLMultiFit.h:77
int TestDelta(double absTol, double relTol) const
test using abs and relative tolerance |dx| < absTol + relTol*|x| for every component
GSLMultiFitFunctionWrapper fFunc
gsl_multifit_fdfsolver * fSolver
const double * Gradient() const
gradient value at the minimum
std::string Name() const
GSLMultiFit(const gsl_multifit_fdfsolver_type *type=nullptr)
Default constructor No need to specify the type so far since only one solver exists so far.
Definition GSLMultiFit.h:61
GSLMultiFit & operator=(const GSLMultiFit &rhs)=delete
const gsl_multifit_fdfsolver_type * fType
const double * CovarMatrix() const
return covariance matrix of the parameters
const double * X() const
parameter values at the minimum
int Set(const std::vector< Func > &funcVec, const double *x)
set the solver parameters
Double_t x[n]
Definition legend1.C:17
Namespace for new Math classes and functions.
tbb::task_arena is an alias of tbb::interface7::task_arena, which doesn't allow to forward declare tb...