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 <cstdlib>
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, ROOT::Fit::EChisquareType type );
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 if (fitOption.PChi2 == 1) {
182 opt.fErrors1 = true; // we are not using errors in chi2, it is like setting = 1
183 } else if (fitOption.PChi2 == 2) {
184 opt.fErrors1 = false; // we need the errors in weighted likelihood fit
185 }
186
187 if (fitOption.BinVolume) {
188 opt.fBinVolume = true; // scale by bin volume
189 if (fitOption.BinVolume == 2) opt.fNormBinVolume = true; // scale by normalized bin volume
190 }
191
192 if (opt.fUseRange) {
193#ifdef DEBUG
194 printf("use range \n" );
195#endif
197 }
198#ifdef DEBUG
199 printf("range size %d\n", range.Size(0) );
200 if (range.Size(0)) {
201 double x1; double x2; range.GetRange(0,x1,x2);
202 printf(" range in x = [%f,%f] \n",x1,x2);
203 }
204#endif
205
206 // fill data
207 std::shared_ptr<ROOT::Fit::BinData> fitdata(new ROOT::Fit::BinData(opt,range) );
208 ROOT::Fit::FillData(*fitdata, h1, f1);
209 if (fitdata->Size() == 0 ) {
210 Warning("Fit","Fit data is empty ");
211 return -1;
212 }
213
214#ifdef DEBUG
215 printf("HFit:: data size is %d \n",fitdata->Size());
216 for (unsigned int i = 0; i < fitdata->Size(); ++i) {
217 if (fitdata->NDim() == 1) printf(" x[%d] = %f - value = %f \n", i,*(fitdata->Coords(i)),fitdata->Value(i) );
218 }
219#endif
220
221 // switch off linear fitting in case data has coordinate errors and the option is set
222 if (fitdata->GetErrorType() == ROOT::Fit::BinData::kCoordError && fitdata->Opt().fCoordErrors ) linear = false;
223 // linear fit cannot be done also in case of asymmetric errors
224 if (fitdata->GetErrorType() == ROOT::Fit::BinData::kAsymError && fitdata->Opt().fAsymErrors ) linear = false;
225
226 // this functions use the TVirtualFitter
227 if (special != 0 && !fitOption.Bound && !linear) {
228 if (special == 100) ROOT::Fit::InitGaus (*fitdata,f1); // gaussian
229 else if (special == 110 || special == 112) ROOT::Fit::Init2DGaus(*fitdata,f1); // 2D gaussians ( xygaus or bigaus)
230 else if (special == 400) ROOT::Fit::InitGaus (*fitdata,f1); // landau (use the same)
231 else if (special == 410) ROOT::Fit::Init2DGaus(*fitdata,f1); // 2D landau (use the same)
232
233 else if (special == 200) ROOT::Fit::InitExpo (*fitdata, f1); // exponential
234
235 }
236
237
238 // set the fit function
239 // if option grad is specified use gradient
240 if ( (linear || fitOption.Gradient) )
241 fitter->SetFunction(ROOT::Math::WrappedMultiTF1(*f1));
242#ifdef R__HAS_VECCORE
243 else if(f1->IsVectorized())
245#endif
246 else
247 fitter->SetFunction(static_cast<const ROOT::Math::IParamMultiFunction &>(ROOT::Math::WrappedMultiTF1(*f1) ) );
248
249 // error normalization in case of zero error in the data
250 if (fitdata->GetErrorType() == ROOT::Fit::BinData::kNoError) fitConfig.SetNormErrors(true);
251 // error normalization also in case of W or WW options (weights = 1)
252 if (fitdata->Opt().fErrors1) fitConfig.SetNormErrors(true);
253 // normalize errors also in case you are fitting a Ndim histo with a N-1 function
254 if (int(fitdata->NDim()) == hdim -1 ) fitConfig.SetNormErrors(true);
255
256
257 // here need to get some static extra information (like max iterations, error def, etc...)
258
259
260 // parameter settings and transfer the parameters values, names and limits from the functions
261 // is done automatically in the Fitter.cxx
262 for (int i = 0; i < npar; ++i) {
263 ROOT::Fit::ParameterSettings & parSettings = fitConfig.ParSettings(i);
264
265 // check limits
266 double plow,pup;
267 f1->GetParLimits(i,plow,pup);
268 if (plow*pup != 0 && plow >= pup) { // this is a limitation - cannot fix a parameter to zero value
269 parSettings.Fix();
270 }
271 else if (plow < pup ) {
272 if (!TMath::Finite(pup) && TMath::Finite(plow) )
273 parSettings.SetLowerLimit(plow);
274 else if (!TMath::Finite(plow) && TMath::Finite(pup) )
275 parSettings.SetUpperLimit(pup);
276 else
277 parSettings.SetLimits(plow,pup);
278 }
279
280 // set the parameter step size (by default are set to 0.3 of value)
281 // if function provides meaningful error values
282 double err = f1->GetParError(i);
283 if ( err > 0)
284 parSettings.SetStepSize(err);
285 else if (plow < pup && TMath::Finite(plow) && TMath::Finite(pup) ) { // in case of limits improve step sizes
286 double step = 0.1 * (pup - plow);
287 // check if value is not too close to limit otherwise trim value
288 if ( parSettings.Value() < pup && pup - parSettings.Value() < 2 * step )
289 step = (pup - parSettings.Value() ) / 2;
290 else if ( parSettings.Value() > plow && parSettings.Value() - plow < 2 * step )
291 step = (parSettings.Value() - plow ) / 2;
292
293 parSettings.SetStepSize(step);
294 }
295
296
297 }
298
299 // needed for setting precision ?
300 // - Compute sum of squares of errors in the bin range
301 // should maybe use stat[1] ??
302 // Double_t ey, sumw2=0;
303// for (i=hxfirst;i<=hxlast;i++) {
304// ey = GetBinError(i);
305// sumw2 += ey*ey;
306// }
307
308
309 // set all default minimizer options (tolerance, max iterations, etc..)
310 fitConfig.SetMinimizerOptions(minOption);
311
312 // specific print level options
313 if (fitOption.Verbose) fitConfig.MinimizerOptions().SetPrintLevel(3);
314 if (fitOption.Quiet) fitConfig.MinimizerOptions().SetPrintLevel(0);
315
316 // specific minimizer options depending on minimizer
317 if (linear) {
318 if (fitOption.Robust ) {
319 // robust fitting
320 std::string type = "Robust";
321 // if an h is specified print out the value adding to the type
322 if (fitOption.hRobust > 0 && fitOption.hRobust < 1.)
323 type += " (h=" + ROOT::Math::Util::ToString(fitOption.hRobust) + ")";
324 fitConfig.SetMinimizer("Linear",type.c_str());
325 fitConfig.MinimizerOptions().SetTolerance(fitOption.hRobust); // use tolerance for passing robust parameter
326 }
327 else
328 fitConfig.SetMinimizer("Linear","");
329 }
330 else {
331 if (fitOption.More) fitConfig.SetMinimizer("Minuit","MigradImproved");
332 }
333
334
335 // check if Error option (run Hesse and Minos) then
336 if (fitOption.Errors) {
337 // run Hesse and Minos
338 fitConfig.SetParabErrors(true);
339 fitConfig.SetMinosErrors(true);
340 }
341
342
343 // do fitting
344
345#ifdef DEBUG
346 if (fitOption.Like) printf("do likelihood fit...\n");
347 if (linear) printf("do linear fit...\n");
348 printf("do now fit...\n");
349#endif
350
351 bool fitok = false;
352
353
354 // check if can use option user
355 TVirtualFitter::FCNFunc_t userFcn = nullptr;
356 // option user is enabled only when running in serial mode
357 if (fitOption.User && TVirtualFitter::GetFitter() ) {
358 userFcn = (TVirtualFitter::GetFitter())->GetFCN();
359 (TVirtualFitter::GetFitter())->SetUserFunc(f1);
360 }
361
362
363 if (fitOption.User && userFcn) // user provided fit objective function
364 fitok = fitter->FitFCN( userFcn );
365 else if (fitOption.Like) {// likelihood fit
366 // perform a weighted likelihood fit by applying weight correction to errors
367 bool weight = ((fitOption.Like & 2) == 2);
368 fitConfig.SetWeightCorrection(weight);
369 bool extended = ((fitOption.Like & 4 ) != 4 );
370 //if (!extended) Info("HFitImpl","Do a not -extended binned fit");
371
372 // pass fitdata as a shared pointer so ownership is shared with Fitter and FitResult class
373 fitok = fitter->LikelihoodFit(fitdata, extended, fitOption.ExecPolicy);
374 }
375 else{ // standard least square fit
376 fitok = fitter->Fit(fitdata, fitOption.ExecPolicy);
377 }
378 if ( !fitok && !fitOption.Quiet )
379 Warning("Fit","Abnormal termination of minimization.");
380 iret |= !fitok;
381
382
383 const ROOT::Fit::FitResult & fitResult = fitter->Result();
384 // one could set directly the fit result in TF1
385 iret = fitResult.Status();
386 if (!fitResult.IsEmpty() ) {
387 // set in f1 the result of the fit
388 f1->SetChisquare(fitResult.Chi2() );
389 f1->SetNDF(fitResult.Ndf() );
390 f1->SetNumberFitPoints(fitdata->Size() );
391
392 assert((Int_t)fitResult.Parameters().size() >= f1->GetNpar() );
393 f1->SetParameters( const_cast<double*>(&(fitResult.Parameters().front())));
394 if ( int( fitResult.Errors().size()) >= f1->GetNpar() )
395 f1->SetParErrors( &(fitResult.Errors().front()) );
396
397
398 }
399
400// - Store fitted function in histogram functions list and draw
401 if (!fitOption.Nostore) {
403 HFit::StoreAndDrawFitFunction(h1, f1, range, !fitOption.Plus, !fitOption.Nograph, goption);
404 }
405
406 // print the result
407 // if using Fitter class must be done here
408 // use old style Minuit for TMinuit and if no corrections have been applied
409 if (!fitOption.Quiet) {
411 if (fitter->GetMinimizer() && fitConfig.MinimizerType() == "Minuit" &&
412 !fitConfig.NormalizeErrors() && fitOption.Like <= 1) {
413 fitter->GetMinimizer()->PrintResults(); // use old style Minuit
414 }
415 else {
416 // print using FitResult class
417 if (fitOption.Verbose) fitResult.PrintCovMatrix(std::cout);
418 fitResult.Print(std::cout);
419 }
420 }
421
422 // store result in the backward compatible TVirtualFitter
423 {
426 TBackCompFitter * bcfitter = new TBackCompFitter(fitter, fitdata);
427 bcfitter->SetFitOption(fitOption);
428 bcfitter->SetObjectFit(h1);
429 bcfitter->SetUserFunc(f1);
431 if (userFcn) {
432 bcfitter->SetFCN(userFcn);
433 // for interpreted FCN functions
434 if (lastFitter->GetMethodCall() ) bcfitter->SetMethodCall(lastFitter->GetMethodCall() );
435 }
436
437 // delete last fitter if it has been created here before
438 if (lastFitter) {
439 TBackCompFitter * lastBCFitter = dynamic_cast<TBackCompFitter *> (lastFitter);
440 if (lastBCFitter && lastBCFitter->TestBit(TBackCompFitter::kCanDeleteLast) )
441 delete lastBCFitter;
442 }
443 //N.B= this might create a memory leak if user does not delete the fitter they create
444 TVirtualFitter::SetFitter( bcfitter );
445 }
446
447 if (fitOption.StoreResult)
448 {
449 TString name = "TFitResult-";
450 name = name + h1->GetName() + "-" + f1->GetName();
451 TString title = "TFitResult-";
452 title += h1->GetTitle();
453 tfr->SetName(name);
454 tfr->SetTitle(title);
455 return TFitResultPtr(tfr);
456 }
457 else
458 return TFitResultPtr(iret);
459}
460
461
463 // get range from histogram and update the DataRange class
464 // if a ranges already exist in that dimension use that one
465
466 Int_t ndim = GetDimension(h1);
467
468 double xmin = 0, xmax = 0, ymin = 0, ymax = 0, zmin = 0, zmax = 0;
469 if (range.Size(0) == 0) {
470 TAxis & xaxis = *(h1->GetXaxis());
471 Int_t hxfirst = xaxis.GetFirst();
472 Int_t hxlast = xaxis.GetLast();
473 Double_t binwidx = xaxis.GetBinWidth(hxlast);
474 xmin = xaxis.GetBinLowEdge(hxfirst);
475 xmax = xaxis.GetBinLowEdge(hxlast) +binwidx;
476 range.AddRange(xmin,xmax);
477 }
478
479 if (ndim > 1) {
480 if (range.Size(1) == 0) {
481 TAxis & yaxis = *(h1->GetYaxis());
482 Int_t hyfirst = yaxis.GetFirst();
483 Int_t hylast = yaxis.GetLast();
484 Double_t binwidy = yaxis.GetBinWidth(hylast);
485 ymin = yaxis.GetBinLowEdge(hyfirst);
486 ymax = yaxis.GetBinLowEdge(hylast) +binwidy;
487 range.AddRange(1,ymin,ymax);
488 }
489 }
490 if (ndim > 2) {
491 if (range.Size(2) == 0) {
492 TAxis & zaxis = *(h1->GetZaxis());
493 Int_t hzfirst = zaxis.GetFirst();
494 Int_t hzlast = zaxis.GetLast();
495 Double_t binwidz = zaxis.GetBinWidth(hzlast);
496 zmin = zaxis.GetBinLowEdge(hzfirst);
497 zmax = zaxis.GetBinLowEdge(hzlast) +binwidz;
498 range.AddRange(2,zmin,zmax);
499 }
500 }
501#ifdef DEBUG
502 std::cout << "xmin,xmax" << xmin << " " << xmax << std::endl;
503#endif
504
505}
506
508 // get range for graph (used sub-set histogram)
509 // N.B. : this is different than in previous implementation of TGraph::Fit where range used was from xmin to xmax.
510 TH1 * h1 = gr->GetHistogram();
511 // an histogram is normally always returned for a TGraph
512 if (h1) HFit::GetDrawingRange(h1, range);
513}
515 // get range for multi-graph (used sub-set histogram)
516 // N.B. : this is different than in previous implementation of TMultiGraph::Fit where range used was from data xmin to xmax.
517 TH1 * h1 = mg->GetHistogram();
518 if (h1) {
520 }
521 else if (range.Size(0) == 0) {
522 // compute range from all the TGraph's belonging to the MultiGraph
523 double xmin = std::numeric_limits<double>::infinity();
524 double xmax = -std::numeric_limits<double>::infinity();
525 TIter next(mg->GetListOfGraphs() );
526 TGraph * g = nullptr;
527 while ( (g = (TGraph*) next() ) ) {
528 double x1 = 0, x2 = 0, y1 = 0, y2 = 0;
529 g->ComputeRange(x1,y1,x2,y2);
530 if (x1 < xmin) xmin = x1;
531 if (x2 > xmax) xmax = x2;
532 }
533 range.AddRange(xmin,xmax);
534 }
535}
537 // get range for graph2D (used sub-set histogram)
538 // N.B. : this is different than in previous implementation of TGraph2D::Fit. There range used was always(0,0)
539 // cannot use TGraph2D::GetHistogram which makes an interpolation
540 //TH1 * h1 = gr->GetHistogram();
541 //if (h1) HFit::GetDrawingRange(h1, range);
542 // not very efficient (t.b.i.)
543 if (range.Size(0) == 0) {
544 double xmin = gr->GetXmin();
545 double xmax = gr->GetXmax();
546 range.AddRange(0,xmin,xmax);
547 }
548 if (range.Size(1) == 0) {
549 double ymin = gr->GetYmin();
550 double ymax = gr->GetYmax();
551 range.AddRange(1,ymin,ymax);
552 }
553}
554
556 // get range from histogram and update the DataRange class
557 // if a ranges already exist in that dimension use that one
558
559 Int_t ndim = GetDimension(s1);
560
561 for ( int i = 0; i < ndim; ++i ) {
562 if ( range.Size(i) == 0 ) {
563 TAxis *axis = s1->GetAxis(i);
564 range.AddRange(i, axis->GetXmin(), axis->GetXmax());
565 }
566 }
567}
568
569template<class FitObject>
570void HFit::StoreAndDrawFitFunction(FitObject * h1, TF1 * f1, const ROOT::Fit::DataRange & range, bool delOldFunction, bool drawFunction, const char *goption) {
571// - Store fitted function in histogram functions list and draw
572// should have separate functions for 1,2,3d ? t.b.d in case
573
574#ifdef DEBUG
575 std::cout <<"draw and store fit function " << f1->GetName() << std::endl;
576#endif
577
578
579 Int_t ndim = GetDimension(h1);
580 double xmin = 0, xmax = 0, ymin = 0, ymax = 0, zmin = 0, zmax = 0;
581 if (range.Size(0) ) range.GetRange(0,xmin,xmax);
582 if (range.Size(1) ) range.GetRange(1,ymin,ymax);
583 if (range.Size(2) ) range.GetRange(2,zmin,zmax);
584
585
586#ifdef DEBUG
587 std::cout <<"draw and store fit function " << f1->GetName()
588 << " Range in x = [ " << xmin << " , " << xmax << " ]" << std::endl;
589#endif
590
591 TList * funcList = h1->GetListOfFunctions();
592 if (funcList == nullptr){
593 Error("StoreAndDrawFitFunction","Function list has not been created - cannot store the fitted function");
594 return;
595 }
596
597 // delete the function in the list only if
598 // the function we are fitting is not in that list
599 // If this is the case we re-use that function object and
600 // we do not create a new one (if delOldFunction is true)
601 bool reuseOldFunction = false;
602 if (delOldFunction) {
603 TIter next(funcList, kIterBackward);
604 TObject *obj;
605 while ((obj = next())) {
606 if (obj->InheritsFrom(TF1::Class())) {
607 if (obj != f1) {
608 funcList->Remove(obj);
609 delete obj;
610 }
611 else {
612 reuseOldFunction = true;
613 }
614 }
615 }
616 }
617
618 TF1 *fnew1 = nullptr;
619 TF2 *fnew2 = nullptr;
620 TF3 *fnew3 = nullptr;
621
622 // copy TF1 using TClass to avoid slicing in case of derived classes
623 if (ndim < 2) {
624 if (!reuseOldFunction) {
625 fnew1 = (TF1*)f1->IsA()->New();
626 R__ASSERT(fnew1);
627 f1->Copy(*fnew1);
628 funcList->Add(fnew1);
629 }
630 else {
631 fnew1 = f1;
632 }
633 fnew1->SetParent( h1 );
634 fnew1->SetRange(xmin,xmax);
635 fnew1->Save(xmin,xmax,0,0,0,0);
636 if (!drawFunction) fnew1->SetBit(TF1::kNotDraw);
637 fnew1->AddToGlobalList(false);
638 } else if (ndim < 3) {
639 if (!reuseOldFunction) {
640 fnew2 = (TF2*)f1->IsA()->New();
641 R__ASSERT(fnew2);
642 f1->Copy(*fnew2);
643 funcList->Add(fnew2);
644 }
645 else {
646 fnew2 = dynamic_cast<TF2*>(f1);
647 R__ASSERT(fnew2);
648 }
649 fnew2->SetRange(xmin,ymin,xmax,ymax);
650 fnew2->SetParent( h1 );
651 fnew2->Save(xmin,xmax,ymin,ymax,0,0);
652 if (!drawFunction) fnew2->SetBit(TF1::kNotDraw);
653 fnew2->AddToGlobalList(false);
654 } else {
655 if (!reuseOldFunction) {
656 fnew3 = (TF3*)f1->IsA()->New();
657 R__ASSERT(fnew3);
658 f1->Copy(*fnew3);
659 funcList->Add(fnew3);
660 }
661 else {
662 fnew3 = dynamic_cast<TF3*>(f1);
663 R__ASSERT(fnew3);
664 }
665 fnew3->SetRange(xmin,ymin,zmin,xmax,ymax,zmax);
666 fnew3->SetParent( h1 );
667 fnew3->Save(xmin,xmax,ymin,ymax,zmin,zmax);
668 if (!drawFunction) fnew3->SetBit(TF1::kNotDraw);
669 fnew3->AddToGlobalList(false);
670 }
671 if (h1->TestBit(kCanDelete)) return;
672 // draw only in case of histograms
673 if (drawFunction && ndim < 3 && h1->InheritsFrom(TH1::Class() ) ) {
674 // no need to re-draw the histogram if the histogram is already in the pad
675 // in that case the function will be just drawn (if option N is not set)
676 if (!gPad || (gPad && gPad->GetListOfPrimitives()->FindObject(h1) == nullptr ) )
677 h1->Draw(goption);
678 }
679 if (gPad) gPad->Modified(); // this is not in TH1 code (needed ??)
680
681 return;
682}
683
684
686 // - Decode list of options into fitOption (used by both TGraph and TH1)
687 // works for both histograms and graph depending on the enum FitObjectType defined in HFit
690 }
691
692 if (option == nullptr) return;
693 if (!option[0]) return;
694
695 TString opt = option;
696 opt.ToUpper();
697
698 // parse firt the specific options
699 if (type == EFitObjectType::kHistogram) {
700
701 if (opt.Contains("WIDTH")) {
702 fitOption.BinVolume = 1; // scale content by the bin width
703 if (opt.Contains("NORMWIDTH")) {
704 // for variable bins: scale content by the bin width normalized by a reference value (typically the minimum bin)
705 // this option is for variable bin widths
706 fitOption.BinVolume = 2;
707 opt.ReplaceAll("NORMWIDTH","");
708 }
709 else
710 opt.ReplaceAll("WIDTH","");
711 }
712
713 // if (opt.Contains("MULTIPROC")) {
714 // fitOption.ExecPolicy = ROOT::Fit::kMultiprocess;
715 // opt.ReplaceAll("MULTIPROC","");
716 // }
717
718 if (opt.Contains("SERIAL")) {
720 opt.ReplaceAll("SERIAL","");
721 }
722
723 if (opt.Contains("MULTITHREAD")) {
725 opt.ReplaceAll("MULTITHREAD","");
726 }
727
728 if (opt.Contains("I")) fitOption.Integral= 1; // integral of function in the bin (no sense for graph)
729 if (opt.Contains("W")) fitOption.W1 = 1; // all non-empty bins or points have weight =1 (for chi2 fit)
730 if (opt.Contains("WW")) fitOption.W1 = 2; //all bins have weight=1, even empty bins
731 if (opt.Contains("L")) fitOption.Like = 1;
732 if (opt.Contains("X")) fitOption.Chi2 = 1;
733 if (opt.Contains("P")) {
734 fitOption.PChi2 = 1;
735 if (fitOption.W1) { // option contains also w is a weighted Pearson chi2 fit
736 fitOption.PChi2 = 2;
737 fitOption.W1 = 0; // does not make sense to have errors=1 in Pearson chi2 fits
738 }
739 }
740
741 // specific likelihood fit options
742 if (fitOption.Like == 1) {
743 //if (opt.Contains("LL")) fitOption.Like = 2;
744 if (opt.Contains("W")){ fitOption.Like = 2; fitOption.W1=0;}// (weighted likelihood)
745 if (opt.Contains("MULTI")) {
746 if (fitOption.Like == 2) fitOption.Like = 6; // weighted multinomial
747 else fitOption.Like = 4; // multinomial likelihood fit instead of Poisson
748 opt.ReplaceAll("MULTI","");
749 }
750 // give precedence for likelihood options
751 if (fitOption.Chi2 || fitOption.PChi2 )
752 Warning("Fit","Cannot use P or X option in combination of L. Ignore the chi2 option and perform a likelihood fit");
753 }
754 }
755 // specific Graph options (need to be parsed before)
756 else if (type == EFitObjectType::kGraph) {
757 opt.ReplaceAll("ROB", "H");
758 opt.ReplaceAll("EX0", "T");
759
760 //for robust fitting, see if # of good points is defined
761 // decode parameters for robust fitting
762 Double_t h=0;
763 if (opt.Contains("H=0.")) {
764 int start = opt.Index("H=0.");
765 int numpos = start + strlen("H=0.");
766 int numlen = 0;
767 int len = opt.Length();
768 while( (numpos+numlen<len) && isdigit(opt[numpos+numlen]) ) numlen++;
769 TString num = opt(numpos,numlen);
770 opt.Remove(start+strlen("H"),strlen("=0.")+numlen);
771 h = atof(num.Data());
772 h*=TMath::Power(10, -numlen);
773 }
774
775 if (opt.Contains("H")) { fitOption.Robust = 1; fitOption.hRobust = h; }
776 if (opt.Contains("T")) fitOption.NoErrX = 1; // no error in X
777 if (opt.Contains("W")) fitOption.W1 = 1; // ignorer all point errors when fitting
778 }
779
780 if (opt.Contains("U")) {
781 // user option can work only when not running in multiple threads
783 fitOption.User = 1;
784 } else {
785 Warning("FitOptionsMake","Cannot use User (U) fit option when running in multi-thread mode. The option is ignored");
786 fitOption.User = 0;
787 }
788 }
789 if (opt.Contains("Q")) fitOption.Quiet = 1;
790 if (opt.Contains("V")) {fitOption.Verbose = 1; fitOption.Quiet = 0;}
791
792
793 if (opt.Contains("E")) fitOption.Errors = 1;
794 if (opt.Contains("R")) fitOption.Range = 1;
795 if (opt.Contains("G")) fitOption.Gradient= 1;
796 if (opt.Contains("M")) fitOption.More = 1;
797 if (opt.Contains("N")) fitOption.Nostore = 1;
798 if (opt.Contains("0")) fitOption.Nograph = 1;
799 if (opt.Contains("+")) fitOption.Plus = 1;
800 if (opt.Contains("B")) fitOption.Bound = 1;
801 if (opt.Contains("C")) fitOption.Nochisq = 1;
802 if (opt.Contains("F")) fitOption.Minuit = 1;
803 if (opt.Contains("S")) fitOption.StoreResult = 1;
804
805}
806
808 if (foption.Like) {
809 Info("CheckGraphFitOptions","L (Log Likelihood fit) is an invalid option when fitting a graph. It is ignored");
810 foption.Like = 0;
811 }
812 if (foption.Integral) {
813 Info("CheckGraphFitOptions","I (use function integral) is an invalid option when fitting a graph. It is ignored");
814 foption.Integral = 0;
815 }
816 return;
817}
818
819// implementation of unbin fit function (defined in HFitInterface)
820
822 // do unbin fit, ownership of fitdata is passed later to the TBackFitter class
823
824 // create a shared pointer to the fit data to managed it
825 std::shared_ptr<ROOT::Fit::UnBinData> fitdata(data);
826
827#ifdef DEBUG
828 printf("tree data size is %d \n",fitdata->Size());
829 for (unsigned int i = 0; i < fitdata->Size(); ++i) {
830 if (fitdata->NDim() == 1) printf(" x[%d] = %f \n", i,*(fitdata->Coords(i) ) );
831 }
832#endif
833 if (fitdata->Size() == 0 ) {
834 Warning("Fit","Fit data is empty ");
835 return -1;
836 }
837
838 // create an empty TFitResult
839 std::shared_ptr<TFitResult> tfr(new TFitResult() );
840 // create the fitter
841 std::shared_ptr<ROOT::Fit::Fitter> fitter(new ROOT::Fit::Fitter(tfr) );
842 ROOT::Fit::FitConfig & fitConfig = fitter->Config();
843
844 // dimension is given by data because TF1 pointer can have wrong one
845 unsigned int dim = fitdata->NDim();
846
847 // set the fit function
848 // if option grad is specified use gradient
849 // need to create a wrapper for an automatic normalized TF1 ???
850 if ( fitOption.Gradient ) {
851 assert ( (int) dim == fitfunc->GetNdim() );
852 fitter->SetFunction(ROOT::Math::WrappedMultiTF1(*fitfunc) );
853 }
854 else
855 fitter->SetFunction(static_cast<const ROOT::Math::IParamMultiFunction &>(ROOT::Math::WrappedMultiTF1(*fitfunc, dim) ) );
856
857 // parameter setting is done automaticaly in the Fitter class
858 // need only to set limits
859 int npar = fitfunc->GetNpar();
860 for (int i = 0; i < npar; ++i) {
861 ROOT::Fit::ParameterSettings & parSettings = fitConfig.ParSettings(i);
862 double plow,pup;
863 fitfunc->GetParLimits(i,plow,pup);
864 // this is a limitation of TF1 interface - cannot fix a parameter to zero value
865 if (plow*pup != 0 && plow >= pup) {
866 parSettings.Fix();
867 }
868 else if (plow < pup ) {
869 if (!TMath::Finite(pup) && TMath::Finite(plow) )
870 parSettings.SetLowerLimit(plow);
871 else if (!TMath::Finite(plow) && TMath::Finite(pup) )
872 parSettings.SetUpperLimit(pup);
873 else
874 parSettings.SetLimits(plow,pup);
875 }
876
877 // set the parameter step size (by default are set to 0.3 of value)
878 // if function provides meaningful error values
879 double err = fitfunc->GetParError(i);
880 if ( err > 0)
881 parSettings.SetStepSize(err);
882 else if (plow < pup && TMath::Finite(plow) && TMath::Finite(pup) ) { // in case of limits improve step sizes
883 double step = 0.1 * (pup - plow);
884 // check if value is not too close to limit otherwise trim value
885 if ( parSettings.Value() < pup && pup - parSettings.Value() < 2 * step )
886 step = (pup - parSettings.Value() ) / 2;
887 else if ( parSettings.Value() > plow && parSettings.Value() - plow < 2 * step )
888 step = (parSettings.Value() - plow ) / 2;
889
890 parSettings.SetStepSize(step);
891 }
892
893 }
894
895 fitConfig.SetMinimizerOptions(minOption);
896
897 if (fitOption.Verbose) fitConfig.MinimizerOptions().SetPrintLevel(3);
898 if (fitOption.Quiet) fitConfig.MinimizerOptions().SetPrintLevel(0);
899
900 // more
901 if (fitOption.More) fitConfig.SetMinimizer("Minuit","MigradImproved");
902
903 // chech if Minos or more options
904 if (fitOption.Errors) {
905 // run Hesse and Minos
906 fitConfig.SetParabErrors(true);
907 fitConfig.SetMinosErrors(true);
908 }
909 // use weight correction
910 if ( (fitOption.Like & 2) == 2)
911 fitConfig.SetWeightCorrection(true);
912
913 bool extended = (fitOption.Like & 1) == 1;
914
915 bool fitok = false;
916 fitok = fitter->LikelihoodFit(fitdata, extended, fitOption.ExecPolicy);
917 if ( !fitok && !fitOption.Quiet )
918 Warning("UnBinFit","Abnormal termination of minimization.");
919
920 const ROOT::Fit::FitResult & fitResult = fitter->Result();
921 // one could set directly the fit result in TF1
922 int iret = fitResult.Status();
923 if (!fitResult.IsEmpty() ) {
924 // set in fitfunc the result of the fit
925 fitfunc->SetNDF(fitResult.Ndf() );
926 fitfunc->SetNumberFitPoints(fitdata->Size() );
927
928 assert( (Int_t)fitResult.Parameters().size() >= fitfunc->GetNpar() );
929 fitfunc->SetParameters( const_cast<double*>(&(fitResult.Parameters().front())));
930 if ( int( fitResult.Errors().size()) >= fitfunc->GetNpar() )
931 fitfunc->SetParErrors( &(fitResult.Errors().front()) );
932
933 }
934
935 // store fitting result in the backward compatible TVirtualFitter object
936 // lock in case running in a multi-thread enabled mode
937 {
940 TBackCompFitter * bcfitter = new TBackCompFitter(fitter, fitdata);
941 // cannot use anymore now fitdata (given away ownership)
942 fitdata = nullptr;
943 bcfitter->SetFitOption(fitOption);
944 bcfitter->SetUserFunc(fitfunc);
945 // delete previous fitter and replace with the new one
946 if (lastFitter) delete lastFitter;
947 TVirtualFitter::SetFitter( bcfitter );
948
949 // print results
950 if (fitOption.Verbose) fitResult.PrintCovMatrix(std::cout);
951 else if (!fitOption.Quiet) fitResult.Print(std::cout);
952 }
953
954
955 if (fitOption.StoreResult)
956 {
957 TString name = "TFitResult-";
958 name = name + "UnBinData-" + fitfunc->GetName();
959 TString title = "TFitResult-";
960 title += name;
961 tfr->SetName(name);
962 tfr->SetTitle(title);
963 return TFitResultPtr(tfr);
964 }
965 else
966 return TFitResultPtr(iret);
967}
968
969
970// implementations of ROOT::Fit::FitObject functions (defined in HFitInterface) in terms of the template HFit::Fit
971
973moption, const char *goption, ROOT::Fit::DataRange & range) {
974 // check fit options
975 // check if have weights in case of weighted likelihood
976 if ( ((foption.Like & 2) == 2) && h1->GetSumw2N() == 0) {
977 Warning("HFit::FitObject","A weighted likelihood fit is requested but histogram is not weighted - do a standard Likelihood fit");
978 foption.Like = 1;
979 }
980 // histogram fitting
981 return HFit::Fit(h1,f1,foption,moption,goption,range);
982}
983
984TFitResultPtr ROOT::Fit::FitObject(TGraph * 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 // TGraph fitting
988 return HFit::Fit(gr,f1,foption,moption,goption,range);
989}
990
991TFitResultPtr ROOT::Fit::FitObject(TMultiGraph * 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 // TMultiGraph fitting
995 return HFit::Fit(gr,f1,foption,moption,goption,range);
996}
997
998TFitResultPtr ROOT::Fit::FitObject(TGraph2D * gr, TF1 *f1 , Foption_t & foption , const ROOT::Math::MinimizerOptions & moption, const char *goption, ROOT::Fit::DataRange & range) {
999 // exclude options not valid for graphs
1001 // TGraph2D fitting
1002 return HFit::Fit(gr,f1,foption,moption,goption,range);
1003}
1004
1005TFitResultPtr ROOT::Fit::FitObject(THnBase * s1, TF1 *f1 , Foption_t & foption , const ROOT::Math::MinimizerOptions & moption, const char *goption, ROOT::Fit::DataRange & range) {
1006 // sparse histogram fitting
1007 return HFit::Fit(s1,f1,foption,moption,goption,range);
1008}
1009
1010
1011
1012// Int_t TGraph2D::DoFit(TF2 *f2 ,Option_t *option ,Option_t *goption) {
1013// // internal graph2D fitting methods
1014// Foption_t fitOption;
1015// ROOT::Fit::FitOptionsMake(option,fitOption);
1016
1017// // create range and minimizer options with default values
1018// ROOT::Fit::DataRange range(2);
1019// ROOT::Math::MinimizerOptions minOption;
1020// return ROOT::Fit::FitObject(this, f2 , fitOption , minOption, goption, range);
1021// }
1022
1023
1024// function to compute the simple chi2 for graphs and histograms
1025
1026
1027double ROOT::Fit::Chisquare(const TH1 & h1, TF1 & f1, bool useRange, ROOT::Fit::EChisquareType type) {
1028 return HFit::ComputeChi2(h1,f1,useRange, type);
1029}
1030
1031double ROOT::Fit::Chisquare(const TGraph & g, TF1 & f1, bool useRange) {
1033}
1034
1035template<class FitObject>
1036double HFit::ComputeChi2(const FitObject & obj, TF1 & f1, bool useRange, ROOT::Fit::EChisquareType type ) {
1037
1038 // implement using the fitting classes
1040 opt.fUseEmpty = (type != ROOT::Fit::EChisquareType::kNeyman); // use empty bin when not using Neyman chisquare (observed error)
1042 opt.fErrors1 = (type == ROOT::Fit::EChisquareType::kPearson); // not using observed errors in Pearson chi2
1043
1045 // get range of function
1046 if (useRange) HFit::GetFunctionRange(f1,range);
1047 // fill the data set
1048 ROOT::Fit::BinData data(opt,range);
1049 ROOT::Fit::FillData(data, &obj, &f1);
1050 if (data.Size() == 0 ) {
1051 Warning("Chisquare","data set is empty - return -1");
1052 return -1;
1053 }
1056 // use the poisson log-lokelihood (Baker-Cousins chi2)
1058 return 2.* nll( f1.GetParameters() ) ;
1059 }
1060 ROOT::Fit::Chi2Function chi2(data, wf1);
1061 return chi2(f1.GetParameters() );
1062
1063}
#define g(i)
Definition RSha256.hxx:105
#define s1(x)
Definition RSha256.hxx:91
#define h(i)
Definition RSha256.hxx:106
constexpr Bool_t kFALSE
Definition RtypesCore.h:94
constexpr Bool_t kTRUE
Definition RtypesCore.h:93
const Bool_t kIterBackward
Definition TCollection.h:43
#define R__ASSERT(e)
Checks condition e and reports a fatal error if it's false.
Definition TError.h:125
void Info(const char *location, const char *msgfmt,...)
Use this function for informational messages.
Definition TError.cxx:218
void Error(const char *location, const char *msgfmt,...)
Use this function in case an error occurred.
Definition TError.cxx:185
void Warning(const char *location, const char *msgfmt,...)
Use this function in warning situations.
Definition TError.cxx:229
Option_t Option_t option
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void data
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 UChar_t len
Option_t Option_t TPoint TPoint const char x2
Option_t Option_t TPoint TPoint const char x1
Option_t Option_t TPoint TPoint const char y2
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
Option_t Option_t TPoint TPoint const char y1
char name[80]
Definition TGX11.cxx:110
float xmin
float ymin
float xmax
float ymax
@ kCanDelete
Definition TObject.h:367
R__EXTERN TVirtualMutex * gROOTMutex
Definition TROOT.h:63
R__EXTERN TVirtualMutex * gGlobalMutex
#define R__LOCKGUARD(mutex)
#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 binned 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 SetMinimizer(const char *type, const char *algo=nullptr)
set minimizer type
Definition FitConfig.h:179
void SetMinosErrors(bool on=true)
set Minos errors computation to be performed after fitting
Definition FitConfig.h:229
void SetNormErrors(bool on=true)
set the option to normalize the error on the result according to chi2/ndf
Definition FitConfig.h:223
bool NormalizeErrors() const
flag to check if resulting errors are be normalized according to chi2/ndf
Definition FitConfig.h:202
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 matrix computation
Definition FitConfig.h:232
void SetParabErrors(bool on=true)
set parabolic errors
Definition FitConfig.h:226
const std::string & MinimizerType() const
return type of minimizer package
Definition FitConfig.h:187
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
class containing 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:108
const std::vector< double > & Errors() const
parameter errors (return st::vector)
Definition FitResult.h:162
const std::vector< double > & Parameters() const
parameter values (return std::vector)
Definition FitResult.h:167
unsigned int Ndf() const
Number of degree of freedom.
Definition FitResult.h:156
double Chi2() const
Return the Chi2 value after fitting In case of unbinned fits (or not defined one, see the documentati...
Definition FitResult.h:153
void Print(std::ostream &os, bool covmat=false) const
print the result and optionally covariance matrix and correlations
void PrintCovMatrix(std::ostream &os) const
print error matrix and correlations
int Status() const
minimizer status code
Definition FitResult.h:128
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
return parameter value
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 un-binned data sets (just x coordinates values) of any dimensions.
Definition UnBinData.h:46
IParamFunction interface (abstract class) describing multi-dimensional parametric functions It is a d...
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:31
Double_t GetXmax() const
Definition TAxis.h:140
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:139
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.
void SetFCN(void(*fcn)(Int_t &, Double_t *, Double_t &f, Double_t *, Int_t)) override
Override setFCN to use the Adapter to Minuit2 FCN interface To set the address of the minimization fu...
void * New(ENewType defConstructor=kClassNew, Bool_t quiet=kFALSE) const
Return a pointer to a newly allocated object of this class.
Definition TClass.cxx:5045
1-Dim function class
Definition TF1.h:233
virtual Int_t GetNumber() const
Definition TF1.h:526
virtual void GetParLimits(Int_t ipar, Double_t &parmin, Double_t &parmax) const
Return limits for parameter ipar.
Definition TF1.cxx:1940
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:3419
virtual Double_t GetParError(Int_t ipar) const
Return value of parameter number ipar.
Definition TF1.cxx:1930
static TClass * Class()
virtual void SetChisquare(Double_t chi2)
Definition TF1.h:640
virtual void SetRange(Double_t xmin, Double_t xmax)
Initialize the upper and lower bounds to draw the function.
Definition TF1.cxx:3528
virtual Int_t GetNpar() const
Definition TF1.h:509
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:3490
virtual void SetParent(TObject *p=nullptr)
Definition TF1.h:706
virtual Double_t * GetParameters() const
Definition TF1.h:548
void Copy(TObject &f1) const override
Copy this F1 to a new F1.
Definition TF1.cxx:1005
virtual void SetNumberFitPoints(Int_t npfits)
Definition TF1.h:652
@ kNotDraw
Definition TF1.h:346
virtual void GetRange(Double_t *xmin, Double_t *xmax) const
Return range of a generic N-D function.
Definition TF1.cxx:2281
virtual Bool_t IsLinear() const
Definition TF1.h:630
bool IsVectorized()
Definition TF1.h:467
virtual void SetParameters(const Double_t *params)
Definition TF1.h:677
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:3161
TClass * IsA() const override
Definition TF1.h:755
virtual Int_t GetNdim() const
Definition TF1.h:513
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:847
A 2-Dim function with parameters.
Definition TF2.h:29
void Save(Double_t xmin, Double_t xmax, Double_t ymin, Double_t ymax, Double_t zmin, Double_t zmax) override
Save values of function in array fSave.
Definition TF2.cxx:796
void SetRange(Double_t xmin, Double_t xmax) override
Initialize the upper and lower bounds to draw the function.
Definition TF2.h:146
A 3-Dim function with parameters.
Definition TF3.h:28
void Save(Double_t xmin, Double_t xmax, Double_t ymin, Double_t ymax, Double_t zmin, Double_t zmax) override
Save values of function in array fSave.
Definition TF3.cxx:551
void SetRange(Double_t xmin, Double_t xmax) override
Initialize the upper and lower bounds to draw the function.
Definition TF3.h:143
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 TH1F * GetHistogram() const
Returns a pointer to the histogram used to draw the axis Takes into account the two following cases.
Definition TGraph.cxx:1428
TH1 is the base class of all histogram classes in ROOT.
Definition TH1.h:59
TAxis * GetZaxis()
Definition TH1.h:326
static TClass * Class()
virtual Int_t GetDimension() const
Definition TH1.h:283
TAxis * GetXaxis()
Definition TH1.h:324
TAxis * GetYaxis()
Definition TH1.h:325
void Draw(Option_t *option="") override
Draw this histogram with options.
Definition TH1.cxx:3066
TList * GetListOfFunctions() const
Definition TH1.h:244
virtual Int_t GetSumw2N() const
Definition TH1.h:315
Multidimensional histogram base.
Definition THnBase.h:43
A doubly linked list.
Definition TList.h:38
void Add(TObject *obj) override
Definition TList.h:83
TObject * Remove(TObject *obj) override
Remove object from the list.
Definition TList.cxx:820
A TMultiGraph is a collection of TGraph (or derived) objects.
Definition TMultiGraph.h:34
TList * GetListOfGraphs() const
Definition TMultiGraph.h:68
TH1F * GetHistogram()
Returns a pointer to the histogram used to draw the axis.
const char * GetName() const override
Returns name of object.
Definition TNamed.h:47
const char * GetTitle() const override
Returns title of object.
Definition TNamed.h:48
Mother of all ROOT objects.
Definition TObject.h:41
R__ALWAYS_INLINE Bool_t TestBit(UInt_t f) const
Definition TObject.h:199
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:786
virtual Bool_t InheritsFrom(const char *classname) const
Returns kTRUE if object inherits from class "classname".
Definition TObject.cxx:530
Basic string class.
Definition TString.h:139
Ssiz_t Length() const
Definition TString.h:417
const char * Data() const
Definition TString.h:376
TString & ReplaceAll(const TString &s1, const TString &s2)
Definition TString.h:704
void ToUpper()
Change string to upper case.
Definition TString.cxx:1195
TString & Remove(Ssiz_t pos)
Definition TString.h:685
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition TString.h:632
Ssiz_t Index(const char *pat, Ssiz_t i=0, ECaseCompare cmp=kExact) const
Definition TString.h:651
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:462
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
double ComputeChi2(const FitObject &h1, TF1 &f1, bool useRange, ROOT::Fit::EChisquareType type)
void FitOptionsMake(const char *option, Foption_t &fitOption)
void CheckGraphFitOptions(Foption_t &fitOption)
Definition HFitImpl.cxx:807
void StoreAndDrawFitFunction(FitObject *h1, TF1 *f1, const ROOT::Fit::DataRange &range, bool, bool, const char *goption)
Definition HFitImpl.cxx:570
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:972
double Chisquare(const TH1 &h1, TF1 &f1, bool useRange, EChisquareType type)
compute the chi2 value for an histogram given a function (see TH1::Chisquare for the documentation)
void FitOptionsMake(EFitObjectType type, const char *option, Foption_t &fitOption)
Decode list of options into fitOption.
Definition HFitImpl.cxx:685
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:821
void FillData(BinData &dv, const TH1 *hist, TF1 *func=nullptr)
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 ...
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:570
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:770
LongDouble_t Power(LongDouble_t x, LongDouble_t y)
Returns x raised to the power y.
Definition TMath.h:721
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 all errors equal to 1, i.e. fit without errors (default is false)
Definition DataOptions.h:52
bool fNormBinVolume
normalize data by a normalized the bin volume (bin volume divided by a reference value)
Definition DataOptions.h:49
bool fUseRange
use the function range when creating the fit data (default is false)
Definition DataOptions.h:51
bool fUseEmpty
use empty bins (default is false) with a fixed error of 1
Definition DataOptions.h:50
bool fIntegral
use integral of bin content instead of bin center (default is false)
Definition DataOptions.h:47
bool fExpErrors
use expected errors from the function and not from the data
Definition DataOptions.h:53
bool fBinVolume
normalize data by the bin volume (it is used in the Poisson likelihood fits)
Definition DataOptions.h:48
bool fCoordErrors
use errors on the x coordinates when available (default is true)
Definition DataOptions.h:54