Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
HFitImpl.cxx
Go to the documentation of this file.
1// new HFit function
2//______________________________________________________________________________
3
4
5#include "TH1.h"
6#include "TF1.h"
7#include "TF2.h"
8#include "TF3.h"
9#include "TError.h"
10#include "TGraph.h"
11#include "TMultiGraph.h"
12#include "TGraph2D.h"
13#include "THnBase.h"
14
15#include "Fit/Fitter.h"
16#include "Fit/FitConfig.h"
17#include "Fit/BinData.h"
18#include "Fit/UnBinData.h"
19#include "Fit/Chi2FCN.h"
21#include "HFitInterface.h"
23#include "Math/Minimizer.h"
24
25#include "Math/WrappedTF1.h"
27
28#include "TList.h"
29#include "TMath.h"
30#include "TROOT.h"
31
32#include "TVirtualPad.h" // for gPad
33
34#include "TBackCompFitter.h"
35#include "TFitResultPtr.h"
36#include "TFitResult.h"
37
38#include <stdlib.h>
39#include <cmath>
40#include <memory>
41#include <limits>
42
43//#define DEBUG
44
45// utility functions used in TH1::Fit
46
47namespace HFit {
48
49
50 int GetDimension(const TH1 * h1) { return h1->GetDimension(); }
51 int GetDimension(const TGraph * ) { return 1; }
52 int GetDimension(const TMultiGraph * ) { return 1; }
53 int GetDimension(const TGraph2D * ) { return 2; }
54 int GetDimension(const THnBase * s1) { return s1->GetNdimensions(); }
55
56 int CheckFitFunction(const TF1 * f1, int hdim);
57
58
59 void GetFunctionRange(const TF1 & f1, ROOT::Fit::DataRange & range);
60
61 void FitOptionsMake(const char *option, Foption_t &fitOption);
62
63 void CheckGraphFitOptions(Foption_t &fitOption);
64
65
71
72
73 template <class FitObject>
74 TFitResultPtr Fit(FitObject * h1, TF1 *f1 , Foption_t & option , const ROOT::Math::MinimizerOptions & moption, const char *goption, ROOT::Fit::DataRange & range);
75
76 template <class FitObject>
77 void StoreAndDrawFitFunction(FitObject * h1, TF1 * f1, const ROOT::Fit::DataRange & range, bool, bool, const char *goption);
78
79 template <class FitObject>
80 double ComputeChi2(const FitObject & h1, TF1 &f1, bool useRange, bool usePL );
81
82
83
84}
85
86int HFit::CheckFitFunction(const TF1 * f1, int dim) {
87 // Check validity of fitted function
88 if (!f1) {
89 Error("Fit", "function may not be null pointer");
90 return -1;
91 }
92 if (f1->IsZombie()) {
93 Error("Fit", "function is zombie");
94 return -2;
95 }
96
97 int npar = f1->GetNpar();
98 if (npar <= 0) {
99 Error("Fit", "function %s has illegal number of parameters = %d", f1->GetName(), npar);
100 return -3;
101 }
102
103 // Check that function has same dimension as histogram
104 if (f1->GetNdim() > dim) {
105 Error("Fit","function %s dimension, %d, is greater than fit object dimension, %d",
106 f1->GetName(), f1->GetNdim(), dim);
107 return -4;
108 }
109 if (f1->GetNdim() < dim-1) {
110 Error("Fit","function %s dimension, %d, is smaller than fit object dimension -1, %d",
111 f1->GetName(), f1->GetNdim(), dim);
112 return -5;
113 }
114
115 return 0;
116
117}
118
119
121 // get the range form the function and fill and return the DataRange object
122 Double_t fxmin, fymin, fzmin, fxmax, fymax, fzmax;
123 f1.GetRange(fxmin, fymin, fzmin, fxmax, fymax, fzmax);
124 // support only one range - so add only if was not set before
125 if (range.Size(0) == 0) range.AddRange(0,fxmin,fxmax);
126 if (range.Size(1) == 0) range.AddRange(1,fymin,fymax);
127 if (range.Size(2) == 0) range.AddRange(2,fzmin,fzmax);
128 return;
129}
130
131
132template<class FitObject>
133TFitResultPtr HFit::Fit(FitObject * h1, TF1 *f1 , Foption_t & fitOption , const ROOT::Math::MinimizerOptions & minOption, const char *goption, ROOT::Fit::DataRange & range)
134{
135 // perform fit of histograms, or graphs using new fitting classes
136 // use same routines for fitting both graphs and histograms
137
138#ifdef DEBUG
139 printf("fit function %s\n",f1->GetName() );
140#endif
141
142 // replacement function using new fitter
143 int hdim = HFit::GetDimension(h1);
144 int iret = HFit::CheckFitFunction(f1, hdim);
145 if (iret != 0) return iret;
146
147
148
149 // integral option is not supported in this case
150 if (f1->GetNdim() < hdim ) {
151 if (fitOption.Integral) Info("Fit","Ignore Integral option. Model function dimension is less than the data object dimension");
152 if (fitOption.Like) Info("Fit","Ignore Likelihood option. Model function dimension is less than the data object dimension");
153 fitOption.Integral = 0;
154 fitOption.Like = 0;
155 }
156
157 Int_t special = f1->GetNumber();
158 Bool_t linear = f1->IsLinear();
159 Int_t npar = f1->GetNpar();
160 if (special==299+npar) linear = kTRUE; // for polynomial functions
161 // do not use linear fitter in these case
162 if (fitOption.Bound || fitOption.Like || fitOption.Errors || fitOption.Gradient || fitOption.More || fitOption.User|| fitOption.Integral || fitOption.Minuit)
163 linear = kFALSE;
164
165 // create an empty TFitResult
166 std::shared_ptr<TFitResult> tfr(new TFitResult() );
167 // create the fitter from an empty fit result
168 std::shared_ptr<ROOT::Fit::Fitter> fitter(new ROOT::Fit::Fitter(std::static_pointer_cast<ROOT::Fit::FitResult>(tfr) ) );
169 ROOT::Fit::FitConfig & fitConfig = fitter->Config();
170
171 // create options
173 opt.fIntegral = fitOption.Integral;
174 opt.fUseRange = fitOption.Range;
175 opt.fExpErrors = fitOption.PChi2; // pearson chi2 with expected errors
176 if (fitOption.Like || fitOption.PChi2) opt.fUseEmpty = true; // use empty bins in log-likelihood fits
177 if (special==300) opt.fCoordErrors = false; // no need to use coordinate errors in a pol0 fit
178 if (fitOption.NoErrX) opt.fCoordErrors = false; // do not use coordinate errors when requested
179 if (fitOption.W1 ) opt.fErrors1 = true;
180 if (fitOption.W1 > 1) opt.fUseEmpty = true; // use empty bins with weight=1
181
182 if (fitOption.BinVolume) {
183 opt.fBinVolume = true; // scale by bin volume
184 if (fitOption.BinVolume == 2) opt.fNormBinVolume = true; // scale by normalized bin volume
185 }
186
187 if (opt.fUseRange) {
188#ifdef DEBUG
189 printf("use range \n" );
190#endif
192 }
193#ifdef DEBUG
194 printf("range size %d\n", range.Size(0) );
195 if (range.Size(0)) {
196 double x1; double x2; range.GetRange(0,x1,x2);
197 printf(" range in x = [%f,%f] \n",x1,x2);
198 }
199#endif
200
201 // fill data
202 std::shared_ptr<ROOT::Fit::BinData> fitdata(new ROOT::Fit::BinData(opt,range) );
203 ROOT::Fit::FillData(*fitdata, h1, f1);
204 if (fitdata->Size() == 0 ) {
205 Warning("Fit","Fit data is empty ");
206 return -1;
207 }
208
209#ifdef DEBUG
210 printf("HFit:: data size is %d \n",fitdata->Size());
211 for (unsigned int i = 0; i < fitdata->Size(); ++i) {
212 if (fitdata->NDim() == 1) printf(" x[%d] = %f - value = %f \n", i,*(fitdata->Coords(i)),fitdata->Value(i) );
213 }
214#endif
215
216 // switch off linear fitting in case data has coordinate errors and the option is set
217 if (fitdata->GetErrorType() == ROOT::Fit::BinData::kCoordError && fitdata->Opt().fCoordErrors ) linear = false;
218 // linear fit cannot be done also in case of asymmetric errors
219 if (fitdata->GetErrorType() == ROOT::Fit::BinData::kAsymError && fitdata->Opt().fAsymErrors ) linear = false;
220
221 // this functions use the TVirtualFitter
222 if (special != 0 && !fitOption.Bound && !linear) {
223 if (special == 100) ROOT::Fit::InitGaus (*fitdata,f1); // gaussian
224 else if (special == 110 || special == 112) ROOT::Fit::Init2DGaus(*fitdata,f1); // 2D gaussians ( xygaus or bigaus)
225 else if (special == 400) ROOT::Fit::InitGaus (*fitdata,f1); // landau (use the same)
226 else if (special == 410) ROOT::Fit::Init2DGaus(*fitdata,f1); // 2D landau (use the same)
227
228 else if (special == 200) ROOT::Fit::InitExpo (*fitdata, f1); // exponential
229
230 }
231
232
233 // set the fit function
234 // if option grad is specified use gradient
235 if ( (linear || fitOption.Gradient) )
236 fitter->SetFunction(ROOT::Math::WrappedMultiTF1(*f1));
237#ifdef R__HAS_VECCORE
238 else if(f1->IsVectorized())
240#endif
241 else
242 fitter->SetFunction(static_cast<const ROOT::Math::IParamMultiFunction &>(ROOT::Math::WrappedMultiTF1(*f1) ) );
243
244 // error normalization in case of zero error in the data
245 if (fitdata->GetErrorType() == ROOT::Fit::BinData::kNoError) fitConfig.SetNormErrors(true);
246 // error normalization also in case of W or WW options (weights = 1)
247 if (fitdata->Opt().fErrors1) fitConfig.SetNormErrors(true);
248 // normalize errors also in case you are fitting a Ndim histo with a N-1 function
249 if (int(fitdata->NDim()) == hdim -1 ) fitConfig.SetNormErrors(true);
250
251
252 // here need to get some static extra information (like max iterations, error def, etc...)
253
254
255 // parameter settings and transfer the parameters values, names and limits from the functions
256 // is done automatically in the Fitter.cxx
257 for (int i = 0; i < npar; ++i) {
258 ROOT::Fit::ParameterSettings & parSettings = fitConfig.ParSettings(i);
259
260 // check limits
261 double plow,pup;
262 f1->GetParLimits(i,plow,pup);
263 if (plow*pup != 0 && plow >= pup) { // this is a limitation - cannot fix a parameter to zero value
264 parSettings.Fix();
265 }
266 else if (plow < pup ) {
267 if (!TMath::Finite(pup) && TMath::Finite(plow) )
268 parSettings.SetLowerLimit(plow);
269 else if (!TMath::Finite(plow) && TMath::Finite(pup) )
270 parSettings.SetUpperLimit(pup);
271 else
272 parSettings.SetLimits(plow,pup);
273 }
274
275 // set the parameter step size (by default are set to 0.3 of value)
276 // if function provides meaningful error values
277 double err = f1->GetParError(i);
278 if ( err > 0)
279 parSettings.SetStepSize(err);
280 else if (plow < pup && TMath::Finite(plow) && TMath::Finite(pup) ) { // in case of limits improve step sizes
281 double step = 0.1 * (pup - plow);
282 // check if value is not too close to limit otherwise trim value
283 if ( parSettings.Value() < pup && pup - parSettings.Value() < 2 * step )
284 step = (pup - parSettings.Value() ) / 2;
285 else if ( parSettings.Value() > plow && parSettings.Value() - plow < 2 * step )
286 step = (parSettings.Value() - plow ) / 2;
287
288 parSettings.SetStepSize(step);
289 }
290
291
292 }
293
294 // needed for setting precision ?
295 // - Compute sum of squares of errors in the bin range
296 // should maybe use stat[1] ??
297 // Double_t ey, sumw2=0;
298// for (i=hxfirst;i<=hxlast;i++) {
299// ey = GetBinError(i);
300// sumw2 += ey*ey;
301// }
302
303
304 // set all default minimizer options (tolerance, max iterations, etc..)
305 fitConfig.SetMinimizerOptions(minOption);
306
307 // specific print level options
308 if (fitOption.Verbose) fitConfig.MinimizerOptions().SetPrintLevel(3);
309 if (fitOption.Quiet) fitConfig.MinimizerOptions().SetPrintLevel(0);
310
311 // specific minimizer options depending on minimizer
312 if (linear) {
313 if (fitOption.Robust ) {
314 // robust fitting
315 std::string type = "Robust";
316 // if an h is specified print out the value adding to the type
317 if (fitOption.hRobust > 0 && fitOption.hRobust < 1.)
318 type += " (h=" + ROOT::Math::Util::ToString(fitOption.hRobust) + ")";
319 fitConfig.SetMinimizer("Linear",type.c_str());
320 fitConfig.MinimizerOptions().SetTolerance(fitOption.hRobust); // use tolerance for passing robust parameter
321 }
322 else
323 fitConfig.SetMinimizer("Linear","");
324 }
325 else {
326 if (fitOption.More) fitConfig.SetMinimizer("Minuit","MigradImproved");
327 }
328
329
330 // check if Error option (run Hesse and Minos) then
331 if (fitOption.Errors) {
332 // run Hesse and Minos
333 fitConfig.SetParabErrors(true);
334 fitConfig.SetMinosErrors(true);
335 }
336
337
338 // do fitting
339
340#ifdef DEBUG
341 if (fitOption.Like) printf("do likelihood fit...\n");
342 if (linear) printf("do linear fit...\n");
343 printf("do now fit...\n");
344#endif
345
346 bool fitok = false;
347
348
349 // check if can use option user
350 //typedef void (* MinuitFCN_t )(int &npar, double *gin, double &f, double *u, int flag);
351 TVirtualFitter::FCNFunc_t userFcn = 0;
352 if (fitOption.User && TVirtualFitter::GetFitter() ) {
353 userFcn = (TVirtualFitter::GetFitter())->GetFCN();
354 (TVirtualFitter::GetFitter())->SetUserFunc(f1);
355 }
356
357
358 if (fitOption.User && userFcn) // user provided fit objective function
359 fitok = fitter->FitFCN( userFcn );
360 else if (fitOption.Like) {// likelihood fit
361 // perform a weighted likelihood fit by applying weight correction to errors
362 bool weight = ((fitOption.Like & 2) == 2);
363 fitConfig.SetWeightCorrection(weight);
364 bool extended = ((fitOption.Like & 4 ) != 4 );
365 //if (!extended) Info("HFitImpl","Do a not -extended binned fit");
366
367 // pass fitdata as a shared pointer so ownership is shared with Fitter and FitResult class
368 fitok = fitter->LikelihoodFit(fitdata, extended, fitOption.ExecPolicy);
369 }
370 else{ // standard least square fit
371 fitok = fitter->Fit(fitdata, fitOption.ExecPolicy);
372 }
373 if ( !fitok && !fitOption.Quiet )
374 Warning("Fit","Abnormal termination of minimization.");
375 iret |= !fitok;
376
377
378 const ROOT::Fit::FitResult & fitResult = fitter->Result();
379 // one could set directly the fit result in TF1
380 iret = fitResult.Status();
381 if (!fitResult.IsEmpty() ) {
382 // set in f1 the result of the fit
383 f1->SetChisquare(fitResult.Chi2() );
384 f1->SetNDF(fitResult.Ndf() );
385 f1->SetNumberFitPoints(fitdata->Size() );
386
387 assert((Int_t)fitResult.Parameters().size() >= f1->GetNpar() );
388 f1->SetParameters( const_cast<double*>(&(fitResult.Parameters().front())));
389 if ( int( fitResult.Errors().size()) >= f1->GetNpar() )
390 f1->SetParErrors( &(fitResult.Errors().front()) );
391
392
393 }
394
395// - Store fitted function in histogram functions list and draw
396 if (!fitOption.Nostore) {
398 HFit::StoreAndDrawFitFunction(h1, f1, range, !fitOption.Plus, !fitOption.Nograph, goption);
399 }
400
401 // print the result
402 // if using Fitter class must be done here
403 // use old style Minuit for TMinuit and if no corrections have been applied
404 if (!fitOption.Quiet) {
405 if (fitter->GetMinimizer() && fitConfig.MinimizerType() == "Minuit" &&
406 !fitConfig.NormalizeErrors() && fitOption.Like <= 1) {
407 fitter->GetMinimizer()->PrintResults(); // use old style Minuit
408 }
409 else {
410 // print using FitResult class
411 if (fitOption.Verbose) fitResult.PrintCovMatrix(std::cout);
412 fitResult.Print(std::cout);
413 }
414 }
415
416
417 // store result in the backward compatible VirtualFitter
418 // in case multi-thread is not enabled
419 if (!gGlobalMutex) {
421 TBackCompFitter * bcfitter = new TBackCompFitter(fitter, fitdata);
422 bcfitter->SetFitOption(fitOption);
423 bcfitter->SetObjectFit(h1);
424 bcfitter->SetUserFunc(f1);
426 if (userFcn) {
427 bcfitter->SetFCN(userFcn);
428 // for interpreted FCN functions
429 if (lastFitter->GetMethodCall() ) bcfitter->SetMethodCall(lastFitter->GetMethodCall() );
430 }
431
432 // delete last fitter if it has been created here before
433 if (lastFitter) {
434 TBackCompFitter * lastBCFitter = dynamic_cast<TBackCompFitter *> (lastFitter);
435 if (lastBCFitter && lastBCFitter->TestBit(TBackCompFitter::kCanDeleteLast) )
436 delete lastBCFitter;
437 }
438 //N.B= this might create a memory leak if user does not delete the fitter they create
439 TVirtualFitter::SetFitter( bcfitter );
440 }
441
442 // use old-style for printing the results
443 // if (fitOption.Verbose) bcfitter->PrintResults(2,0.);
444 // else if (!fitOption.Quiet) bcfitter->PrintResults(1,0.);
445
446 if (fitOption.StoreResult)
447 {
448 TString name = "TFitResult-";
449 name = name + h1->GetName() + "-" + f1->GetName();
450 TString title = "TFitResult-";
451 title += h1->GetTitle();
452 tfr->SetName(name);
453 tfr->SetTitle(title);
454 return TFitResultPtr(tfr);
455 }
456 else
457 return TFitResultPtr(iret);
458}
459
460
462 // get range from histogram and update the DataRange class
463 // if a ranges already exist in that dimension use that one
464
465 Int_t ndim = GetDimension(h1);
466
467 double xmin = 0, xmax = 0, ymin = 0, ymax = 0, zmin = 0, zmax = 0;
468 if (range.Size(0) == 0) {
469 TAxis & xaxis = *(h1->GetXaxis());
470 Int_t hxfirst = xaxis.GetFirst();
471 Int_t hxlast = xaxis.GetLast();
472 Double_t binwidx = xaxis.GetBinWidth(hxlast);
473 xmin = xaxis.GetBinLowEdge(hxfirst);
474 xmax = xaxis.GetBinLowEdge(hxlast) +binwidx;
475 range.AddRange(xmin,xmax);
476 }
477
478 if (ndim > 1) {
479 if (range.Size(1) == 0) {
480 TAxis & yaxis = *(h1->GetYaxis());
481 Int_t hyfirst = yaxis.GetFirst();
482 Int_t hylast = yaxis.GetLast();
483 Double_t binwidy = yaxis.GetBinWidth(hylast);
484 ymin = yaxis.GetBinLowEdge(hyfirst);
485 ymax = yaxis.GetBinLowEdge(hylast) +binwidy;
486 range.AddRange(1,ymin,ymax);
487 }
488 }
489 if (ndim > 2) {
490 if (range.Size(2) == 0) {
491 TAxis & zaxis = *(h1->GetZaxis());
492 Int_t hzfirst = zaxis.GetFirst();
493 Int_t hzlast = zaxis.GetLast();
494 Double_t binwidz = zaxis.GetBinWidth(hzlast);
495 zmin = zaxis.GetBinLowEdge(hzfirst);
496 zmax = zaxis.GetBinLowEdge(hzlast) +binwidz;
497 range.AddRange(2,zmin,zmax);
498 }
499 }
500#ifdef DEBUG
501 std::cout << "xmin,xmax" << xmin << " " << xmax << std::endl;
502#endif
503
504}
505
507 // get range for graph (used sub-set histogram)
508 // N.B. : this is different than in previous implementation of TGraph::Fit where range used was from xmin to xmax.
509 TH1 * h1 = gr->GetHistogram();
510 // an histogram is normally always returned for a TGraph
511 if (h1) HFit::GetDrawingRange(h1, range);
512}
514 // get range for multi-graph (used sub-set histogram)
515 // N.B. : this is different than in previous implementation of TMultiGraph::Fit where range used was from data xmin to xmax.
516 TH1 * h1 = mg->GetHistogram();
517 if (h1) {
519 }
520 else if (range.Size(0) == 0) {
521 // compute range from all the TGraph's belonging to the MultiGraph
522 double xmin = std::numeric_limits<double>::infinity();
523 double xmax = -std::numeric_limits<double>::infinity();
524 TIter next(mg->GetListOfGraphs() );
525 TGraph * g = 0;
526 while ( (g = (TGraph*) next() ) ) {
527 double x1 = 0, x2 = 0, y1 = 0, y2 = 0;
528 g->ComputeRange(x1,y1,x2,y2);
529 if (x1 < xmin) xmin = x1;
530 if (x2 > xmax) xmax = x2;
531 }
532 range.AddRange(xmin,xmax);
533 }
534}
536 // get range for graph2D (used sub-set histogram)
537 // N.B. : this is different than in previous implementation of TGraph2D::Fit. There range used was always(0,0)
538 // cannot use TGraph2D::GetHistogram which makes an interpolation
539 //TH1 * h1 = gr->GetHistogram();
540 //if (h1) HFit::GetDrawingRange(h1, range);
541 // not very efficient (t.b.i.)
542 if (range.Size(0) == 0) {
543 double xmin = gr->GetXmin();
544 double xmax = gr->GetXmax();
545 range.AddRange(0,xmin,xmax);
546 }
547 if (range.Size(1) == 0) {
548 double ymin = gr->GetYmin();
549 double ymax = gr->GetYmax();
550 range.AddRange(1,ymin,ymax);
551 }
552}
553
555 // get range from histogram and update the DataRange class
556 // if a ranges already exist in that dimension use that one
557
558 Int_t ndim = GetDimension(s1);
559
560 for ( int i = 0; i < ndim; ++i ) {
561 if ( range.Size(i) == 0 ) {
562 TAxis *axis = s1->GetAxis(i);
563 range.AddRange(i, axis->GetXmin(), axis->GetXmax());
564 }
565 }
566}
567
568template<class FitObject>
569void HFit::StoreAndDrawFitFunction(FitObject * h1, TF1 * f1, const ROOT::Fit::DataRange & range, bool delOldFunction, bool drawFunction, const char *goption) {
570// - Store fitted function in histogram functions list and draw
571// should have separate functions for 1,2,3d ? t.b.d in case
572
573#ifdef DEBUG
574 std::cout <<"draw and store fit function " << f1->GetName() << std::endl;
575#endif
576
577
578 Int_t ndim = GetDimension(h1);
579 double xmin = 0, xmax = 0, ymin = 0, ymax = 0, zmin = 0, zmax = 0;
580 if (range.Size(0) ) range.GetRange(0,xmin,xmax);
581 if (range.Size(1) ) range.GetRange(1,ymin,ymax);
582 if (range.Size(2) ) range.GetRange(2,zmin,zmax);
583
584
585#ifdef DEBUG
586 std::cout <<"draw and store fit function " << f1->GetName()
587 << " Range in x = [ " << xmin << " , " << xmax << " ]" << std::endl;
588#endif
589
590 TList * funcList = h1->GetListOfFunctions();
591 if (funcList == 0){
592 Error("StoreAndDrawFitFunction","Function list has not been created - cannot store the fitted function");
593 return;
594 }
595
596 // delete the function in the list only if
597 // the function we are fitting is not in that list
598 // If this is the case we re-use that function object and
599 // we do not create a new one (if delOldFunction is true)
600 bool reuseOldFunction = false;
601 if (delOldFunction) {
602 TIter next(funcList, kIterBackward);
603 TObject *obj;
604 while ((obj = next())) {
605 if (obj->InheritsFrom(TF1::Class())) {
606 if (obj != f1) {
607 funcList->Remove(obj);
608 delete obj;
609 }
610 else {
611 reuseOldFunction = true;
612 }
613 }
614 }
615 }
616
617 TF1 *fnew1 = 0;
618 TF2 *fnew2 = 0;
619 TF3 *fnew3 = 0;
620
621 // copy TF1 using TClass to avoid slicing in case of derived classes
622 if (ndim < 2) {
623 if (!reuseOldFunction) {
624 fnew1 = (TF1*)f1->IsA()->New();
625 R__ASSERT(fnew1);
626 f1->Copy(*fnew1);
627 funcList->Add(fnew1);
628 }
629 else {
630 fnew1 = f1;
631 }
632 fnew1->SetParent( h1 );
633 fnew1->SetRange(xmin,xmax);
634 fnew1->Save(xmin,xmax,0,0,0,0);
635 if (!drawFunction) fnew1->SetBit(TF1::kNotDraw);
636 fnew1->AddToGlobalList(false);
637 } else if (ndim < 3) {
638 if (!reuseOldFunction) {
639 fnew2 = (TF2*)f1->IsA()->New();
640 R__ASSERT(fnew2);
641 f1->Copy(*fnew2);
642 funcList->Add(fnew2);
643 }
644 else {
645 fnew2 = dynamic_cast<TF2*>(f1);
646 R__ASSERT(fnew2);
647 }
648 fnew2->SetRange(xmin,ymin,xmax,ymax);
649 fnew2->SetParent( h1 );
650 fnew2->Save(xmin,xmax,ymin,ymax,0,0);
651 if (!drawFunction) fnew2->SetBit(TF1::kNotDraw);
652 fnew2->AddToGlobalList(false);
653 } else {
654 if (!reuseOldFunction) {
655 fnew3 = (TF3*)f1->IsA()->New();
656 R__ASSERT(fnew3);
657 f1->Copy(*fnew3);
658 funcList->Add(fnew3);
659 }
660 else {
661 fnew3 = dynamic_cast<TF3*>(f1);
662 R__ASSERT(fnew3);
663 }
664 fnew3->SetRange(xmin,ymin,zmin,xmax,ymax,zmax);
665 fnew3->SetParent( h1 );
666 fnew3->Save(xmin,xmax,ymin,ymax,zmin,zmax);
667 if (!drawFunction) fnew3->SetBit(TF1::kNotDraw);
668 fnew3->AddToGlobalList(false);
669 }
670 if (h1->TestBit(kCanDelete)) return;
671 // draw only in case of histograms
672 if (drawFunction && ndim < 3 && h1->InheritsFrom(TH1::Class() ) ) {
673 // no need to re-draw the histogram if the histogram is already in the pad
674 // in that case the function will be just drawn (if option N is not set)
675 if (!gPad || (gPad && gPad->GetListOfPrimitives()->FindObject(h1) == NULL ) )
676 h1->Draw(goption);
677 }
678 if (gPad) gPad->Modified(); // this is not in TH1 code (needed ??)
679
680 return;
681}
682
683
684void ROOT::Fit::FitOptionsMake(EFitObjectType type, const char *option, Foption_t &fitOption) {
685 // - Decode list of options into fitOption (used by both TGraph and TH1)
686 // works for both histograms and graph depending on the enum FitObjectType defined in HFit
689 }
690
691 if (option == 0) return;
692 if (!option[0]) return;
693
694 TString opt = option;
695 opt.ToUpper();
696
697 // parse firt the specific options
698 if (type == kHistogram) {
699
700 if (opt.Contains("WIDTH")) {
701 fitOption.BinVolume = 1; // scale content by the bin width
702 if (opt.Contains("NORMWIDTH")) {
703 // for variable bins: scale content by the bin width normalized by a reference value (typically the minimum bin)
704 // this option is for variable bin widths
705 fitOption.BinVolume = 2;
706 opt.ReplaceAll("NORMWIDTH","");
707 }
708 else
709 opt.ReplaceAll("WIDTH","");
710 }
711
712 // if (opt.Contains("MULTIPROC")) {
713 // fitOption.ExecPolicy = ROOT::Fit::kMultiprocess;
714 // opt.ReplaceAll("MULTIPROC","");
715 // }
716
717 if (opt.Contains("SERIAL")) {
719 opt.ReplaceAll("SERIAL","");
720 }
721
722 if (opt.Contains("MULTITHREAD")) {
724 opt.ReplaceAll("MULTITHREAD","");
725 }
726
727 if (opt.Contains("I")) fitOption.Integral= 1; // integral of function in the bin (no sense for graph)
728 if (opt.Contains("W")) fitOption.W1 = 1; // all non-empty bins or points have weight =1 (for chi2 fit)
729 if (opt.Contains("WW")) fitOption.W1 = 2; //all bins have weight=1, even empty bins
730 if (opt.Contains("L")) fitOption.Like = 1;
731 if (opt.Contains("X")) fitOption.Chi2 = 1;
732 if (opt.Contains("P")) fitOption.PChi2 = 1;
733
734 // specific likelihood fit options
735 if (fitOption.Like == 1) {
736 //if (opt.Contains("LL")) fitOption.Like = 2;
737 if (opt.Contains("W")){ fitOption.Like = 2; fitOption.W1=0;}// (weighted likelihood)
738 if (opt.Contains("MULTI")) {
739 if (fitOption.Like == 2) fitOption.Like = 6; // weighted multinomial
740 else fitOption.Like = 4; // multinomial likelihood fit instead of Poisson
741 opt.ReplaceAll("MULTI","");
742 }
743 // give precedence for likelihood options
744 if (fitOption.Chi2 == 1 || fitOption.PChi2 == 1)
745 Warning("Fit","Cannot use P or X option in combination of L. Ignore the chi2 option and perform a likelihood fit");
746 }
747 if (fitOption.PChi2 && fitOption.W1) {
748 Warning("FitOptionsMake", "Ignore option W or WW when used together with option P (Pearson chi2)");
749 fitOption.W1 = 0; // with Pearson chi2 W option is ignored
750 }
751 }
752 // specific Graph options (need to be parsed before)
753 else if (type == kGraph) {
754 opt.ReplaceAll("ROB", "H");
755 opt.ReplaceAll("EX0", "T");
756
757 //for robust fitting, see if # of good points is defined
758 // decode parameters for robust fitting
759 Double_t h=0;
760 if (opt.Contains("H=0.")) {
761 int start = opt.Index("H=0.");
762 int numpos = start + strlen("H=0.");
763 int numlen = 0;
764 int len = opt.Length();
765 while( (numpos+numlen<len) && isdigit(opt[numpos+numlen]) ) numlen++;
766 TString num = opt(numpos,numlen);
767 opt.Remove(start+strlen("H"),strlen("=0.")+numlen);
768 h = atof(num.Data());
769 h*=TMath::Power(10, -numlen);
770 }
771
772 if (opt.Contains("H")) { fitOption.Robust = 1; fitOption.hRobust = h; }
773 if (opt.Contains("T")) fitOption.NoErrX = 1; // no error in X
774 if (opt.Contains("W")) fitOption.W1 = 1; // ignorer all point errors when fitting
775 }
776
777 if (opt.Contains("U")) fitOption.User = 1;
778 if (opt.Contains("Q")) fitOption.Quiet = 1;
779 if (opt.Contains("V")) {fitOption.Verbose = 1; fitOption.Quiet = 0;}
780
781
782 if (opt.Contains("E")) fitOption.Errors = 1;
783 if (opt.Contains("R")) fitOption.Range = 1;
784 if (opt.Contains("G")) fitOption.Gradient= 1;
785 if (opt.Contains("M")) fitOption.More = 1;
786 if (opt.Contains("N")) fitOption.Nostore = 1;
787 if (opt.Contains("0")) fitOption.Nograph = 1;
788 if (opt.Contains("+")) fitOption.Plus = 1;
789 if (opt.Contains("B")) fitOption.Bound = 1;
790 if (opt.Contains("C")) fitOption.Nochisq = 1;
791 if (opt.Contains("F")) fitOption.Minuit = 1;
792 if (opt.Contains("S")) fitOption.StoreResult = 1;
793
794}
795
797 if (foption.Like) {
798 Info("CheckGraphFitOptions","L (Log Likelihood fit) is an invalid option when fitting a graph. It is ignored");
799 foption.Like = 0;
800 }
801 if (foption.Integral) {
802 Info("CheckGraphFitOptions","I (use function integral) is an invalid option when fitting a graph. It is ignored");
803 foption.Integral = 0;
804 }
805 return;
806}
807
808// implementation of unbin fit function (defined in HFitInterface)
809
811 // do unbin fit, ownership of fitdata is passed later to the TBackFitter class
812
813 // create a shared pointer to the fit data to managed it
814 std::shared_ptr<ROOT::Fit::UnBinData> fitdata(data);
815
816#ifdef DEBUG
817 printf("tree data size is %d \n",fitdata->Size());
818 for (unsigned int i = 0; i < fitdata->Size(); ++i) {
819 if (fitdata->NDim() == 1) printf(" x[%d] = %f \n", i,*(fitdata->Coords(i) ) );
820 }
821#endif
822 if (fitdata->Size() == 0 ) {
823 Warning("Fit","Fit data is empty ");
824 return -1;
825 }
826
827 // create an empty TFitResult
828 std::shared_ptr<TFitResult> tfr(new TFitResult() );
829 // create the fitter
830 std::shared_ptr<ROOT::Fit::Fitter> fitter(new ROOT::Fit::Fitter(tfr) );
831 ROOT::Fit::FitConfig & fitConfig = fitter->Config();
832
833 // dimension is given by data because TF1 pointer can have wrong one
834 unsigned int dim = fitdata->NDim();
835
836 // set the fit function
837 // if option grad is specified use gradient
838 // need to create a wrapper for an automatic normalized TF1 ???
839 if ( fitOption.Gradient ) {
840 assert ( (int) dim == fitfunc->GetNdim() );
841 fitter->SetFunction(ROOT::Math::WrappedMultiTF1(*fitfunc) );
842 }
843 else
844 fitter->SetFunction(static_cast<const ROOT::Math::IParamMultiFunction &>(ROOT::Math::WrappedMultiTF1(*fitfunc, dim) ) );
845
846 // parameter setting is done automaticaly in the Fitter class
847 // need only to set limits
848 int npar = fitfunc->GetNpar();
849 for (int i = 0; i < npar; ++i) {
850 ROOT::Fit::ParameterSettings & parSettings = fitConfig.ParSettings(i);
851 double plow,pup;
852 fitfunc->GetParLimits(i,plow,pup);
853 // this is a limitation of TF1 interface - cannot fix a parameter to zero value
854 if (plow*pup != 0 && plow >= pup) {
855 parSettings.Fix();
856 }
857 else if (plow < pup ) {
858 if (!TMath::Finite(pup) && TMath::Finite(plow) )
859 parSettings.SetLowerLimit(plow);
860 else if (!TMath::Finite(plow) && TMath::Finite(pup) )
861 parSettings.SetUpperLimit(pup);
862 else
863 parSettings.SetLimits(plow,pup);
864 }
865
866 // set the parameter step size (by default are set to 0.3 of value)
867 // if function provides meaningful error values
868 double err = fitfunc->GetParError(i);
869 if ( err > 0)
870 parSettings.SetStepSize(err);
871 else if (plow < pup && TMath::Finite(plow) && TMath::Finite(pup) ) { // in case of limits improve step sizes
872 double step = 0.1 * (pup - plow);
873 // check if value is not too close to limit otherwise trim value
874 if ( parSettings.Value() < pup && pup - parSettings.Value() < 2 * step )
875 step = (pup - parSettings.Value() ) / 2;
876 else if ( parSettings.Value() > plow && parSettings.Value() - plow < 2 * step )
877 step = (parSettings.Value() - plow ) / 2;
878
879 parSettings.SetStepSize(step);
880 }
881
882 }
883
884 fitConfig.SetMinimizerOptions(minOption);
885
886 if (fitOption.Verbose) fitConfig.MinimizerOptions().SetPrintLevel(3);
887 if (fitOption.Quiet) fitConfig.MinimizerOptions().SetPrintLevel(0);
888
889 // more
890 if (fitOption.More) fitConfig.SetMinimizer("Minuit","MigradImproved");
891
892 // chech if Minos or more options
893 if (fitOption.Errors) {
894 // run Hesse and Minos
895 fitConfig.SetParabErrors(true);
896 fitConfig.SetMinosErrors(true);
897 }
898 // use weight correction
899 if ( (fitOption.Like & 2) == 2)
900 fitConfig.SetWeightCorrection(true);
901
902 bool extended = (fitOption.Like & 1) == 1;
903
904 bool fitok = false;
905 fitok = fitter->LikelihoodFit(fitdata, extended, fitOption.ExecPolicy);
906 if ( !fitok && !fitOption.Quiet )
907 Warning("UnBinFit","Abnormal termination of minimization.");
908
909 const ROOT::Fit::FitResult & fitResult = fitter->Result();
910 // one could set directly the fit result in TF1
911 int iret = fitResult.Status();
912 if (!fitResult.IsEmpty() ) {
913 // set in fitfunc the result of the fit
914 fitfunc->SetNDF(fitResult.Ndf() );
915 fitfunc->SetNumberFitPoints(fitdata->Size() );
916
917 assert( (Int_t)fitResult.Parameters().size() >= fitfunc->GetNpar() );
918 fitfunc->SetParameters( const_cast<double*>(&(fitResult.Parameters().front())));
919 if ( int( fitResult.Errors().size()) >= fitfunc->GetNpar() )
920 fitfunc->SetParErrors( &(fitResult.Errors().front()) );
921
922 }
923
924 // store result in the backward compatible VirtualFitter
925 // in case not running in a multi-thread enabled mode
926 if (gGlobalMutex) {
928 // pass ownership of Fitter and Fitdata to TBackCompFitter (fitter pointer cannot be used afterwards)
929 TBackCompFitter * bcfitter = new TBackCompFitter(fitter, fitdata);
930 // cannot use anymore now fitdata (given away ownership)
931 fitdata = 0;
932 bcfitter->SetFitOption(fitOption);
933 //bcfitter->SetObjectFit(fTree);
934 bcfitter->SetUserFunc(fitfunc);
935
936 if (lastFitter) delete lastFitter;
937 TVirtualFitter::SetFitter( bcfitter );
938
939 // use old-style for printing the results
940 // if (fitOption.Verbose) bcfitter->PrintResults(2,0.);
941 // else if (!fitOption.Quiet) bcfitter->PrintResults(1,0.);
942
943 }
944 // print results
945 if (fitOption.Verbose) fitResult.PrintCovMatrix(std::cout);
946 else if (!fitOption.Quiet) fitResult.Print(std::cout);
947
948 if (fitOption.StoreResult)
949 {
950 TString name = "TFitResult-";
951 name = name + "UnBinData-" + fitfunc->GetName();
952 TString title = "TFitResult-";
953 title += name;
954 tfr->SetName(name);
955 tfr->SetTitle(title);
956 return TFitResultPtr(tfr);
957 }
958 else
959 return TFitResultPtr(iret);
960}
961
962
963// implementations of ROOT::Fit::FitObject functions (defined in HFitInterface) in terms of the template HFit::Fit
964
966moption, const char *goption, ROOT::Fit::DataRange & range) {
967 // check fit options
968 // check if have weights in case of weighted likelihood
969 if ( ((foption.Like & 2) == 2) && h1->GetSumw2N() == 0) {
970 Warning("HFit::FitObject","A weighted likelihood fit is requested but histogram is not weighted - do a standard Likelihood fit");
971 foption.Like = 1;
972 }
973 // histogram fitting
974 return HFit::Fit(h1,f1,foption,moption,goption,range);
975}
976
977TFitResultPtr ROOT::Fit::FitObject(TGraph * gr, TF1 *f1 , Foption_t & foption , const ROOT::Math::MinimizerOptions & moption, const char *goption, ROOT::Fit::DataRange & range) {
978 // exclude options not valid for graphs
980 // TGraph fitting
981 return HFit::Fit(gr,f1,foption,moption,goption,range);
982}
983
984TFitResultPtr ROOT::Fit::FitObject(TMultiGraph * gr, TF1 *f1 , Foption_t & foption , const ROOT::Math::MinimizerOptions & moption, const char *goption, ROOT::Fit::DataRange & range) {
985 // exclude options not valid for graphs
987 // TMultiGraph fitting
988 return HFit::Fit(gr,f1,foption,moption,goption,range);
989}
990
991TFitResultPtr ROOT::Fit::FitObject(TGraph2D * gr, TF1 *f1 , Foption_t & foption , const ROOT::Math::MinimizerOptions & moption, const char *goption, ROOT::Fit::DataRange & range) {
992 // exclude options not valid for graphs
994 // TGraph2D fitting
995 return HFit::Fit(gr,f1,foption,moption,goption,range);
996}
997
998TFitResultPtr ROOT::Fit::FitObject(THnBase * s1, TF1 *f1 , Foption_t & foption , const ROOT::Math::MinimizerOptions & moption, const char *goption, ROOT::Fit::DataRange & range) {
999 // sparse histogram fitting
1000 return HFit::Fit(s1,f1,foption,moption,goption,range);
1001}
1002
1003
1004
1005// Int_t TGraph2D::DoFit(TF2 *f2 ,Option_t *option ,Option_t *goption) {
1006// // internal graph2D fitting methods
1007// Foption_t fitOption;
1008// ROOT::Fit::FitOptionsMake(option,fitOption);
1009
1010// // create range and minimizer options with default values
1011// ROOT::Fit::DataRange range(2);
1012// ROOT::Math::MinimizerOptions minOption;
1013// return ROOT::Fit::FitObject(this, f2 , fitOption , minOption, goption, range);
1014// }
1015
1016
1017// function to compute the simple chi2 for graphs and histograms
1018
1019double ROOT::Fit::Chisquare(const TH1 & h1, TF1 & f1, bool useRange, bool usePL) {
1020 return HFit::ComputeChi2(h1,f1,useRange, usePL);
1021}
1022
1023double ROOT::Fit::Chisquare(const TGraph & g, TF1 & f1, bool useRange) {
1024 return HFit::ComputeChi2(g,f1, useRange, false);
1025}
1026
1027template<class FitObject>
1028double HFit::ComputeChi2(const FitObject & obj, TF1 & f1, bool useRange, bool usePL ) {
1029
1030 // implement using the fitting classes
1032 if (usePL) opt.fUseEmpty=true;
1034 // get range of function
1035 if (useRange) HFit::GetFunctionRange(f1,range);
1036 // fill the data set
1037 ROOT::Fit::BinData data(opt,range);
1038 ROOT::Fit::FillData(data, &obj, &f1);
1039 if (data.Size() == 0 ) {
1040 Warning("Chisquare","data set is empty - return -1");
1041 return -1;
1042 }
1044 if (usePL) {
1045 // use the poisson log-lokelihood (Baker-Cousins chi2)
1046 ROOT::Fit::PoissonLLFunction nll(data, wf1);
1047 return 2.* nll( f1.GetParameters() ) ;
1048 }
1049 ROOT::Fit::Chi2Function chi2(data, wf1);
1050 return chi2(f1.GetParameters() );
1051
1052}
@ kGraph
Definition Buttons.h:34
#define g(i)
Definition RSha256.hxx:105
#define s1(x)
Definition RSha256.hxx:91
#define h(i)
Definition RSha256.hxx:106
static const double x2[5]
static const double x1[5]
const Bool_t kFALSE
Definition RtypesCore.h:101
double Double_t
Definition RtypesCore.h:59
const Bool_t kTRUE
Definition RtypesCore.h:100
const Bool_t kIterBackward
Definition TCollection.h:43
#define R__ASSERT(e)
Definition TError.h:118
void Info(const char *location, const char *msgfmt,...)
Use this function for informational messages.
Definition TError.cxx:220
void Error(const char *location, const char *msgfmt,...)
Use this function in case an error occurred.
Definition TError.cxx:187
void Warning(const char *location, const char *msgfmt,...)
Use this function in warning situations.
Definition TError.cxx:231
char name[80]
Definition TGX11.cxx:110
int type
Definition TGX11.cxx:121
float xmin
float ymin
float xmax
float ymax
@ kCanDelete
Definition TObject.h:369
R__EXTERN TVirtualMutex * gGlobalMutex
#define gPad
Class describing the binned data sets : vectors of x coordinates, y values and optionally error on y ...
Definition BinData.h:52
Chi2FCN class for binnned fits using the least square methods.
Definition Chi2FCN.h:46
class describing the range in the coordinates it supports multiple range in a coordinate.
Definition DataRange.h:35
void AddRange(unsigned int icoord, double xmin, double xmax)
add a range [xmin,xmax] for the new coordinate icoord Adding a range does not delete existing one,...
Definition DataRange.cxx:94
unsigned int Size(unsigned int icoord=0) const
return range size for coordinate icoord (starts from zero) Size == 0 indicates no range is present [-...
Definition DataRange.h:71
void GetRange(unsigned int irange, unsigned int icoord, double &xmin, double &xmax) const
get the i-th range for given coordinate.
Definition DataRange.h:104
Class describing the configuration of the fit, options and parameter settings using the ROOT::Fit::Pa...
Definition FitConfig.h:47
void SetMinosErrors(bool on=true)
set Minos erros computation to be performed after fitting
Definition FitConfig.h:231
void SetNormErrors(bool on=true)
set the option to normalize the error on the result according to chi2/ndf
Definition FitConfig.h:225
void SetMinimizer(const char *type, const char *algo=0)
set minimizer type
Definition FitConfig.h:181
bool NormalizeErrors() const
flag to check if resulting errors are be normalized according to chi2/ndf
Definition FitConfig.h:204
void SetMinimizerOptions(const ROOT::Math::MinimizerOptions &minopt)
set all the minimizer options using class MinimizerOptions
void SetWeightCorrection(bool on=true)
apply the weight correction for error matric computation
Definition FitConfig.h:234
void SetParabErrors(bool on=true)
set parabolic erros
Definition FitConfig.h:228
const std::string & MinimizerType() const
return type of minimizer package
Definition FitConfig.h:189
const ParameterSettings & ParSettings(unsigned int i) const
get the parameter settings for the i-th parameter (const method)
Definition FitConfig.h:76
ROOT::Math::MinimizerOptions & MinimizerOptions()
access to the minimizer control parameter (non const method)
Definition FitConfig.h:167
unsigned int Size() const
return number of fit points
Definition FitData.h:303
class containg the result of the fit and all the related information (fitted parameter values,...
Definition FitResult.h:47
bool IsEmpty() const
True if a fit result does not exist (even invalid) with parameter values.
Definition FitResult.h:117
const std::vector< double > & Errors() const
parameter errors (return st::vector)
Definition FitResult.h:169
const std::vector< double > & Parameters() const
parameter values (return std::vector)
Definition FitResult.h:174
unsigned int Ndf() const
Number of degree of freedom.
Definition FitResult.h:163
double Chi2() const
Chi2 fit value in case of likelihood must be computed ?
Definition FitResult.h:160
void Print(std::ostream &os, bool covmat=false) const
print the result and optionaly covariance matrix and correlations
void PrintCovMatrix(std::ostream &os) const
print error matrix and correlations
int Status() const
minimizer status code
Definition FitResult.h:137
Fitter class, entry point for performing all type of fits.
Definition Fitter.h:77
Class, describing value, limits and step size of the parameters Provides functionality also to set/re...
void SetStepSize(double err)
set the step size
void SetLimits(double low, double up)
set a double side limit, if low == up the parameter is fixed if low > up the limits are removed The c...
double Value() const
copy constructor and assignment operators (leave them to the compiler)
void SetUpperLimit(double up)
set a single upper limit
void Fix()
fix the parameter
void SetLowerLimit(double low)
set a single lower limit
class evaluating the log likelihood for binned Poisson likelihood fits it is template to distinguish ...
Class describing the unbinned data sets (just x coordinates values) of any dimensions.
Definition UnBinData.h:42
IParamFunction interface (abstract class) describing multi-dimensional parameteric functions It is a ...
void SetPrintLevel(int level)
set print level
void SetTolerance(double tol)
set the tolerance
Class to Wrap a ROOT Function class (like TF1) in a IParamMultiFunction interface of multi-dimensions...
Class to manage histogram axis.
Definition TAxis.h:30
Double_t GetXmax() const
Definition TAxis.h:134
virtual Double_t GetBinLowEdge(Int_t bin) const
Return low edge of bin.
Definition TAxis.cxx:518
Int_t GetLast() const
Return last bin on the axis i.e.
Definition TAxis.cxx:469
Double_t GetXmin() const
Definition TAxis.h:133
virtual Double_t GetBinWidth(Int_t bin) const
Return bin width.
Definition TAxis.cxx:540
Int_t GetFirst() const
Return first bin on the axis i.e.
Definition TAxis.cxx:458
Backward compatible implementation of TVirtualFitter.
virtual void SetMethodCall(TMethodCall *m)
For using interpreted function passed by the user.
virtual void SetFCN(void(*fcn)(Int_t &, Double_t *, Double_t &f, Double_t *, Int_t))
Override setFCN to use the Adapter to Minuit2 FCN interface To set the address of the minimization fu...
1-Dim function class
Definition TF1.h:213
virtual Int_t GetNumber() const
Definition TF1.h:498
virtual void GetParLimits(Int_t ipar, Double_t &parmin, Double_t &parmax) const
Return limits for parameter ipar.
Definition TF1.cxx:1936
virtual void SetNDF(Int_t ndf)
Set the number of degrees of freedom ndf should be the number of points used in a fit - the number of...
Definition TF1.cxx:3439
virtual Double_t GetParError(Int_t ipar) const
Return value of parameter number ipar.
Definition TF1.cxx:1926
virtual void SetChisquare(Double_t chi2)
Definition TF1.h:612
virtual void Copy(TObject &f1) const
Copy this F1 to a new F1.
Definition TF1.cxx:1000
virtual void SetRange(Double_t xmin, Double_t xmax)
Initialize the upper and lower bounds to draw the function.
Definition TF1.cxx:3539
virtual void SetParent(TObject *p=0)
Definition TF1.h:665
virtual Int_t GetNpar() const
Definition TF1.h:481
virtual void SetParErrors(const Double_t *errors)
Set errors for all active parameters when calling this function, the array errors must have at least ...
Definition TF1.cxx:3504
virtual Double_t * GetParameters() const
Definition TF1.h:520
virtual void SetNumberFitPoints(Int_t npfits)
Definition TF1.h:624
@ kNotDraw
Definition TF1.h:326
virtual void GetRange(Double_t *xmin, Double_t *xmax) const
Return range of a generic N-D function.
Definition TF1.cxx:2275
virtual Bool_t IsLinear() const
Definition TF1.h:602
bool IsVectorized()
Definition TF1.h:440
virtual void SetParameters(const Double_t *params)
Definition TF1.h:644
virtual void Save(Double_t xmin, Double_t xmax, Double_t ymin, Double_t ymax, Double_t zmin, Double_t zmax)
Save values of function in array fSave.
Definition TF1.cxx:3154
virtual Int_t GetNdim() const
Definition TF1.h:485
virtual Bool_t AddToGlobalList(Bool_t on=kTRUE)
Add to global list of functions (gROOT->GetListOfFunctions() ) return previous status (true if the fu...
Definition TF1.cxx:841
A 2-Dim function with parameters.
Definition TF2.h:29
virtual void Save(Double_t xmin, Double_t xmax, Double_t ymin, Double_t ymax, Double_t zmin, Double_t zmax)
Save values of function in array fSave.
Definition TF2.cxx:799
virtual void SetRange(Double_t xmin, Double_t xmax)
Initialize the upper and lower bounds to draw the function.
Definition TF2.h:148
A 3-Dim function with parameters.
Definition TF3.h:28
virtual void SetRange(Double_t xmin, Double_t xmax)
Initialize the upper and lower bounds to draw the function.
Definition TF3.h:145
virtual void Save(Double_t xmin, Double_t xmax, Double_t ymin, Double_t ymax, Double_t zmin, Double_t zmax)
Save values of function in array fSave.
Definition TF3.cxx:535
Provides an indirection to the TFitResult class and with a semantics identical to a TFitResult pointe...
Extends the ROOT::Fit::Result class with a TNamed inheritance providing easy possibility for I/O.
Definition TFitResult.h:34
Graphics object made of three arrays X, Y and Z with the same number of points each.
Definition TGraph2D.h:41
A TGraph is an object made of two arrays X and Y with npoints each.
Definition TGraph.h:41
virtual void ComputeRange(Double_t &xmin, Double_t &ymin, Double_t &xmax, Double_t &ymax) const
Compute the x/y range of the points in this graph.
Definition TGraph.cxx:669
TH1F * GetHistogram() const
Returns a pointer to the histogram used to draw the axis Takes into account the two following cases.
Definition TGraph.cxx:1485
TH1 is the base class of all histogram classes in ROOT.
Definition TH1.h:58
TAxis * GetZaxis()
Definition TH1.h:322
virtual Int_t GetDimension() const
Definition TH1.h:282
TAxis * GetXaxis()
Get the behaviour adopted by the object about the statoverflows. See EStatOverflows for more informat...
Definition TH1.h:320
TAxis * GetYaxis()
Definition TH1.h:321
TList * GetListOfFunctions() const
Definition TH1.h:243
virtual void Draw(Option_t *option="")
Draw this histogram with options.
Definition TH1.cxx:3074
virtual Int_t GetSumw2N() const
Definition TH1.h:314
Multidimensional histogram base.
Definition THnBase.h:43
A doubly linked list.
Definition TList.h:38
virtual void Add(TObject *obj)
Definition TList.h:81
virtual TObject * Remove(TObject *obj)
Remove object from the list.
Definition TList.cxx:822
A TMultiGraph is a collection of TGraph (or derived) objects.
Definition TMultiGraph.h:36
TList * GetListOfGraphs() const
Definition TMultiGraph.h:70
TH1F * GetHistogram()
Returns a pointer to the histogram used to draw the axis.
virtual const char * GetTitle() const
Returns title of object.
Definition TNamed.h:48
virtual const char * GetName() const
Returns name of object.
Definition TNamed.h:47
Mother of all ROOT objects.
Definition TObject.h:41
R__ALWAYS_INLINE Bool_t TestBit(UInt_t f) const
Definition TObject.h:201
R__ALWAYS_INLINE Bool_t IsZombie() const
Definition TObject.h:153
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition TObject.cxx:766
virtual Bool_t InheritsFrom(const char *classname) const
Returns kTRUE if object inherits from class "classname".
Definition TObject.cxx:515
Basic string class.
Definition TString.h:136
Ssiz_t Length() const
Definition TString.h:410
const char * Data() const
Definition TString.h:369
TString & ReplaceAll(const TString &s1, const TString &s2)
Definition TString.h:692
void ToUpper()
Change string to upper case.
Definition TString.cxx:1163
TString & Remove(Ssiz_t pos)
Definition TString.h:673
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition TString.h:624
Ssiz_t Index(const char *pat, Ssiz_t i=0, ECaseCompare cmp=kExact) const
Definition TString.h:639
Abstract Base Class for Fitting.
virtual void SetFitOption(Foption_t option)
virtual void SetObjectFit(TObject *obj)
TMethodCall * GetMethodCall() const
void(* FCNFunc_t)(Int_t &npar, Double_t *gin, Double_t &f, Double_t *u, Int_t flag)
virtual void SetUserFunc(TObject *userfunc)
static TVirtualFitter * GetFitter()
static: return the current Fitter
static void SetFitter(TVirtualFitter *fitter, Int_t maxpar=25)
Static function to set an alternative fitter.
TGraphErrors * gr
Definition legend1.C:25
TH1F * h1
Definition legend1.C:5
TF1 * f1
Definition legend1.C:11
void GetDrawingRange(TH1 *h1, ROOT::Fit::DataRange &range)
Definition HFitImpl.cxx:461
void GetFunctionRange(const TF1 &f1, ROOT::Fit::DataRange &range)
Definition HFitImpl.cxx:120
int CheckFitFunction(const TF1 *f1, int hdim)
Definition HFitImpl.cxx:86
TFitResultPtr Fit(FitObject *h1, TF1 *f1, Foption_t &option, const ROOT::Math::MinimizerOptions &moption, const char *goption, ROOT::Fit::DataRange &range)
Definition HFitImpl.cxx:133
void FitOptionsMake(const char *option, Foption_t &fitOption)
void CheckGraphFitOptions(Foption_t &fitOption)
Definition HFitImpl.cxx:796
void StoreAndDrawFitFunction(FitObject *h1, TF1 *f1, const ROOT::Fit::DataRange &range, bool, bool, const char *goption)
Definition HFitImpl.cxx:569
double ComputeChi2(const FitObject &h1, TF1 &f1, bool useRange, bool usePL)
int GetDimension(const TH1 *h1)
Definition HFitImpl.cxx:50
TFitResultPtr FitObject(TH1 *h1, TF1 *f1, Foption_t &option, const ROOT::Math::MinimizerOptions &moption, const char *goption, ROOT::Fit::DataRange &range)
fitting function for a TH1 (called from TH1::Fit)
Definition HFitImpl.cxx:965
void FitOptionsMake(EFitObjectType type, const char *option, Foption_t &fitOption)
Decode list of options into fitOption.
Definition HFitImpl.cxx:684
void Init2DGaus(const ROOT::Fit::BinData &data, TF1 *f1)
compute initial parameter for 2D gaussian function given the fit data Set the sigma limits for zero t...
TFitResultPtr UnBinFit(ROOT::Fit::UnBinData *data, TF1 *f1, Foption_t &option, const ROOT::Math::MinimizerOptions &moption)
fit an unbin data set (from tree or from histogram buffer) using a TF1 pointer and fit options.
Definition HFitImpl.cxx:810
void FillData(BinData &dv, const TH1 *hist, TF1 *func=0)
fill the data vector from a TH1.
void InitExpo(const ROOT::Fit::BinData &data, TF1 *f1)
compute initial parameter for an exponential function given the fit data Set the constant and slope a...
void InitGaus(const ROOT::Fit::BinData &data, TF1 *f1)
compute initial parameter for gaussian function given the fit data Set the sigma limits for zero top ...
double Chisquare(const TH1 &h1, TF1 &f1, bool useRange, bool usePL=false)
compute the chi2 value for an histogram given a function (see TH1::Chisquare for the documentation)
std::string ToString(const T &val)
Utility function for conversion to strings.
Definition Util.h:50
Bool_t IsImplicitMTEnabled()
Returns true if the implicit multi-threading in ROOT is enabled.
Definition TROOT.cxx:558
Int_t Finite(Double_t x)
Check if it is finite with a mask in order to be consistent in presence of fast math.
Definition TMath.h:721
LongDouble_t Power(LongDouble_t x, LongDouble_t y)
Definition TMath.h:685
int Range
Definition Foption.h:39
int Nograph
Definition Foption.h:42
int Quiet
Definition Foption.h:29
int Like
Definition Foption.h:34
int W1
Definition Foption.h:36
int Gradient
Definition Foption.h:40
ROOT::EExecutionPolicy ExecPolicy
Definition Foption.h:52
int StoreResult
Definition Foption.h:49
int Nochisq
Definition Foption.h:45
int Robust
Definition Foption.h:48
double hRobust
Definition Foption.h:51
int Plus
Definition Foption.h:43
int Integral
Definition Foption.h:44
int Bound
Definition Foption.h:31
int Nostore
Definition Foption.h:41
int More
Definition Foption.h:38
int PChi2
Definition Foption.h:33
int Chi2
Definition Foption.h:32
int Minuit
Definition Foption.h:46
int Errors
Definition Foption.h:37
int NoErrX
Definition Foption.h:47
int Verbose
Definition Foption.h:30
int User
Definition Foption.h:35
int BinVolume
Definition Foption.h:50
DataOptions : simple structure holding the options on how the data are filled.
Definition DataOptions.h:28
bool fErrors1
use the function range when creating the fit data (default is false)
Definition DataOptions.h:52
bool fNormBinVolume
normalize data by the bin volume (it is used in the Poisson likelihood fits)
Definition DataOptions.h:49
bool fUseRange
use empty bins (default is false) with a fixed error of 1
Definition DataOptions.h:51
bool fUseEmpty
normalize data by a normalized the bin volume (bin volume divided by a reference value)
Definition DataOptions.h:50
bool fExpErrors
use all errors equal to 1, i.e. fit without errors (default is false)
Definition DataOptions.h:53
bool fBinVolume
use integral of bin content instead of bin center (default is false)
Definition DataOptions.h:48
bool fCoordErrors
use expected errors from the function and not from the data
Definition DataOptions.h:54