ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
TFitter.cxx
Go to the documentation of this file.
1 // @(#)root/minuit:$Id$
2 // Author: Rene Brun 31/08/99
3 /*************************************************************************
4  * Copyright (C) 1995-2000, Rene Brun and Fons Rademakers. *
5  * All rights reserved. *
6  * *
7  * For the licensing terms see $ROOTSYS/LICENSE. *
8  * For the list of contributors see $ROOTSYS/README/CREDITS. *
9  *************************************************************************/
10 
11 #include "TMinuit.h"
12 #include "TFitter.h"
13 #include "TH1.h"
14 #include "TF1.h"
15 #include "TF2.h"
16 #include "TF3.h"
17 #include "TList.h"
18 #include "TGraph.h"
19 #include "TGraph2D.h"
20 #include "TMultiGraph.h"
21 #include "TMath.h"
22 
23 extern void H1FitChisquare(Int_t &npar, Double_t *gin, Double_t &f, Double_t *u, Int_t flag);
24 extern void H1FitLikelihood(Int_t &npar, Double_t *gin, Double_t &f, Double_t *u, Int_t flag);
25 extern void GraphFitChisquare(Int_t &npar, Double_t *gin, Double_t &f, Double_t *u, Int_t flag);
26 extern void Graph2DFitChisquare(Int_t &npar, Double_t *gin, Double_t &f, Double_t *u, Int_t flag);
27 extern void MultiGraphFitChisquare(Int_t &npar, Double_t *gin, Double_t &f, Double_t *u, Int_t flag);
28 extern void F2Fit(Int_t &npar, Double_t *gin, Double_t &f, Double_t *u, Int_t flag);
29 extern void F3Fit(Int_t &npar, Double_t *gin, Double_t &f, Double_t *u, Int_t flag);
30 
32 
33 ////////////////////////////////////////////////////////////////////////////////
34 ///*-*-*-*-*-*-*-*-*-*-*default constructor*-*-*-*-*-*-*-*-*-*-*-*-*
35 ///*-* ===================
36 
37 TFitter::TFitter(Int_t maxpar)
38 {
39  fMinuit = new TMinuit(maxpar);
40  fNlog = 0;
41  fSumLog = 0;
42  fCovar = 0;
43  SetName("MinuitFitter");
44 }
45 
46 ////////////////////////////////////////////////////////////////////////////////
47 ///*-*-*-*-*-*-*-*-*-*-*default destructor*-*-*-*-*-*-*-*-*-*-*-*-*-*
48 ///*-* ==================
49 
51 {
52  if (fCovar) delete [] fCovar;
53  if (fSumLog) delete [] fSumLog;
54  delete fMinuit;
55 }
56 
57 ////////////////////////////////////////////////////////////////////////////////
58 /// return a chisquare equivalent
59 
61 {
62  Double_t amin = 0;
63  H1FitChisquare(npar,params,amin,params,1);
64  return amin;
65 }
66 
67 ////////////////////////////////////////////////////////////////////////////////
68 /// reset the fitter environment
69 
71 {
72  if (fCovar) {delete [] fCovar; fCovar = 0;}
73  fMinuit->mncler();
74 
75  //reset the internal Minuit random generator to its initial state
76  Double_t val = 3;
77  Int_t inseed = 12345;
78  fMinuit->mnrn15(val,inseed);
79 }
80 
81 ////////////////////////////////////////////////////////////////////////////////
82 /// Execute a fitter command;
83 /// command : command string
84 /// args : list of nargs command arguments
85 
86 Int_t TFitter::ExecuteCommand(const char *command, Double_t *args, Int_t nargs)
87 {
88  if (fCovar) {delete [] fCovar; fCovar = 0;}
89  Int_t ierr = 0;
90  fMinuit->mnexcm(command,args,nargs,ierr);
91  return ierr;
92 }
93 
94 ////////////////////////////////////////////////////////////////////////////////
95 /// Fix parameter ipar.
96 
98 {
99  if (fCovar) {delete [] fCovar; fCovar = 0;}
100  fMinuit->FixParameter(ipar);
101 }
102 
103 ////////////////////////////////////////////////////////////////////////////////
104 ///Computes point-by-point confidence intervals for the fitted function
105 ///Parameters:
106 ///n - number of points
107 ///ndim - dimensions of points
108 ///x - points, at which to compute the intervals, for ndim > 1
109 /// should be in order: (x0,y0, x1, y1, ... xn, yn)
110 ///ci - computed intervals are returned in this array
111 ///cl - confidence level, default=0.95
112 ///NOTE, that the intervals are approximate for nonlinear(in parameters) models
113 
115 {
116  TF1 *f = (TF1*)fUserFunc;
117  Int_t npar = f->GetNumberFreeParameters();
118  Int_t npar_real = f->GetNpar();
119  Double_t *grad = new Double_t[npar_real];
120  Double_t *sum_vector = new Double_t[npar];
121  Bool_t *fixed=0;
122  Double_t al, bl;
123  if (npar_real != npar){
124  fixed = new Bool_t[npar_real];
125  memset(fixed,0,npar_real*sizeof(Bool_t));
126 
127  for (Int_t ipar=0; ipar<npar_real; ipar++){
128  fixed[ipar]=0;
129  f->GetParLimits(ipar,al,bl);
130  if (al*bl != 0 && al >= bl) {
131  //this parameter is fixed
132  fixed[ipar]=1;
133  }
134  }
135  }
136  Double_t c=0;
137 
138  Double_t *matr = GetCovarianceMatrix();
139  if (!matr)
140  return;
141  Double_t t = TMath::StudentQuantile(0.5 + cl/2, f->GetNDF());
142  Double_t chidf = TMath::Sqrt(f->GetChisquare()/f->GetNDF());
143  Int_t igrad, ifree=0;
144  for (Int_t ipoint=0; ipoint<n; ipoint++){
145  c=0;
146  f->GradientPar(x+ndim*ipoint, grad);
147  //multiply the covariance matrix by gradient
148  for (Int_t irow=0; irow<npar; irow++){
149  sum_vector[irow]=0;
150  igrad = 0;
151  for (Int_t icol=0; icol<npar; icol++){
152  igrad = 0;
153  ifree=0;
154  if (fixed) {
155  //find the free parameter #icol
156  while (ifree<icol+1){
157  if (fixed[igrad]==0) ifree++;
158  igrad++;
159  }
160  igrad--;
161  //now the [igrad] element of gradient corresponds to [icol] element of cov.matrix
162  } else {
163  igrad = icol;
164  }
165  sum_vector[irow]+=matr[irow*npar_real+icol]*grad[igrad];
166  }
167  }
168  igrad = 0;
169  for (Int_t i=0; i<npar; i++){
170  igrad = 0; ifree=0;
171  if (fixed) {
172  //find the free parameter #icol
173  while (ifree<i+1){
174  if (fixed[igrad]==0) ifree++;
175  igrad++;
176  }
177  igrad--;
178  } else {
179  igrad = i;
180  }
181  c+=grad[igrad]*sum_vector[i];
182  }
183 
184  c=TMath::Sqrt(c);
185  ci[ipoint]=c*t*chidf;
186  }
187 
188  delete [] grad;
189  delete [] sum_vector;
190  if (fixed)
191  delete [] fixed;
192 }
193 
194 ////////////////////////////////////////////////////////////////////////////////
195 ///Computes confidence intervals at level cl. Default is 0.95
196 ///The TObject parameter can be a TGraphErrors, a TGraph2DErrors or a TH1,2,3.
197 ///For Graphs, confidence intervals are computed for each point,
198 ///the value of the graph at that point is set to the function value at that
199 ///point, and the graph y-errors (or z-errors) are set to the value of
200 ///the confidence interval at that point.
201 ///For Histograms, confidence intervals are computed for each bin center
202 ///The bin content of this bin is then set to the function value at the bin
203 ///center, and the bin error is set to the confidence interval value.
204 ///NOTE: confidence intervals are approximate for nonlinear models!
205 ///
206 ///Allowed combinations:
207 ///Fitted object Passed object
208 ///TGraph TGraphErrors, TH1
209 ///TGraphErrors, AsymmErrors TGraphErrors, TH1
210 ///TH1 TGraphErrors, TH1
211 ///TGraph2D TGraph2DErrors, TH2
212 ///TGraph2DErrors TGraph2DErrors, TH2
213 ///TH2 TGraph2DErrors, TH2
214 ///TH3 TH3
215 
217 {
218  if (obj->InheritsFrom(TGraph::Class())) {
219  TGraph *gr = (TGraph*)obj;
220  if (!gr->GetEY()){
221  Error("GetConfidenceIntervals", "A TGraphErrors should be passed instead of a graph");
222  return;
223  }
225  Error("GetConfidenceIntervals", "A TGraph2DErrors should be passed instead of a graph");
226  return;
227  }
229  if (((TH1*)(fObjectFit))->GetDimension()>1){
230  Error("GetConfidenceIntervals", "A TGraph2DErrors or a TH23 should be passed instead of a graph");
231  return;
232  }
233  }
234  GetConfidenceIntervals(gr->GetN(),1,gr->GetX(), gr->GetEY(), cl);
235  for (Int_t i=0; i<gr->GetN(); i++)
236  gr->SetPoint(i, gr->GetX()[i], ((TF1*)(fUserFunc))->Eval(gr->GetX()[i]));
237  }
238 
239  //TGraph2D/////////////////
240  else if (obj->InheritsFrom(TGraph2D::Class())) {
241  TGraph2D *gr2 = (TGraph2D*)obj;
242  if (!gr2->GetEZ()){
243  Error("GetConfidenceIntervals", "A TGraph2DErrors should be passed instead of a TGraph2D");
244  return;
245  }
247  Error("GetConfidenceIntervals", "A TGraphErrors should be passed instead of a TGraph2D");
248  return;
249  }
251  if (((TH1*)(fObjectFit))->GetDimension()==1){
252  Error("GetConfidenceIntervals", "A TGraphErrors or a TH1 should be passed instead of a graph");
253  return;
254  }
255  }
256  TF2 *f=(TF2*)fUserFunc;
257  Double_t xy[2];
258  Int_t np = gr2->GetN();
259  Int_t npar = f->GetNpar();
260  Double_t *grad = new Double_t[npar];
261  Double_t *sum_vector = new Double_t[npar];
262  Double_t *x = gr2->GetX();
263  Double_t *y = gr2->GetY();
264  Double_t t = TMath::StudentQuantile(0.5 + cl/2, f->GetNDF());
265  Double_t chidf = TMath::Sqrt(f->GetChisquare()/f->GetNDF());
267  Double_t c = 0;
268  for (Int_t ipoint=0; ipoint<np; ipoint++){
269  xy[0]=x[ipoint];
270  xy[1]=y[ipoint];
271  f->GradientPar(xy, grad);
272  for (Int_t irow=0; irow<f->GetNpar(); irow++){
273  sum_vector[irow]=0;
274  for (Int_t icol=0; icol<npar; icol++)
275  sum_vector[irow]+=matr[irow*npar+icol]*grad[icol];
276  }
277  c = 0;
278  for (Int_t i=0; i<npar; i++)
279  c+=grad[i]*sum_vector[i];
280  c=TMath::Sqrt(c);
281  gr2->SetPoint(ipoint, xy[0], xy[1], f->EvalPar(xy));
282  gr2->GetEZ()[ipoint]=c*t*chidf;
283 
284  }
285  delete [] grad;
286  delete [] sum_vector;
287  }
288 
289  //TH1/////////////////////////
290  else if (obj->InheritsFrom(TH1::Class())) {
292  if (((TH1*)obj)->GetDimension()>1){
293  Error("GetConfidenceIntervals", "Fitted graph and passed histogram have different number of dimensions");
294  return;
295  }
296  }
298  if (((TH1*)obj)->GetDimension()!=2){
299  Error("GetConfidenceIntervals", "Fitted graph and passed histogram have different number of dimensions");
300  return;
301  }
302  }
304  if (((TH1*)(fObjectFit))->GetDimension()!=((TH1*)(obj))->GetDimension()){
305  Error("GetConfidenceIntervals", "Fitted and passed histograms have different number of dimensions");
306  return;
307  }
308  }
309 
310 
311  TH1 *hfit = (TH1*)obj;
312  TF1 *f = (TF1*)GetUserFunc();
313  Int_t npar = f->GetNpar();
314  Double_t *grad = new Double_t[npar];
315  Double_t *sum_vector = new Double_t[npar];
316  Double_t x[3];
317 
318  Int_t hxfirst = hfit->GetXaxis()->GetFirst();
319  Int_t hxlast = hfit->GetXaxis()->GetLast();
320  Int_t hyfirst = hfit->GetYaxis()->GetFirst();
321  Int_t hylast = hfit->GetYaxis()->GetLast();
322  Int_t hzfirst = hfit->GetZaxis()->GetFirst();
323  Int_t hzlast = hfit->GetZaxis()->GetLast();
324 
325  TAxis *xaxis = hfit->GetXaxis();
326  TAxis *yaxis = hfit->GetYaxis();
327  TAxis *zaxis = hfit->GetZaxis();
328  Double_t t = TMath::StudentQuantile(0.5 + cl/2, f->GetNDF());
329  Double_t chidf = TMath::Sqrt(f->GetChisquare()/f->GetNDF());
331  Double_t c=0;
332  for (Int_t binz=hzfirst; binz<=hzlast; binz++){
333  x[2]=zaxis->GetBinCenter(binz);
334  for (Int_t biny=hyfirst; biny<=hylast; biny++) {
335  x[1]=yaxis->GetBinCenter(biny);
336  for (Int_t binx=hxfirst; binx<=hxlast; binx++) {
337  x[0]=xaxis->GetBinCenter(binx);
338  f->GradientPar(x, grad);
339  for (Int_t irow=0; irow<npar; irow++){
340  sum_vector[irow]=0;
341  for (Int_t icol=0; icol<npar; icol++)
342  sum_vector[irow]+=matr[irow*npar+icol]*grad[icol];
343  }
344  c = 0;
345  for (Int_t i=0; i<npar; i++)
346  c+=grad[i]*sum_vector[i];
347  c=TMath::Sqrt(c);
348  hfit->SetBinContent(binx, biny, binz, f->EvalPar(x));
349  hfit->SetBinError(binx, biny, binz, c*t*chidf);
350  }
351  }
352  }
353  delete [] grad;
354  delete [] sum_vector;
355  }
356  else {
357  Error("GetConfidenceIntervals", "This object type is not supported");
358  return;
359  }
360 
361 }
362 
363 ////////////////////////////////////////////////////////////////////////////////
364 /// return a pointer to the covariance matrix
365 
367 {
368  if (fCovar) return fCovar;
369  Int_t npars = fMinuit->GetNumPars();
370  ((TFitter*)this)->fCovar = new Double_t[npars*npars];
371  fMinuit->mnemat(fCovar,npars);
372  return fCovar;
373 }
374 
375 ////////////////////////////////////////////////////////////////////////////////
376 /// return element i,j from the covariance matrix
377 
379 {
381  Int_t npars = fMinuit->GetNumPars();
382  if (i < 0 || i >= npars || j < 0 || j >= npars) {
383  Error("GetCovarianceMatrixElement","Illegal arguments i=%d, j=%d",i,j);
384  return 0;
385  }
386  return fCovar[j+npars*i];
387 }
388 
389 ////////////////////////////////////////////////////////////////////////////////
390 /// return current errors for a parameter
391 /// ipar : parameter number
392 /// eplus : upper error
393 /// eminus : lower error
394 /// eparab : parabolic error
395 /// globcc : global correlation coefficient
396 
397 Int_t TFitter::GetErrors(Int_t ipar,Double_t &eplus, Double_t &eminus, Double_t &eparab, Double_t &globcc) const
398 {
399 
400  Int_t ierr = 0;
401  fMinuit->mnerrs(ipar, eplus,eminus,eparab,globcc);
402  return ierr;
403 }
404 
405 
406 ////////////////////////////////////////////////////////////////////////////////
407 /// return the total number of parameters (free + fixed)
408 
410 {
411  return fMinuit->fNpar + fMinuit->fNpfix;
412 }
413 
414 ////////////////////////////////////////////////////////////////////////////////
415 /// return the number of free parameters
416 
418 {
419  return fMinuit->fNpar;
420 }
421 
422 
423 ////////////////////////////////////////////////////////////////////////////////
424 /// return error of parameter ipar
425 
427 {
428  Int_t ierr = 0;
429  TString pname;
430  Double_t value,verr,vlow,vhigh;
431 
432  fMinuit->mnpout(ipar, pname,value,verr,vlow,vhigh,ierr);
433  return verr;
434 }
435 
436 
437 ////////////////////////////////////////////////////////////////////////////////
438 /// return current value of parameter ipar
439 
441 {
442  Int_t ierr = 0;
443  TString pname;
444  Double_t value,verr,vlow,vhigh;
445 
446  fMinuit->mnpout(ipar, pname,value,verr,vlow,vhigh,ierr);
447  return value;
448 }
449 
450 ////////////////////////////////////////////////////////////////////////////////
451 /// return current values for a parameter
452 /// ipar : parameter number
453 /// parname : parameter name
454 /// value : initial parameter value
455 /// verr : initial error for this parameter
456 /// vlow : lower value for the parameter
457 /// vhigh : upper value for the parameter
458 /// WARNING! parname must be suitably dimensionned in the calling function.
459 
460 Int_t TFitter::GetParameter(Int_t ipar, char *parname,Double_t &value,Double_t &verr,Double_t &vlow, Double_t &vhigh) const
461 {
462  Int_t ierr = 0;
463  TString pname;
464  fMinuit->mnpout(ipar, pname,value,verr,vlow,vhigh,ierr);
465  strcpy(parname,pname.Data());
466  return ierr;
467 }
468 
469 ////////////////////////////////////////////////////////////////////////////////
470 /// return name of parameter ipar
471 
472 const char *TFitter::GetParName(Int_t ipar) const
473 {
474  if (!fMinuit || ipar < 0 || ipar > fMinuit->fNu) return "";
475  return fMinuit->fCpnam[ipar];
476 }
477 
478 ////////////////////////////////////////////////////////////////////////////////
479 /// return global fit parameters
480 /// amin : chisquare
481 /// edm : estimated distance to minimum
482 /// errdef
483 /// nvpar : number of variable parameters
484 /// nparx : total number of parameters
485 
486 Int_t TFitter::GetStats(Double_t &amin, Double_t &edm, Double_t &errdef, Int_t &nvpar, Int_t &nparx) const
487 {
488  Int_t ierr = 0;
489  fMinuit->mnstat(amin,edm,errdef,nvpar,nparx,ierr);
490  return ierr;
491 }
492 
493 ////////////////////////////////////////////////////////////////////////////////
494 /// return Sum(log(i) i=0,n
495 /// used by log likelihood fits
496 
498 {
499  if (n < 0) return 0;
500  if (n > fNlog) {
501  if (fSumLog) delete [] fSumLog;
502  fNlog = 2*n+1000;
503  fSumLog = new Double_t[fNlog+1];
504  Double_t fobs = 0;
505  for (Int_t j=0;j<=fNlog;j++) {
506  if (j > 1) fobs += TMath::Log(j);
507  fSumLog[j] = fobs;
508  }
509  }
510  if (fSumLog) return fSumLog[n];
511  return 0;
512 }
513 
514 
515 ////////////////////////////////////////////////////////////////////////////////
516 ///return kTRUE if parameter ipar is fixed, kFALSE othersise)
517 
519 {
520  if (fMinuit->fNiofex[ipar] == 0 ) return kTRUE;
521  return kFALSE;
522 }
523 
524 
525 ////////////////////////////////////////////////////////////////////////////////
526 /// Print fit results
527 
528 void TFitter::PrintResults(Int_t level, Double_t amin) const
529 {
530  fMinuit->mnprin(level,amin);
531 }
532 
533 ////////////////////////////////////////////////////////////////////////////////
534 /// Release parameter ipar.
535 
537 {
538  if (fCovar) {delete [] fCovar; fCovar = 0;}
539  fMinuit->Release(ipar);
540 }
541 
542 ////////////////////////////////////////////////////////////////////////////////
543 /// Specify the address of the fitting algorithm (from the interpreter)
544 
545 void TFitter::SetFCN(void *fcn)
546 {
547  if (fCovar) {delete [] fCovar; fCovar = 0;}
549  fMinuit->SetFCN(fcn);
550 
551 }
552 
553 ////////////////////////////////////////////////////////////////////////////////
554 /// Specify the address of the fitting algorithm
555 
557 {
558  if (fCovar) {delete [] fCovar; fCovar = 0;}
560  fMinuit->SetFCN(fcn);
561 }
562 
563 ////////////////////////////////////////////////////////////////////////////////
564 /// ret fit method (chisquare or loglikelihood)
565 
566 void TFitter::SetFitMethod(const char *name)
567 {
568  if (fCovar) {delete [] fCovar; fCovar = 0;}
569  if (!strcmp(name,"H1FitChisquare")) SetFCN(H1FitChisquare);
570  if (!strcmp(name,"H1FitLikelihood")) SetFCN(H1FitLikelihood);
571  if (!strcmp(name,"GraphFitChisquare")) SetFCN(GraphFitChisquare);
572  if (!strcmp(name,"Graph2DFitChisquare")) SetFCN(Graph2DFitChisquare);
573  if (!strcmp(name,"MultiGraphFitChisquare")) SetFCN(MultiGraphFitChisquare);
574  if (!strcmp(name,"F2Minimizer")) SetFCN(F2Fit);
575  if (!strcmp(name,"F3Minimizer")) SetFCN(F3Fit);
576 }
577 
578 ////////////////////////////////////////////////////////////////////////////////
579 /// set initial values for a parameter
580 /// ipar : parameter number
581 /// parname : parameter name
582 /// value : initial parameter value
583 /// verr : initial error for this parameter
584 /// vlow : lower value for the parameter
585 /// vhigh : upper value for the parameter
586 
587 Int_t TFitter::SetParameter(Int_t ipar,const char *parname,Double_t value,Double_t verr,Double_t vlow, Double_t vhigh)
588 {
589  if (fCovar) {delete [] fCovar; fCovar = 0;}
590  Int_t ierr = 0;
591  fMinuit->mnparm(ipar,parname,value,verr,vlow,vhigh,ierr);
592  return ierr;
593 }
594 
595 ////////////////////////////////////////////////////////////////////////////////
596 /// Minimization function for H1s using a Chisquare method
597 /// Default method (function evaluated at center of bin)
598 /// for each point the cache contains the following info
599 /// -1D : bc,e, xc (bin content, error, x of center of bin)
600 /// -2D : bc,e, xc,yc
601 /// -3D : bc,e, xc,yc,zc
602 
604 {
605  Foption_t fitOption = GetFitOption();
606  if (fitOption.Integral) {
607  FitChisquareI(npar,gin,f,u,flag);
608  return;
609  }
610  Double_t cu,eu,fu,fsum;
611  Double_t dersum[100], grad[100];
612  memset(grad,0,800);
613  Double_t x[3];
614 
615  TH1 *hfit = (TH1*)GetObjectFit();
616  TF1 *f1 = (TF1*)GetUserFunc();
617  Int_t nd = hfit->GetDimension();
618  Int_t j;
619 
620  f1->InitArgs(x,u);
621  npar = f1->GetNpar();
622  if (flag == 2) for (j=0;j<npar;j++) dersum[j] = gin[j] = 0;
623  f = 0;
624 
625  Int_t npfit = 0;
626  Double_t *cache = fCache;
627  for (Int_t i=0;i<fNpoints;i++) {
628  if (nd > 2) x[2] = cache[4];
629  if (nd > 1) x[1] = cache[3];
630  x[0] = cache[2];
631  cu = cache[0];
633  fu = f1->EvalPar(x,u);
634  if (TF1::RejectedPoint()) {cache += fPointSize; continue;}
635  eu = cache[1];
636  if (flag == 2) {
637  for (j=0;j<npar;j++) dersum[j] += 1; //should be the derivative
638  for (j=0;j<npar;j++) grad[j] += dersum[j]*(fu-cu)/eu; dersum[j] = 0;
639  }
640  fsum = (cu-fu)/eu;
641  f += fsum*fsum;
642  npfit++;
643  cache += fPointSize;
644  }
645  f1->SetNumberFitPoints(npfit);
646 }
647 
648 ////////////////////////////////////////////////////////////////////////////////
649 /// Minimization function for H1s using a Chisquare method
650 /// The "I"ntegral method is used
651 /// for each point the cache contains the following info
652 /// -1D : bc,e, xc,xw (bin content, error, x of center of bin, x bin width of bin)
653 /// -2D : bc,e, xc,xw,yc,yw
654 /// -3D : bc,e, xc,xw,yc,yw,zc,zw
655 
657 {
658  Double_t cu,eu,fu,fsum;
659  Double_t dersum[100], grad[100];
660  memset(grad,0,800);
661  Double_t x[3];
662 
663  TH1 *hfit = (TH1*)GetObjectFit();
664  TF1 *f1 = (TF1*)GetUserFunc();
665  Int_t nd = hfit->GetDimension();
666  Int_t j;
667 
668  f1->InitArgs(x,u);
669  npar = f1->GetNpar();
670  if (flag == 2) for (j=0;j<npar;j++) dersum[j] = gin[j] = 0;
671  f = 0;
672 
673  Int_t npfit = 0;
674  Double_t *cache = fCache;
675  for (Int_t i=0;i<fNpoints;i++) {
676  cu = cache[0];
678  f1->SetParameters(u);
679  if (nd < 2) {
680  fu = f1->Integral(cache[2] - 0.5*cache[3],cache[2] + 0.5*cache[3])/cache[3];
681  } else if (nd < 3) {
682  fu = ((TF2*)f1)->Integral(cache[2] - 0.5*cache[3],cache[2] + 0.5*cache[3],cache[4] - 0.5*cache[5],cache[4] + 0.5*cache[5])/(cache[3]*cache[5]);
683  } else {
684  fu = ((TF3*)f1)->Integral(cache[2] - 0.5*cache[3],cache[2] + 0.5*cache[3],cache[4] - 0.5*cache[5],cache[4] + 0.5*cache[5],cache[6] - 0.5*cache[7],cache[6] + 0.5*cache[7])/(cache[3]*cache[5]*cache[7]);
685  }
686  if (TF1::RejectedPoint()) {cache += fPointSize; continue;}
687  eu = cache[1];
688  if (flag == 2) {
689  for (j=0;j<npar;j++) dersum[j] += 1; //should be the derivative
690  for (j=0;j<npar;j++) grad[j] += dersum[j]*(fu-cu)/eu; dersum[j] = 0;
691  }
692  fsum = (cu-fu)/eu;
693  f += fsum*fsum;
694  npfit++;
695  cache += fPointSize;
696  }
697  f1->SetNumberFitPoints(npfit);
698 }
699 
700 
701 ////////////////////////////////////////////////////////////////////////////////
702 /// Minimization function for H1s using a Likelihood method*-*-*-*-*-*
703 /// Basically, it forms the likelihood by determining the Poisson
704 /// probability that given a number of entries in a particular bin,
705 /// the fit would predict it's value. This is then done for each bin,
706 /// and the sum of the logs is taken as the likelihood.
707 /// Default method (function evaluated at center of bin)
708 /// for each point the cache contains the following info
709 /// -1D : bc,e, xc (bin content, error, x of center of bin)
710 /// -2D : bc,e, xc,yc
711 /// -3D : bc,e, xc,yc,zc
712 
714 {
715  Foption_t fitOption = GetFitOption();
716  if (fitOption.Integral) {
717  FitLikelihoodI(npar,gin,f,u,flag);
718  return;
719  }
720  Double_t cu,fu,fobs,fsub;
721  Double_t dersum[100];
722  Double_t x[3];
723  Int_t icu;
724 
725  TH1 *hfit = (TH1*)GetObjectFit();
726  TF1 *f1 = (TF1*)GetUserFunc();
727  Int_t nd = hfit->GetDimension();
728  Int_t j;
729 
730  f1->InitArgs(x,u);
731  npar = f1->GetNpar();
732  if (flag == 2) for (j=0;j<npar;j++) dersum[j] = gin[j] = 0;
733  f = 0;
734 
735  Int_t npfit = 0;
736  Double_t *cache = fCache;
737  for (Int_t i=0;i<fNpoints;i++) {
738  if (nd > 2) x[2] = cache[4];
739  if (nd > 1) x[1] = cache[3];
740  x[0] = cache[2];
741  cu = cache[0];
743  fu = f1->EvalPar(x,u);
744  if (TF1::RejectedPoint()) {cache += fPointSize; continue;}
745  if (flag == 2) {
746  for (j=0;j<npar;j++) {
747  dersum[j] += 1; //should be the derivative
748  //grad[j] += dersum[j]*(fu-cu)/eu; dersum[j] = 0;
749  }
750  }
751  if (fu < 1.e-9) fu = 1.e-9;
752  if (fitOption.Like == 1) {
753  icu = Int_t(cu);
754  fsub = -fu +cu*TMath::Log(fu);
755  if (icu <10000) fobs = GetSumLog(icu);
756  else fobs = TMath::LnGamma(cu+1);
757  } else {
758  fsub = -fu +cu*TMath::Log(fu);
759  fobs = TMath::LnGamma(cu+1);
760  }
761  fsub -= fobs;
762  f -= fsub;
763  npfit++;
764  cache += fPointSize;
765  }
766  f *= 2;
767  f1->SetNumberFitPoints(npfit);
768 }
769 
770 
771 ////////////////////////////////////////////////////////////////////////////////
772 /// Minimization function for H1s using a Likelihood method*-*-*-*-*-*
773 /// Basically, it forms the likelihood by determining the Poisson
774 /// probability that given a number of entries in a particular bin,
775 /// the fit would predict it's value. This is then done for each bin,
776 /// and the sum of the logs is taken as the likelihood.
777 /// The "I"ntegral method is used
778 /// for each point the cache contains the following info
779 /// -1D : bc,e, xc,xw (bin content, error, x of center of bin, x bin width of bin)
780 /// -2D : bc,e, xc,xw,yc,yw
781 /// -3D : bc,e, xc,xw,yc,yw,zc,zw
782 
784 {
785  Double_t cu,fu,fobs,fsub;
786  Double_t dersum[100];
787  Double_t x[3];
788  Int_t icu;
789 
790  TH1 *hfit = (TH1*)GetObjectFit();
791  TF1 *f1 = (TF1*)GetUserFunc();
792  Foption_t fitOption = GetFitOption();
793  Int_t nd = hfit->GetDimension();
794  Int_t j;
795 
796  f1->InitArgs(x,u);
797  npar = f1->GetNpar();
798  if (flag == 2) for (j=0;j<npar;j++) dersum[j] = gin[j] = 0;
799  f = 0;
800 
801  Int_t npfit = 0;
802  Double_t *cache = fCache;
803  for (Int_t i=0;i<fNpoints;i++) {
804  if (nd > 2) x[2] = cache[6];
805  if (nd > 1) x[1] = cache[4];
806  x[0] = cache[2];
807  cu = cache[0];
809  f1->SetParameters(u);
810  if (nd < 2) {
811  fu = f1->Integral(cache[2] - 0.5*cache[3],cache[2] + 0.5*cache[3])/cache[3];
812  } else if (nd < 3) {
813  fu = ((TF2*)f1)->Integral(cache[2] - 0.5*cache[3],cache[2] + 0.5*cache[3],cache[4] - 0.5*cache[5],cache[4] + 0.5*cache[5])/(cache[3]*cache[5]);
814  } else {
815  fu = ((TF3*)f1)->Integral(cache[2] - 0.5*cache[3],cache[2] + 0.5*cache[3],cache[4] - 0.5*cache[5],cache[4] + 0.5*cache[5],cache[6] - 0.5*cache[7],cache[6] + 0.5*cache[7])/(cache[3]*cache[5]*cache[7]);
816  }
817  if (TF1::RejectedPoint()) {cache += fPointSize; continue;}
818  if (flag == 2) {
819  for (j=0;j<npar;j++) {
820  dersum[j] += 1; //should be the derivative
821  //grad[j] += dersum[j]*(fu-cu)/eu; dersum[j] = 0;
822  }
823  }
824  if (fu < 1.e-9) fu = 1.e-9;
825  if (fitOption.Like == 1) {
826  icu = Int_t(cu);
827  fsub = -fu +cu*TMath::Log(fu);
828  if (icu <10000) fobs = GetSumLog(icu);
829  else fobs = TMath::LnGamma(cu+1);
830  } else {
831  fsub = -fu +cu*TMath::Log(fu);
832  fobs = TMath::LnGamma(cu+1);
833  }
834  fsub -= fobs;
835  f -= fsub;
836  npfit++;
837  cache += fPointSize;
838  }
839  f *= 2;
840  f1->SetNumberFitPoints(npfit);
841 }
842 
843 
844 
845 ////////////////////////////////////////////////////////////////////////////////
846 /// Minimization function for H1s using a Chisquare method
847 /// ======================================================
848 
849 void H1FitChisquare(Int_t &npar, Double_t *gin, Double_t &f, Double_t *u, Int_t flag)
850 {
852  hFitter->FitChisquare(npar, gin, f, u, flag);
853 }
854 
855 ////////////////////////////////////////////////////////////////////////////////
856 /// -*-*-*-*Minimization function for H1s using a Likelihood method*-*-*-*-*-*
857 /// =======================================================
858 /// Basically, it forms the likelihood by determining the Poisson
859 /// probability that given a number of entries in a particular bin,
860 /// the fit would predict it's value. This is then done for each bin,
861 /// and the sum of the logs is taken as the likelihood.
862 
863 void H1FitLikelihood(Int_t &npar, Double_t *gin, Double_t &f, Double_t *u, Int_t flag)
864 {
866  hFitter->FitLikelihood(npar, gin, f, u, flag);
867 }
868 
869 ////////////////////////////////////////////////////////////////////////////////
870 ///*-*-*-*-*-*Minimization function for Graphs using a Chisquare method*-*-*-*-*
871 ///*-* =========================================================
872 ///
873 /// In case of a TGraphErrors object, ex, the error along x, is projected
874 /// along the y-direction by calculating the function at the points x-exlow and
875 /// x+exhigh.
876 ///
877 /// The chisquare is computed as the sum of the quantity below at each point:
878 ///
879 /// (y - f(x))**2
880 /// -----------------------------------
881 /// ey**2 + (0.5*(exl + exh)*f'(x))**2
882 ///
883 /// where x and y are the point coordinates and f'(x) is the derivative of function f(x).
884 /// This method to approximate the uncertainty in y because of the errors in x, is called
885 /// "effective variance" method.
886 /// The improvement, compared to the previously used method (f(x+ exhigh) - f(x-exlow))/2
887 /// is of (error of x)**2 order.
888 ///
889 /// In case the function lies below (above) the data point, ey is ey_low (ey_high).
890 ///
891 /// thanks to Andy Haas (haas@yahoo.com) for adding the case with TGraphAsymmErrors
892 /// University of Washington
893 /// February 3, 2004
894 ///
895 /// NOTE:
896 /// 1) By using the "effective variance" method a simple linear regression
897 /// becomes a non-linear case , which takes several iterations
898 /// instead of 0 as in the linear case .
899 ///
900 /// 2) The effective variance technique assumes that there is no correlation
901 /// between the x and y coordinate .
902 ///
903 /// The book by Sigmund Brandt (Data Analysis) contains an interesting
904 /// section how to solve the problem when correclations do exist .
905 
906 void GraphFitChisquare(Int_t &npar, Double_t * /*gin*/, Double_t &f,
907  Double_t *u, Int_t /*flag*/)
908 {
909  Double_t cu,eu,exh,exl,ey,eux,fu,fsum;
910  Double_t x[1];
911  //Double_t xm,xp;
912  Int_t bin, npfits=0;
913 
915  TGraph *gr = (TGraph*)grFitter->GetObjectFit();
916  TF1 *f1 = (TF1*)grFitter->GetUserFunc();
917  Foption_t fitOption = grFitter->GetFitOption();
918  Int_t n = gr->GetN();
919  Double_t *gx = gr->GetX();
920  Double_t *gy = gr->GetY();
921  //Double_t fxmin = f1->GetXmin();
922  //Double_t fxmax = f1->GetXmax();
923  npar = f1->GetNpar();
924 
925  f = 0;
926  for (bin=0;bin<n;bin++) {
927  f1->InitArgs(x,u); //must be inside the loop because of TF1::Derivative calling InitArgs
928  x[0] = gx[bin];
929  if (!f1->IsInside(x)) continue;
930  cu = gy[bin];
932  fu = f1->EvalPar(x,u);
933  if (TF1::RejectedPoint()) continue;
934  fsum = (cu-fu);
935  npfits++;
936  if (fitOption.W1) {
937  f += fsum*fsum;
938  continue;
939  }
940 
941  exh = gr->GetErrorXhigh(bin);
942  exl = gr->GetErrorXlow(bin);
943  if (fsum < 0)
944  ey = gr->GetErrorYhigh(bin);
945  else
946  ey = gr->GetErrorYlow(bin);
947  if (exl < 0) exl = 0;
948  if (exh < 0) exh = 0;
949  if (ey < 0) ey = 0;
950  if (exh > 0 || exl > 0) {
951  //Without the "variance method", we had the 6 next lines instead
952  // of the line above.
953  //xm = x[0] - exl; if (xm < fxmin) xm = fxmin;
954  //xp = x[0] + exh; if (xp > fxmax) xp = fxmax;
955  //Double_t fm,fp;
956  //x[0] = xm; fm = f1->EvalPar(x,u);
957  //x[0] = xp; fp = f1->EvalPar(x,u);
958  //eux = 0.5*(fp-fm);
959 
960  //"Effective Variance" method introduced by Anna Kreshuk
961  // in version 4.00/08.
962  eux = 0.5*(exl + exh)*f1->Derivative(x[0], u);
963  } else
964  eux = 0.;
965  eu = ey*ey+eux*eux;
966  if (eu <= 0) eu = 1;
967  f += fsum*fsum/eu;
968  }
969  f1->SetNumberFitPoints(npfits);
970  //
971 }
972 
973 
974 ////////////////////////////////////////////////////////////////////////////////
975 ///*-*-*-*-*Minimization function for 2D Graphs using a Chisquare method*-*-*-*-*
976 ///*-* ============================================================
977 
978 void Graph2DFitChisquare(Int_t &npar, Double_t * /*gin*/, Double_t &f,
979  Double_t *u, Int_t /*flag*/)
980 {
981  Double_t cu,eu,ex,ey,ez,eux,euy,fu,fsum,fm,fp;
982  Double_t x[2];
983  Double_t xm,xp,ym,yp;
984  Int_t bin, npfits=0;
985 
987  TGraph2D *gr = (TGraph2D*)grFitter->GetObjectFit();
988  TF2 *f2 = (TF2*)grFitter->GetUserFunc();
989  Foption_t fitOption = grFitter->GetFitOption();
990 
991  Int_t n = gr->GetN();
992  Double_t *gx = gr->GetX();
993  Double_t *gy = gr->GetY();
994  Double_t *gz = gr->GetZ();
995  Double_t fxmin = f2->GetXmin();
996  Double_t fxmax = f2->GetXmax();
997  Double_t fymin = f2->GetYmin();
998  Double_t fymax = f2->GetYmax();
999  npar = f2->GetNpar();
1000  f = 0;
1001  for (bin=0;bin<n;bin++) {
1002  f2->InitArgs(x,u);
1003  x[0] = gx[bin];
1004  x[1] = gy[bin];
1005  if (!f2->IsInside(x)) continue;
1006  cu = gz[bin];
1008  fu = f2->EvalPar(x,u);
1009  if (TF2::RejectedPoint()) continue;
1010  fsum = (cu-fu);
1011  npfits++;
1012  if (fitOption.W1) {
1013  f += fsum*fsum;
1014  continue;
1015  }
1016  ex = gr->GetErrorX(bin);
1017  ey = gr->GetErrorY(bin);
1018  ez = gr->GetErrorZ(bin);
1019  if (ex < 0) ex = 0;
1020  if (ey < 0) ey = 0;
1021  if (ez < 0) ez = 0;
1022  eux = euy = 0;
1023  if (ex > 0) {
1024  xm = x[0] - ex; if (xm < fxmin) xm = fxmin;
1025  xp = x[0] + ex; if (xp > fxmax) xp = fxmax;
1026  x[0] = xm; fm = f2->EvalPar(x,u);
1027  x[0] = xp; fp = f2->EvalPar(x,u);
1028  eux = fp-fm;
1029  }
1030  if (ey > 0) {
1031  x[0] = gx[bin];
1032  ym = x[1] - ey; if (ym < fymin) ym = fxmin;
1033  yp = x[1] + ey; if (yp > fymax) yp = fymax;
1034  x[1] = ym; fm = f2->EvalPar(x,u);
1035  x[1] = yp; fp = f2->EvalPar(x,u);
1036  euy = fp-fm;
1037  }
1038  eu = ez*ez+eux*eux+euy*euy;
1039  if (eu <= 0) eu = 1;
1040  f += fsum*fsum/eu;
1041  }
1042  f2->SetNumberFitPoints(npfits);
1043 }
1044 
1045 ////////////////////////////////////////////////////////////////////////////////
1046 
1048  Double_t *u, Int_t /*flag*/)
1049 {
1050  Double_t cu,eu,exh,exl,ey,eux,fu,fsum;
1051  Double_t x[1];
1052  //Double_t xm,xp;
1053  Int_t bin, npfits=0;
1054 
1056  TMultiGraph *mg = (TMultiGraph*)grFitter->GetObjectFit();
1057  TF1 *f1 = (TF1*)grFitter->GetUserFunc();
1058  Foption_t fitOption = grFitter->GetFitOption();
1059  TGraph *gr;
1060  TIter next(mg->GetListOfGraphs());
1061 
1062  Int_t n;
1063  Double_t *gx;
1064  Double_t *gy;
1065  //Double_t fxmin = f1->GetXmin();
1066  //Double_t fxmax = f1->GetXmax();
1067  npar = f1->GetNpar();
1068 
1069  f = 0;
1070 
1071  while ((gr = (TGraph*) next())) {
1072  n = gr->GetN();
1073  gx = gr->GetX();
1074  gy = gr->GetY();
1075  for (bin=0;bin<n;bin++) {
1076  f1->InitArgs(x,u); //must be inside the loop because of TF1::Derivative calling InitArgs
1077  x[0] = gx[bin];
1078  if (!f1->IsInside(x)) continue;
1079  cu = gy[bin];
1081  fu = f1->EvalPar(x,u);
1082  if (TF1::RejectedPoint()) continue;
1083  fsum = (cu-fu);
1084  npfits++;
1085  if (fitOption.W1) {
1086  f += fsum*fsum;
1087  continue;
1088  }
1089  exh = gr->GetErrorXhigh(bin);
1090  exl = gr->GetErrorXlow(bin);
1091  ey = gr->GetErrorY(bin);
1092  if (exl < 0) exl = 0;
1093  if (exh < 0) exh = 0;
1094  if (ey < 0) ey = 0;
1095  if (exh > 0 && exl > 0) {
1096  //Without the "variance method", we had the 6 next lines instead
1097  // of the line above.
1098  //xm = x[0] - exl; if (xm < fxmin) xm = fxmin;
1099  //xp = x[0] + exh; if (xp > fxmax) xp = fxmax;
1100  //Double_t fm,fp;
1101  //x[0] = xm; fm = f1->EvalPar(x,u);
1102  //x[0] = xp; fp = f1->EvalPar(x,u);
1103  //eux = 0.5*(fp-fm);
1104 
1105  //"Effective Variance" method introduced by Anna Kreshuk
1106  // in version 4.00/08.
1107  eux = 0.5*(exl + exh)*f1->Derivative(x[0], u);
1108  } else
1109  eux = 0.;
1110  eu = ey*ey+eux*eux;
1111  if (eu <= 0) eu = 1;
1112  f += fsum*fsum/eu;
1113  }
1114  }
1115  f1->SetNumberFitPoints(npfits);
1116 }
1117 
1118 ////////////////////////////////////////////////////////////////////////////////
1119 
1120 void F2Fit(Int_t &/*npar*/, Double_t * /*gin*/, Double_t &f,Double_t *u, Int_t /*flag*/)
1121 {
1123  TF2 *f2 = (TF2*)fitter->GetObjectFit();
1124  f2->InitArgs(u, f2->GetParameters() );
1125  f = f2->EvalPar(u);
1126 }
1127 
1128 ////////////////////////////////////////////////////////////////////////////////
1129 
1130 void F3Fit(Int_t &/*npar*/, Double_t * /*gin*/, Double_t &f,Double_t *u, Int_t /*flag*/)
1131 {
1133  TF3 *f3 = (TF3*)fitter->GetObjectFit();
1134  f3->InitArgs(u, f3->GetParameters() );
1135  f = f3->EvalPar(u);
1136 }
void GraphFitChisquare(Int_t &npar, Double_t *gin, Double_t &f, Double_t *u, Int_t flag)
*-*-*-*-*-*Minimization function for Graphs using a Chisquare method*-*-*-*-* *-* ===================...
Definition: TFitter.cxx:906
Int_t GetFirst() const
Return first bin on the axis i.e.
Definition: TAxis.cxx:429
virtual void mnstat(Double_t &fmin, Double_t &fedm, Double_t &errdef, Int_t &npari, Int_t &nparx, Int_t &istat)
*-*-*-*-*Returns concerning the current status of the minimization*-*-*-*-* *-* =====================...
Definition: TMinuit.cxx:7699
Double_t * fCache
virtual void mnprin(Int_t inkode, Double_t fval)
*-*-*-*Prints the values of the parameters at the time of the call*-*-*-*-* *-* =====================...
Definition: TMinuit.cxx:6373
virtual void FitChisquareI(Int_t &npar, Double_t *gin, Double_t &f, Double_t *u, Int_t flag)
Minimization function for H1s using a Chisquare method The "I"ntegral method is used for each point t...
Definition: TFitter.cxx:656
virtual void SetParameters(const Double_t *params)
Definition: TF1.h:432
ClassImp(TFitter) TFitter
*-*-*-*-*-*-*-*-*-*-*default constructor*-*-*-*-*-*-*-*-*-*-*-*-* *-* =================== ...
Definition: TFitter.cxx:31
void grad()
Definition: grad.C:17
static float fu(float x)
Definition: main.cpp:53
virtual Double_t GetErrorYhigh(Int_t bin) const
This function is called by GraphFitChisquare.
Definition: TGraph.cxx:1410
void fcn(Int_t &npar, Double_t *gin, Double_t &f, Double_t *par, Int_t iflag)
Definition: Ifit.C:26
Implementation in C++ of the Minuit package written by F.
Definition: TMinuit.h:34
Double_t Log(Double_t x)
Definition: TMath.h:526
virtual Bool_t InheritsFrom(const char *classname) const
Returns kTRUE if object inherits from class "classname".
Definition: TObject.cxx:487
virtual void FitLikelihood(Int_t &npar, Double_t *gin, Double_t &f, Double_t *u, Int_t flag)
Minimization function for H1s using a Likelihood method*-*-*-*-*-* Basically, it forms the likelihood...
Definition: TFitter.cxx:713
virtual Double_t GetYmax() const
Definition: TF2.h:120
ClassImp(TSeqCollection) Int_t TSeqCollection TIter next(this)
Return index of object in collection.
virtual void mnemat(Double_t *emat, Int_t ndim)
Calculates the external error matrix from the internal matrix.
Definition: TMinuit.cxx:2523
TList * GetListOfGraphs() const
Definition: TMultiGraph.h:71
virtual Int_t GetNumberFreeParameters() const
Return the number of free parameters.
Definition: TF1.cxx:1579
return c
const char Option_t
Definition: RtypesCore.h:62
virtual Foption_t GetFitOption() const
Double_t StudentQuantile(Double_t p, Double_t ndf, Bool_t lower_tail=kTRUE)
Computes quantiles of the Student's t-distribution 1st argument is the probability, at which the quantile is computed 2nd argument - the number of degrees of freedom of the Student distribution When the 3rd argument lower_tail is kTRUE (default)- the algorithm returns such x0, that P(x < x0)=p upper tail (lower_tail is kFALSE)- the algorithm returns such x0, that P(x > x0)=p the algorithm was taken from G.W.Hill, "Algorithm 396, Student's t-quantiles" "Communications of the ACM", 13(10), October 1970.
Definition: TMath.cxx:2553
virtual Int_t GetDimension() const
Definition: TH1.h:283
Double_t * GetX() const
Definition: TGraph2D.h:128
tuple f2
Definition: surfaces.py:24
virtual void FitChisquare(Int_t &npar, Double_t *gin, Double_t &f, Double_t *u, Int_t flag)
Minimization function for H1s using a Chisquare method Default method (function evaluated at center o...
Definition: TFitter.cxx:603
virtual void PrintResults(Int_t level, Double_t amin) const
Print fit results.
Definition: TFitter.cxx:528
tuple pname
Definition: tree.py:131
void f3()
Definition: na49.C:50
A TMultiGraph is a collection of TGraph (or derived) objects.
Definition: TMultiGraph.h:37
virtual Double_t GetParError(Int_t ipar) const
return error of parameter ipar
Definition: TFitter.cxx:426
Basic string class.
Definition: TString.h:137
virtual TObject * GetUserFunc() const
virtual void SetNumberFitPoints(Int_t npfits)
Definition: TF1.h:421
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
const Bool_t kFALSE
Definition: Rtypes.h:92
virtual Int_t GetNumberTotalParameters() const
return the total number of parameters (free + fixed)
Definition: TFitter.cxx:409
virtual Int_t ExecuteCommand(const char *command, Double_t *args, Int_t nargs)
Execute a fitter command; command : command string args : list of nargs command arguments.
Definition: TFitter.cxx:86
virtual void SetFCN(void *fcn)
To set the address of the minimization objective function.
virtual Double_t Integral(Double_t a, Double_t b, Double_t epsrel=1.e-12)
IntegralOneDim or analytical integral.
Definition: TF1.cxx:2241
Int_t GetN() const
Definition: TGraph.h:132
virtual Int_t GetNDF() const
Return the number of degrees of freedom in the fit the fNDF parameter has been previously computed du...
Definition: TF1.cxx:1568
virtual Double_t GetXmin() const
Definition: TF1.h:381
TFile * f
virtual Int_t GetNumPars() const
returns the total number of parameters that have been defined as fixed or free.
Definition: TMinuit.cxx:866
virtual Double_t * GetEY() const
Definition: TGraph.h:142
virtual Double_t GetErrorXhigh(Int_t bin) const
This function is called by GraphFitChisquare.
Definition: TGraph.cxx:1390
virtual Int_t FixParameter(Int_t parNo)
fix a parameter
Definition: TMinuit.cxx:821
virtual ~TFitter()
*-*-*-*-*-*-*-*-*-*-*default destructor*-*-*-*-*-*-*-*-*-*-*-*-*-* *-* ================== ...
Definition: TFitter.cxx:50
virtual Double_t GetErrorY(Int_t bin) const
This function is called by Graph2DFitChisquare.
Definition: TGraph2D.cxx:1024
const char * Data() const
Definition: TString.h:349
Double_t * GetY() const
Definition: TGraph.h:140
virtual void SetFCN(void *fcn)
*-*-*-*-*-*-*To set the address of the minimization function*-*-*-*-*-*-*-* *-* =====================...
Definition: TMinuit.cxx:947
virtual Double_t GetParameter(Int_t ipar) const
return current value of parameter ipar
Definition: TFitter.cxx:440
Double_t x[n]
Definition: legend1.C:17
Double_t GetChisquare() const
Definition: TF1.h:329
int GetDimension(const TH1 *h1)
Definition: HFitImpl.cxx:50
Int_t fNpar
Definition: TMinuit.h:48
virtual void ReleaseParameter(Int_t ipar)
Release parameter ipar.
Definition: TFitter.cxx:536
void Class()
Definition: Class.C:29
TObject * fObjectFit
virtual Int_t GetErrors(Int_t ipar, Double_t &eplus, Double_t &eminus, Double_t &eparab, Double_t &globcc) const
return current errors for a parameter ipar : parameter number eplus : upper error eminus : lower erro...
Definition: TFitter.cxx:397
virtual void mnrn15(Double_t &val, Int_t &inseed)
*-*-*-*-*-*-*This is a super-portable random number generator*-*-*-*-*-*-* *-* ======================...
Definition: TMinuit.cxx:6688
virtual Int_t GetStats(Double_t &amin, Double_t &edm, Double_t &errdef, Int_t &nvpar, Int_t &nparx) const
return global fit parameters amin : chisquare edm : estimated distance to minimum errdef nvpar : numb...
Definition: TFitter.cxx:486
virtual Double_t Chisquare(Int_t npar, Double_t *params) const
return a chisquare equivalent
Definition: TFitter.cxx:60
tuple np
Definition: multifit.py:30
Double_t * fSumLog
Definition: TFitter.h:35
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition: TObject.cxx:918
virtual void SetBinError(Int_t bin, Double_t error)
see convention for numbering bins in TH1::GetBin
Definition: TH1.cxx:8528
virtual Double_t * GetEZ() const
Definition: TGraph2D.h:133
void H1FitChisquare(Int_t &npar, Double_t *gin, Double_t &f, Double_t *u, Int_t flag)
Minimization function for H1s using a Chisquare method
Definition: TFitter.cxx:849
virtual Double_t GetXmax() const
Definition: TF1.h:382
virtual void GetParLimits(Int_t ipar, Double_t &parmin, Double_t &parmax) const
Return limits for parameter ipar.
Definition: TF1.cxx:1618
virtual Double_t GradientPar(Int_t ipar, const Double_t *x, Double_t eps=0.01)
Compute the gradient (derivative) wrt a parameter ipar.
Definition: TF1.cxx:2117
virtual Double_t GetErrorXlow(Int_t bin) const
This function is called by GraphFitChisquare.
Definition: TGraph.cxx:1400
TObject * fUserFunc
virtual Double_t GetCovarianceMatrixElement(Int_t i, Int_t j) const
return element i,j from the covariance matrix
Definition: TFitter.cxx:378
virtual Double_t GetErrorX(Int_t bin) const
This function is called by Graph2DFitChisquare.
Definition: TGraph2D.cxx:1014
TThread * t[5]
Definition: threadsh1.C:13
Double_t * GetX() const
Definition: TGraph.h:139
Int_t fNpfix
Definition: TMinuit.h:44
TString * fCpnam
Character to be plotted at the X,Y contour positions.
Definition: TMinuit.h:172
Class to manage histogram axis.
Definition: TAxis.h:36
XPoint xy[kMAXMK]
Definition: TGX11.cxx:122
virtual void mncler()
*-*-*-*-*-*-*-*-*-*-*Resets the parameter list to UNDEFINED*-*-*-*-*-*-*-* *-* ======================...
Definition: TMinuit.cxx:1121
virtual void mnparm(Int_t k, TString cnamj, Double_t uk, Double_t wk, Double_t a, Double_t b, Int_t &ierflg)
*-*-*-*-*-*-*-*-*Implements one parameter definition*-*-*-*-*-*-*-*-*-*-*-* *-* =====================...
Definition: TMinuit.cxx:5728
A 3-Dim function with parameters.
Definition: TF3.h:30
virtual void SetBinContent(Int_t bin, Double_t content)
Set bin content see convention for numbering bins in TH1::GetBin In case the bin number is greater th...
Definition: TH1.cxx:8543
Double_t * GetZ() const
Definition: TGraph2D.h:130
Double_t * fCovar
Definition: TFitter.h:34
static TVirtualFitter * GetFitter()
static: return the current Fitter
static void RejectPoint(Bool_t reject=kTRUE)
Static function to set the global flag to reject points the fgRejectPoint global flag is tested by al...
Definition: TF1.cxx:3367
int Like
Definition: Foption.h:34
TAxis * GetYaxis()
Definition: TH1.h:320
int Integral
Definition: Foption.h:44
virtual Bool_t IsFixed(Int_t ipar) const
return kTRUE if parameter ipar is fixed, kFALSE othersise)
Definition: TFitter.cxx:518
A 2-Dim function with parameters.
Definition: TF2.h:33
virtual void SetFCN(void *fcn)
Specify the address of the fitting algorithm (from the interpreter)
Definition: TFitter.cxx:545
virtual void mnerrs(Int_t number, Double_t &eplus, Double_t &eminus, Double_t &eparab, Double_t &gcc)
*-*-*-*-*-*-*-*-*-*Utility routine to get MINOS errors*-*-*-*-*-*-*-*-*-*-* *-* =====================...
Definition: TMinuit.cxx:2599
TGraphErrors * gr
Definition: legend1.C:25
int W1
Definition: Foption.h:36
Int_t * fNiofex
Definition: TMinuit.h:134
void MultiGraphFitChisquare(Int_t &npar, Double_t *gin, Double_t &f, Double_t *u, Int_t flag)
Definition: TFitter.cxx:1047
virtual const char * GetParName(Int_t ipar) const
return name of parameter ipar
Definition: TFitter.cxx:472
virtual Double_t Derivative(Double_t x, Double_t *params=0, Double_t epsilon=0.001) const
Returns the first derivative of the function at point x, computed by Richardson's extrapolation metho...
Definition: TF1.cxx:829
virtual void SetFitMethod(const char *name)
ret fit method (chisquare or loglikelihood)
Definition: TFitter.cxx:566
virtual Double_t GetErrorYlow(Int_t bin) const
This function is called by GraphFitChisquare.
Definition: TGraph.cxx:1420
virtual Double_t GetYmin() const
Definition: TF2.h:119
virtual Int_t GetNumberFreeParameters() const
return the number of free parameters
Definition: TFitter.cxx:417
virtual Double_t GetSumLog(Int_t i)
return Sum(log(i) i=0,n used by log likelihood fits
Definition: TFitter.cxx:497
virtual void SetPoint(Int_t point, Double_t x, Double_t y, Double_t z)
Sets point number n.
Definition: TGraph2D.cxx:1707
double Double_t
Definition: RtypesCore.h:55
void F2Fit(Int_t &npar, Double_t *gin, Double_t &f, Double_t *u, Int_t flag)
Definition: TFitter.cxx:1120
virtual Int_t Release(Int_t parNo)
release a parameter
Definition: TMinuit.cxx:888
virtual Double_t * GetCovarianceMatrix() const
return a pointer to the covariance matrix
Definition: TFitter.cxx:366
virtual void GetConfidenceIntervals(Int_t n, Int_t ndim, const Double_t *x, Double_t *ci, Double_t cl=0.95)
Computes point-by-point confidence intervals for the fitted function Parameters: n - number of points...
Definition: TFitter.cxx:114
Int_t fNlog
Definition: TFitter.h:33
Double_t y[n]
Definition: legend1.C:17
Double_t ey[n]
Definition: legend1.C:17
The TH1 histogram class.
Definition: TH1.h:80
virtual void mnexcm(const char *comand, Double_t *plist, Int_t llist, Int_t &ierflg)
*-*-*-*-*-*Interprets a command and takes appropriate action*-*-*-*-*-*-*-* *-* =====================...
Definition: TMinuit.cxx:2688
virtual void FixParameter(Int_t ipar)
Fix parameter ipar.
Definition: TFitter.cxx:97
static const double eu
Definition: Vavilov.cxx:52
TMinuit * fMinuit
Definition: TFitter.h:36
virtual void InitArgs(const Double_t *x, const Double_t *params)
Initialize parameters addresses.
Definition: TF1.cxx:2198
#define name(a, b)
Definition: linkTestLib0.cpp:5
Abstract Base Class for Fitting.
Int_t GetLast() const
Return last bin on the axis i.e.
Definition: TAxis.cxx:440
TAxis * GetZaxis()
Definition: TH1.h:321
virtual Double_t GetBinCenter(Int_t bin) const
Return center of bin.
Definition: TAxis.cxx:449
Mother of all ROOT objects.
Definition: TObject.h:58
virtual Double_t GetErrorZ(Int_t bin) const
This function is called by Graph2DFitChisquare.
Definition: TGraph2D.cxx:1034
virtual Double_t * GetParameters() const
Definition: TF1.h:358
virtual Double_t GetErrorY(Int_t bin) const
This function is called by GraphFitChisquare.
Definition: TGraph.cxx:1380
void Graph2DFitChisquare(Int_t &npar, Double_t *gin, Double_t &f, Double_t *u, Int_t flag)
*-*-*-*-*Minimization function for 2D Graphs using a Chisquare method*-*-*-*-* *-* ==================...
Definition: TFitter.cxx:978
virtual void SetPoint(Int_t i, Double_t x, Double_t y)
Set x and y values for point number i.
Definition: TGraph.cxx:2127
static Bool_t RejectedPoint()
See TF1::RejectPoint above.
Definition: TF1.cxx:3376
virtual void Clear(Option_t *option="")
reset the fitter environment
Definition: TFitter.cxx:70
Int_t GetN() const
Definition: TGraph2D.h:127
1-Dim function class
Definition: TF1.h:149
virtual TObject * GetObjectFit() const
A Graph is a graphics object made of two arrays X and Y with npoints each.
Definition: TGraph.h:53
TF1 * f1
Definition: legend1.C:11
Double_t * GetY() const
Definition: TGraph2D.h:129
Double_t LnGamma(Double_t z)
Computation of ln[gamma(z)] for all z.
Definition: TMath.cxx:490
void H1FitLikelihood(Int_t &npar, Double_t *gin, Double_t &f, Double_t *u, Int_t flag)
-*-*-*-*Minimization function for H1s using a Likelihood method*-*-*-*-*-* Basically, it forms the likelihood by determining the Poisson probability that given a number of entries in a particular bin, the fit would predict it's value.
Definition: TFitter.cxx:863
virtual Bool_t IsInside(const Double_t *x) const
Return kTRUE is the point is inside the function range.
Definition: TF2.cxx:640
Double_t Sqrt(Double_t x)
Definition: TMath.h:464
void F3Fit(Int_t &npar, Double_t *gin, Double_t &f, Double_t *u, Int_t flag)
Definition: TFitter.cxx:1130
Graphics object made of three arrays X, Y and Z with the same number of points each.
Definition: TGraph2D.h:50
virtual Bool_t IsInside(const Double_t *x) const
return kTRUE if the point is inside the function range
Definition: TF1.h:404
const Bool_t kTRUE
Definition: Rtypes.h:91
virtual Double_t EvalPar(const Double_t *x, const Double_t *params=0)
Evaluate function with given coordinates and parameters.
Definition: TF1.cxx:1192
TObject * obj
float value
Definition: math.cpp:443
Int_t fNu
Definition: TMinuit.h:137
Double_t ex[n]
Definition: legend1.C:17
const Int_t n
Definition: legend1.C:16
virtual Int_t GetNpar() const
Definition: TF1.h:342
virtual void mnpout(Int_t iuext, TString &chnam, Double_t &val, Double_t &err, Double_t &xlolim, Double_t &xuplim, Int_t &iuint) const
*-*-*-*Provides the user with information concerning the current status*-*-* *-* ====================...
Definition: TMinuit.cxx:6317
gr SetName("gr")
virtual void FitLikelihoodI(Int_t &npar, Double_t *gin, Double_t &f, Double_t *u, Int_t flag)
Minimization function for H1s using a Likelihood method*-*-*-*-*-* Basically, it forms the likelihood...
Definition: TFitter.cxx:783
Int_t npfits
Definition: fit2dHist.C:46
TAxis * GetXaxis()
Definition: TH1.h:319
virtual Int_t SetParameter(Int_t ipar, const char *parname, Double_t value, Double_t verr, Double_t vlow, Double_t vhigh)
set initial values for a parameter ipar : parameter number parname : parameter name value : initial p...
Definition: TFitter.cxx:587