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