Logo ROOT   6.13/01
Reference Guide
TFumiliMinimizer.cxx
Go to the documentation of this file.
1 // @(#)root/fumili:$Id$
2 // Author: L. Moneta Wed Oct 25 16:28:55 2006
3 
4 /**********************************************************************
5  * *
6  * Copyright (c) 2006 LCG ROOT Math Team, CERN/PH-SFT *
7  * *
8  * *
9  **********************************************************************/
10 
11 // Implementation file for class TFumiliMinimizer
12 
13 #include "TFumiliMinimizer.h"
14 #include "Math/IFunction.h"
15 #include "Math/Util.h"
16 #include "TError.h"
17 
18 #include "TFumili.h"
19 
20 #include <iostream>
21 #include <cassert>
22 #include <algorithm>
23 #include <functional>
24 
25 
26 // setting USE_FUMILI_FUNCTION will use the Derivatives provided by Fumili
27 // instead of what proided in FitUtil::EvalChi2Residual
28 // t.d.: use still standard Chi2 but replace model function
29 // with a gradient function where gradient is computed by TFumili
30 // since TFumili knows the step size can calculate it better
31 // Derivative in FUmili are very fast (1 extra call for each parameter)
32 // + 1 function evaluation
33 //
34 //#define USE_FUMILI_FUNCTION
35 #ifdef USE_FUMILI_FUNCTION
36 bool gUseFumiliFunction = true;
37 //#include "FumiliFunction.h"
38 // fit method function used in TFumiliMinimizer
39 
41 #include "Fit/LogLikelihoodFCN.h"
42 #include "Fit/Chi2FCN.h"
43 #include "TF1.h"
44 #include "TFumili.h"
45 
46 template<class MethodFunc>
47 class FumiliFunction : public ROOT::Math::FitMethodFunction {
48 
50 
51 public:
52  FumiliFunction(TFumili * fumili, const ROOT::Math::FitMethodFunction * func) :
53  ROOT::Math::FitMethodFunction(func->NDim(), func->NPoints() ),
54  fFumili(fumili),
55  fObjFunc(0)
56  {
57  fObjFunc = dynamic_cast<const MethodFunc *>(func);
58  assert(fObjFunc != 0);
59 
60  // create TF1 class from model function
61  fModFunc = new TF1("modfunc",ROOT::Math::ParamFunctor( &fObjFunc->ModelFunction() ) );
62  fFumili->SetUserFunc(fModFunc);
63  }
64 
65  ROOT::Math::FitMethodFunction::Type_t Type() const { return fObjFunc->Type(); }
66 
67  FumiliFunction * Clone() const { return new FumiliFunction(fFumili, fObjFunc); }
68 
69 
70  // recalculate data elemet using Fumili stuff
71  double DataElement(const double * /*par */, unsigned int i, double * g) const {
72 
73  // parameter values are inside TFumili
74 
75  // suppose type is bin likelihood
76  unsigned int npar = fObjFunc->NDim();
77  double y = 0;
78  double invError = 0;
79  const double *x = fObjFunc->Data().GetPoint(i,y,invError);
80  double fval = fFumili->EvalTFN(g,const_cast<double *>( x));
81  fFumili->Derivatives(g, const_cast<double *>( x));
82 
83  if ( fObjFunc->Type() == ROOT::Math::FitMethodFunction::kLogLikelihood) {
84  double logPdf = y * ROOT::Math::Util::EvalLog( fval) - fval;
85  for (unsigned int k = 0; k < npar; ++k) {
86  g[k] *= ( y/fval - 1.) ;//* pdfval;
87  }
88 
89  // std::cout << "x = " << x[0] << " logPdf = " << logPdf << " grad";
90 // for (unsigned int ipar = 0; ipar < npar; ++ipar)
91 // std::cout << g[ipar] << "\t";
92 // std::cout << std::endl;
93 
94  return logPdf;
95  }
96  else if (fObjFunc->Type() == ROOT::Math::FitMethodFunction::kLeastSquare ) {
97  double resVal = (y-fval)*invError;
98  for (unsigned int k = 0; k < npar; ++k) {
99  g[k] *= -invError;
100  }
101  return resVal;
102  }
103 
104  return 0;
105  }
106 
107 
108 private:
109 
110  double DoEval(const double *x ) const {
111  return (*fObjFunc)(x);
112  }
113 
114  TFumili * fFumili;
115  const MethodFunc * fObjFunc;
116  TF1 * fModFunc;
117 
118 };
119 #else
120 bool gUseFumiliFunction = false;
121 #endif
122 //______________________________________________________________________________
123 //
124 // TFumiliMinimizer class implementing the ROOT::Math::Minimizer interface using
125 // TFumili.
126 // This class is normally instantiates using the plug-in manager
127 // (plug-in with name Fumili or TFumili)
128 // In addition the user can choose the minimizer algorithm: Migrad (the default one), Simplex, or Minimize (combined Migrad + Simplex)
129 //
130 //__________________________________________________________________________________________
131 
132 // initialize the static instances
133 
137 
138 
139 ClassImp(TFumiliMinimizer);
140 
141 
143  fDim(0),
144  fNFree(0),
145  fMinVal(0),
146  fEdm(-1),
147  fFumili(0)
148 {
149  // Constructor for TFumiliMinimier class
150 
151  // construct with npar = 0 (by default a value of 25 is used in TFumili for allocating the arrays)
152 #ifdef USE_STATIC_TMINUIT
153  // allocate here only the first time
154  if (fgFumili == 0) fgFumili = new TFumili(0);
155  fFumili = fgFumili;
156 #else
157  if (fFumili) delete fFumili;
158  fFumili = new TFumili(0);
159  fgFumili = fFumili;
160 #endif
161 
162 }
163 
164 
166 {
167  // Destructor implementation.
168  if (fFumili) delete fFumili;
169 }
170 
172  Minimizer()
173 {
174  // Implementation of copy constructor (it is private).
175 }
176 
178 {
179  // Implementation of assignment operator (private)
180  if (this == &rhs) return *this; // time saving self-test
181  return *this;
182 }
183 
184 
185 
187  // Set the objective function to be minimized, by passing a function object implement the
188  // basic multi-dim Function interface. In this case the derivatives will be
189  // calculated by Fumili
190 
191  // Here a TFumili instance is created since only at this point we know the number of parameters
192  // needed to create TFumili
193  fDim = func.NDim();
195 
196  // for Fumili the fit method function interface is required
197  const ROOT::Math::FitMethodFunction * fcnfunc = dynamic_cast<const ROOT::Math::FitMethodFunction *>(&func);
198  if (!fcnfunc) {
199  Error("SetFunction","Wrong Fit method function type used for Fumili");
200  return;
201  }
202  // assign to the static pointer (NO Thread safety here)
203  fgFunc = const_cast<ROOT::Math::FitMethodFunction *>(fcnfunc);
204  fgGradFunc = 0;
205  fFumili->SetFCN(&TFumiliMinimizer::Fcn);
206 
207 #ifdef USE_FUMILI_FUNCTION
208  if (gUseFumiliFunction) {
211  else if (fcnfunc->Type() == ROOT::Math::FitMethodFunction::kLeastSquare)
213  }
214 #endif
215 
216 }
217 
219  // Set the objective function to be minimized, by passing a function object implement the
220  // multi-dim gradient Function interface. In this case the function derivatives are provided
221  // by the user via this interface and there not calculated by Fumili.
222 
223  fDim = func.NDim();
225 
226  // for Fumili the fit method function interface is required
227  const ROOT::Math::FitMethodGradFunction * fcnfunc = dynamic_cast<const ROOT::Math::FitMethodGradFunction *>(&func);
228  if (!fcnfunc) {
229  Error("SetFunction","Wrong Fit method function type used for Fumili");
230  return;
231  }
232  // assign to the static pointer (NO Thread safety here)
233  fgFunc = 0;
234  fgGradFunc = const_cast<ROOT::Math::FitMethodGradFunction *>(fcnfunc);
235  fFumili->SetFCN(&TFumiliMinimizer::Fcn);
236 
237 }
238 
239 void TFumiliMinimizer::Fcn( int & , double * g , double & f, double * x , int /* iflag */) {
240  // implementation of FCN static function used internally by TFumili.
241  // Adapt IMultiGenFunction interface to TFumili FCN static function
242  f = TFumiliMinimizer::EvaluateFCN(const_cast<double*>(x),g);
243 }
244 
245 // void TFumiliMinimizer::FcnGrad( int &, double * g, double & f, double * x , int iflag ) {
246 // // implementation of FCN static function used internally by TFumili.
247 // // Adapt IMultiGradFunction interface to TFumili FCN static function in the case of user
248 // // provided gradient.
249 // ROOT::Math::IMultiGradFunction * gFunc = dynamic_cast<ROOT::Math::IMultiGradFunction *> ( fgFunc);
250 
251 // assert(gFunc != 0);
252 // f = gFunc->operator()(x);
253 
254 // // calculates also derivatives
255 // if (iflag == 2) gFunc->Gradient(x,g);
256 // }
257 
258 double TFumiliMinimizer::EvaluateFCN(const double * x, double * grad) {
259  // function callaed to evaluate the FCN at the value x
260  // calculates also the matrices of the second derivatives of the objective function needed by FUMILI
261 
262 
263  //typedef FumiliFCNAdapter::Function Function;
264 
265 
266 
267  // reset
268 // assert(grad.size() == npar);
269 // grad.assign( npar, 0.0);
270 // hess.assign( hess.size(), 0.0);
271 
272  double sum = 0;
273  unsigned int ndata = 0;
274  unsigned int npar = 0;
275  if (fgFunc) {
276  ndata = fgFunc->NPoints();
277  npar = fgFunc->NDim();
278  fgFunc->UpdateNCalls();
279  }
280  else if (fgGradFunc) {
281  ndata = fgGradFunc->NPoints();
282  npar = fgGradFunc->NDim();
284  }
285 
286  // eventually store this matrix as static member to optimize speed
287  std::vector<double> gf(npar);
288  std::vector<double> hess(npar*(npar+1)/2);
289 
290  // reset gradients
291  for (unsigned int ipar = 0; ipar < npar; ++ipar)
292  grad[ipar] = 0;
293 
294 
295  //loop on the data points
296 //#define DEBUG
297 #ifdef DEBUG
298  std::cout << "=============================================";
299  std::cout << "par = ";
300  for (unsigned int ipar = 0; ipar < npar; ++ipar)
301  std::cout << x[ipar] << "\t";
302  std::cout << std::endl;
303  if (fgFunc) std::cout << "type " << fgFunc->Type() << std::endl;
304 #endif
305 
306 
307  // assume for now least-square
308  // since TFumili doet not use errodef I must diveide chi2 by 2
311 
312  double fval = 0;
313  for (unsigned int i = 0; i < ndata; ++i) {
314  // calculate data element and gradient
315  // DataElement returns (f-y)/s and gf is derivatives of model function multiplied by (-1/sigma)
316  if (gUseFumiliFunction) {
317  fval = fgFunc->DataElement( x, i, &gf[0]);
318  }
319  else {
320  if (fgFunc != 0)
321  fval = fgFunc->DataElement(x, i, &gf[0]);
322  else
323  fval = fgGradFunc->DataElement(x, i, &gf[0]);
324  }
325 
326  // t.b.d should protect for bad values of fval
327  sum += fval*fval;
328 
329  // to be check (TFumili uses a factor of 1/2 for chi2)
330 
331  for (unsigned int j = 0; j < npar; ++j) {
332  grad[j] += fval * gf[j];
333  for (unsigned int k = j; k < npar; ++ k) {
334  int idx = j + k*(k+1)/2;
335  hess[idx] += gf[j] * gf[k];
336  }
337  }
338  }
339  }
342 
343 
344 
345  double fval = 0;
346 
347  //std::cout << "\t x " << x[0] << " " << x[1] << " " << x[2] << std::endl;
348 
349  for (unsigned int i = 0; i < ndata; ++i) {
350 
351  if (gUseFumiliFunction) {
352  fval = fgFunc->DataElement( x, i, &gf[0]);
353  }
354  else {
355  // calculate data element and gradient
356  if (fgFunc != 0)
357  fval = fgFunc->DataElement(x, i, &gf[0]);
358  else
359  fval = fgGradFunc->DataElement(x, i, &gf[0]);
360  }
361 
362  // protect for small values of fval
363  // std::cout << i << " " << fval << " log " << " grad " << gf[0] << " " << gf[1] << " " << gf[2] << std::endl;
364 // sum -= ROOT::Math::Util::EvalLog(fval);
365  sum -= fval;
366 
367  for (unsigned int j = 0; j < npar; ++j) {
368  double gfj = gf[j];// / fval;
369  grad[j] -= gfj;
370  for (unsigned int k = j; k < npar; ++ k) {
371  int idx = j + k*(k+1)/2;
372  hess[idx] += gfj * gf[k];// / (fval );
373  }
374  }
375  }
376  }
377  else {
378  Error("EvaluateFCN"," type of fit method is not supported, it must be chi2 or log-likelihood");
379  }
380 
381  // now TFumili excludes fixed prameter in second-derivative matrix
382  // ned to get them using the static instance of TFumili
383  double * zmatrix = fgFumili->GetZ();
384  double * pl0 = fgFumili->GetPL0(); // parameter limits
385  assert(zmatrix != 0);
386  assert(pl0 != 0);
387  unsigned int k = 0;
388  unsigned int l = 0;
389  for (unsigned int i = 0; i < npar; ++i) {
390  for (unsigned int j = 0; j <= i; ++j) {
391  if (pl0[i] > 0 && pl0[j] > 0) { // only for non-fixed parameters
392  zmatrix[l++] = hess[k];
393  }
394  k++;
395  }
396  }
397 
398 #ifdef DEBUG
399  std::cout << "FCN value " << sum << " grad ";
400  for (unsigned int ipar = 0; ipar < npar; ++ipar)
401  std::cout << grad[ipar] << "\t";
402  std::cout << std::endl << std::endl;
403 #endif
404 
405 
406  return 0.5*sum; // fumili multiply then by 2
407 
408 }
409 
410 
411 
412 bool TFumiliMinimizer::SetVariable(unsigned int ivar, const std::string & name, double val, double step) {
413  // set a free variable.
414  if (fFumili == 0) {
415  Error("SetVariableValue","invalid TFumili pointer. Set function first ");
416  return false;
417  }
418 #ifdef DEBUG
419  std::cout << "set variable " << ivar << " " << name << " value " << val << " step " << step << std::endl;
420 #endif
421 
422  int ierr = fFumili->SetParameter(ivar , name.c_str(), val, step, 0., 0. );
423  if (ierr) {
424  Error("SetVariable","Error for parameter %d ",ivar);
425  return false;
426  }
427  return true;
428 }
429 
430 bool TFumiliMinimizer::SetLimitedVariable(unsigned int ivar, const std::string & name, double val, double step, double lower, double upper) {
431  // set a limited variable.
432  if (fFumili == 0) {
433  Error("SetVariableValue","invalid TFumili pointer. Set function first ");
434  return false;
435  }
436 #ifdef DEBUG
437  std::cout << "set limited variable " << ivar << " " << name << " value " << val << " step " << step << std::endl;
438 #endif
439  int ierr = fFumili->SetParameter(ivar, name.c_str(), val, step, lower, upper );
440  if (ierr) {
441  Error("SetLimitedVariable","Error for parameter %d ",ivar);
442  return false;
443  }
444  return true;
445 }
446 #ifdef LATER
447 bool Fumili2Minimizer::SetLowerLimitedVariable(unsigned int ivar , const std::string & name , double val , double step , double lower ) {
448  // add a lower bounded variable as a double bound one, using a very large number for the upper limit
449  double s = val-lower;
450  double upper = s*1.0E15;
451  if (s != 0) upper = 1.0E15;
452  return SetLimitedVariable(ivar, name, val, step, lower,upper);
453 }
454 #endif
455 
456 
457 bool TFumiliMinimizer::SetFixedVariable(unsigned int ivar, const std::string & name, double val) {
458  // set a fixed variable.
459  if (fFumili == 0) {
460  Error("SetVariableValue","invalid TFumili pointer. Set function first ");
461  return false;
462  }
463 
464 
465  int ierr = fFumili->SetParameter(ivar, name.c_str(), val, 0., val, val );
466  fFumili->FixParameter(ivar);
467 
468 #ifdef DEBUG
469  std::cout << "Fix variable " << ivar << " " << name << " value " << std::endl;
470 #endif
471 
472  if (ierr) {
473  Error("SetFixedVariable","Error for parameter %d ",ivar);
474  return false;
475  }
476  return true;
477 }
478 
479 bool TFumiliMinimizer::SetVariableValue(unsigned int ivar, double val) {
480  // set the variable value
481  if (fFumili == 0) {
482  Error("SetVariableValue","invalid TFumili pointer. Set function first ");
483  return false;
484  }
485  TString name = fFumili->GetParName(ivar);
486  double oldval, verr, vlow, vhigh = 0;
487  int ierr = fFumili->GetParameter( ivar, &name[0], oldval, verr, vlow, vhigh);
488  if (ierr) {
489  Error("SetVariableValue","Error for parameter %d ",ivar);
490  return false;
491  }
492 #ifdef DEBUG
493  std::cout << "set variable " << ivar << " " << name << " value "
494  << val << " step " << verr << std::endl;
495 #endif
496 
497  ierr = fFumili->SetParameter(ivar , name , val, verr, vlow, vhigh );
498  if (ierr) {
499  Error("SetVariableValue","Error for parameter %d ",ivar);
500  return false;
501  }
502  return true;
503 }
504 
506  // perform the minimization using the algorithm chosen previously by the user
507  // By default Migrad is used.
508  // Return true if the found minimum is valid and update internal chached values of
509  // minimum values, errors and covariance matrix.
510 
511  if (fFumili == 0) {
512  Error("SetVariableValue","invalid TFumili pointer. Set function first ");
513  return false;
514  }
515 
516  // need to set static instance to be used when calling FCN
517  fgFumili = fFumili;
518 
519 
520  double arglist[10];
521 
522  // error cannot be set in TFumili (always the same)
523 // arglist[0] = ErrorUp();
524 // fFumili->ExecuteCommand("SET Err",arglist,1);
525 
526  int printlevel = PrintLevel();
527  // not implemented in TFumili yet
528  //arglist[0] = printlevel - 1;
529  //fFumili->ExecuteCommand("SET PRINT",arglist,1,ierr);
530 
531  // suppress warning in case Printlevel() == 0
532  if (printlevel == 0) fFumili->ExecuteCommand("SET NOW",arglist,0);
533  else fFumili->ExecuteCommand("SET WAR",arglist,0);
534 
535 
536  // minimize: use ExecuteCommand instead of Minimize to set tolerance and maxiter
537 
538  arglist[0] = MaxFunctionCalls();
539  arglist[1] = Tolerance();
540 
541  if (printlevel > 0)
542  std::cout << "Minimize using TFumili with tolerance = " << Tolerance()
543  << " max calls " << MaxFunctionCalls() << std::endl;
544 
545  int iret = fFumili->ExecuteCommand("MIGRAD",arglist,2);
546  fStatus = iret;
547  //int iret = fgFumili->Minimize();
548 
549  // Hesse and IMP not implemented
550 // // run improved if needed
551 // if (ierr == 0 && fType == ROOT::Fumili::kMigradImproved)
552 // fFumili->mnexcm("IMPROVE",arglist,1,ierr);
553 
554 // // check if Hesse needs to be run
555 // if (ierr == 0 && IsValidError() ) {
556 // fFumili->mnexcm("HESSE",arglist,1,ierr);
557 // }
558 
559 
560  int ntot;
561  int nfree;
562  double errdef = 0; // err def is not used by Fumili
563  fFumili->GetStats(fMinVal,fEdm,errdef,nfree,ntot);
564 
565  if (printlevel > 0)
566  fFumili->PrintResults(printlevel,fMinVal);
567 
568 
569  assert (static_cast<unsigned int>(ntot) == fDim);
570  assert( nfree == fFumili->GetNumberFreeParameters() );
571  fNFree = nfree;
572 
573 
574  // get parameter values and correlation matrix
575  // fumili stores only lower part of diagonal matrix of the free parameters
576  fParams.resize( fDim);
577  fErrors.resize( fDim);
578  fCovar.resize(fDim*fDim);
579  const double * cv = fFumili->GetCovarianceMatrix();
580  unsigned int l = 0;
581  for (unsigned int i = 0; i < fDim; ++i) {
582  fParams[i] = fFumili->GetParameter( i );
583  fErrors[i] = fFumili->GetParError( i );
584 
585  if ( !fFumili->IsFixed(i) ) {
586  for (unsigned int j = 0; j <=i ; ++j) {
587  if ( !fFumili->IsFixed(j) ) {
588  fCovar[i*fDim + j] = cv[l];
589  fCovar[j*fDim + i] = fCovar[i*fDim + j];
590  l++;
591  }
592  }
593  }
594  }
595 
596  return (iret==0) ? true : false;
597 }
598 
599 
600 // } // end namespace Fit
601 
602 // } // end namespace ROOT
603 
virtual bool SetVariable(unsigned int ivar, const std::string &name, double val, double step)
set free variable
virtual Int_t GetStats(Double_t &amin, Double_t &edm, Double_t &errdef, Int_t &nvpar, Int_t &nparx) const
return global fit parameters amin : chisquare edm : estimated distance to minimum errdef nvpar : numb...
Definition: TFumili.cxx:912
virtual void UpdateNCalls() const
update number of calls
Param Functor class for Multidimensional functions.
Definition: ParamFunctor.h:274
Type_t
enumeration specyfing the possible fit method types
Interface (abstract class) for multi-dimensional functions providing a gradient calculation.
Definition: IFunction.h:326
virtual Double_t * GetCovarianceMatrix() const
return a pointer to the covariance matrix
Definition: TFumili.cxx:791
Namespace for new ROOT classes and functions.
Definition: TFoamSampler.h:19
bool gUseFumiliFunction
virtual void FixParameter(Int_t ipar)
Fixes parameter number ipar.
Definition: TFumili.cxx:781
virtual Double_t GetParameter(Int_t ipar) const
return current value of parameter ipar
Definition: TFumili.cxx:843
virtual unsigned int NPoints() const
return the number of data points used in evaluating the function
static double EvaluateFCN(const double *x, double *g)
implementation of FCN for Fumili when user provided gradient is used
virtual void PrintResults(Int_t k, Double_t p) const
Prints fit results.
Definition: TFumili.cxx:1509
Double_t * GetPL0() const
Definition: TFumili.h:103
TFumiliMinimizer class: minimizer implementation based on TFumili.
Minimizer()
Default constructor.
Definition: Minimizer.h:85
static ROOT::Math::FitMethodFunction * fgFunc
you should not use this method at all Int_t y
Definition: TRolke.cxx:630
int PrintLevel() const
minimizer configuration parameters
Definition: Minimizer.h:411
TFumiliMinimizer(int dummy=0)
Default constructor (an argument is needed by plug-in manager)
virtual Int_t SetParameter(Int_t ipar, const char *parname, Double_t value, Double_t verr, Double_t vlow, Double_t vhigh)
Sets for prameter number ipar initial parameter value, name parname, initial error verr and limits vl...
Definition: TFumili.cxx:1672
class evaluating the log likelihood for binned Poisson likelihood fits it is template to distinguish ...
virtual const char * GetParName(Int_t ipar) const
return name of parameter ipar
Definition: TFumili.cxx:880
std::vector< double > fParams
virtual Type_t Type() const
return the type of method, override if needed
virtual bool SetFixedVariable(unsigned int, const std::string &, double)
set fixed variable (override if minimizer supports them )
virtual bool Minimize()
method to perform the minimization
T EvalLog(T x)
safe evaluation of log(x) with a protections against negative or zero argument to the log smooth line...
Definition: Util.h:56
double Tolerance() const
absolute tolerance
Definition: Minimizer.h:420
virtual Double_t GetParError(Int_t ipar) const
return error of parameter ipar
Definition: TFumili.cxx:834
Documentation for the abstract class IBaseFunctionMultiDim.
Definition: IFunction.h:62
Chi2FCN class for binnned fits using the least square methods.
Definition: Chi2FCN.h:49
virtual bool SetLimitedVariable(unsigned int ivar, const std::string &name, double val, double step, double, double)
set upper/lower limited variable (override if minimizer supports them )
std::vector< double > fCovar
virtual Bool_t IsFixed(Int_t ipar) const
Return kTRUE if parameter ipar is fixed, kFALSE othersise)
Definition: TFumili.cxx:1076
virtual Int_t GetNumberFreeParameters() const
return the number of free parameters
Definition: TFumili.cxx:822
FunctionType::BaseFunc BaseFunction
virtual unsigned int NDim() const
Number of dimension (parameters) .
virtual Int_t ExecuteCommand(const char *command, Double_t *args, Int_t nargs)
Execute MINUIT commands.
Definition: TFumili.cxx:394
~TFumiliMinimizer()
Destructor (no operations)
* x
Deprecated and error prone model selection interface.
Definition: TRolke.cxx:630
Double_t * GetZ() const
Definition: TFumili.h:110
FitMethodFunction class Interface for objective functions (like chi2 and likelihood used in the fit) ...
Definition: Fitter.h:40
virtual void SetFunction(const ROOT::Math::IMultiGenFunction &func)
set the function to minimize
Namespace for new Math classes and functions.
unsigned int MaxFunctionCalls() const
max number of function calls
Definition: Minimizer.h:414
static ROOT::Math::FitMethodGradFunction * fgGradFunc
static TFumili * fgFumili
std::vector< double > fErrors
TFumiliMinimizer & operator=(const TFumiliMinimizer &rhs)
Assignment operator.
void SetParNumber(Int_t ParNum)
Definition: TFumili.h:120
unsigned int fNFree
virtual bool SetVariableValue(unsigned int ivar, double val)
set the value of an existing variable
static void Fcn(int &, double *, double &f, double *, int)
implementation of FCN for Fumili
virtual double DataElement(const double *x, unsigned int i, double *g=0) const =0
method returning the data i-th contribution to the fit objective function For example the residual fo...
virtual unsigned int NDim() const =0
Retrieve the dimension of the function.
BasicFitMethodFunction< ROOT::Math::IMultiGenFunction > FitMethodFunction
Definition: Fitter.h:40