ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
TwoHistoFit2D.C
Go to the documentation of this file.
1 //+ Example to fit two histograms at the same time
2 //Author: Rene Brun
3 
4 #include "TH2D.h"
5 #include "TF2.h"
6 #include "TCanvas.h"
7 #include "TStyle.h"
8 #include "TRandom3.h"
9 #include "TVirtualFitter.h"
10 #include "TList.h"
11 
12 #include <vector>
13 #include <map>
14 #include <iostream>
15 
16 double gauss2D(double *x, double *par) {
17  double z1 = double((x[0]-par[1])/par[2]);
18  double z2 = double((x[1]-par[3])/par[4]);
19  return par[0]*exp(-0.5*(z1*z1+z2*z2));
20 }
21 double my2Dfunc(double *x, double *par) {
22  double *p1 = &par[0];
23  double *p2 = &par[5];
24  return gauss2D(x,p1) + gauss2D(x,p2);
25 }
26 
27 
28 
29 // data need to be globals to be visible by fcn
30 
31 std::vector<std::pair<double, double> > coords;
32 std::vector<double > values;
33 std::vector<double > errors;
34 
35 void myFcn(Int_t & /*nPar*/, Double_t * /*grad*/ , Double_t &fval, Double_t *p, Int_t /*iflag */ )
36 {
37  int n = coords.size();
38  double chi2 = 0;
39  double tmp,x[2];
40  for (int i = 0; i <n; ++i ) {
41  x[0] = coords[i].first;
42  x[1] = coords[i].second;
43  tmp = ( values[i] - my2Dfunc(x,p))/errors[i];
44  chi2 += tmp*tmp;
45  }
46  fval = chi2;
47 }
49 void FillHisto(TH2D * h, int n, double * p) {
50 
51 
52  const double mx1 = p[1];
53  const double my1 = p[3];
54  const double sx1 = p[2];
55  const double sy1 = p[4];
56  const double mx2 = p[6];
57  const double my2 = p[8];
58  const double sx2 = p[7];
59  const double sy2 = p[9];
60  //const double w1 = p[0]*sx1*sy1/(p[5]*sx2*sy2);
61  const double w1 = 0.5;
62 
63  double x, y;
64  for (int i = 0; i < n; ++i) {
65  // generate randoms with larger gaussians
66  rndm.Rannor(x,y);
67 
68  double r = rndm.Rndm(1);
69  if (r < w1) {
70  x = x*sx1 + mx1;
71  y = y*sy1 + my1;
72  }
73  else {
74  x = x*sx2 + mx2;
75  y = y*sy2 + my2;
76  }
77  h->Fill(x,y);
78 
79  }
80 }
81 
82 
83 
84 
85 int TwoHistoFit2D(bool global = true) {
86 
87  // create two histograms
88 
89  int nbx1 = 50;
90  int nby1 = 50;
91  int nbx2 = 50;
92  int nby2 = 50;
93  double xlow1 = 0.;
94  double ylow1 = 0.;
95  double xup1 = 10.;
96  double yup1 = 10.;
97  double xlow2 = 5.;
98  double ylow2 = 5.;
99  double xup2 = 20.;
100  double yup2 = 20.;
101 
102  TH2D * h1 = new TH2D("h1","core",nbx1,xlow1,xup1,nby1,ylow1,yup1);
103  TH2D * h2 = new TH2D("h2","tails",nbx2,xlow2,xup2,nby2,ylow2,yup2);
104 
105  double iniParams[10] = { 100, 6., 2., 7., 3, 100, 12., 3., 11., 2. };
106  // create fit function
107  TF2 * func = new TF2("func",my2Dfunc,xlow2,xup2,ylow2,yup2, 10);
108  func->SetParameters(iniParams);
109 
110  // fill Histos
111  int n1 = 1000000;
112  int n2 = 1000000;
113  // h1->FillRandom("func", n1);
114  //h2->FillRandom("func",n2);
115  FillHisto(h1,n1,iniParams);
116  FillHisto(h2,n2,iniParams);
117 
118  // scale histograms to same heights (for fitting)
119  double dx1 = (xup1-xlow1)/double(nbx1);
120  double dy1 = (yup1-ylow1)/double(nby1);
121  double dx2 = (xup2-xlow2)/double(nbx2);
122  double dy2 = (yup2-ylow2)/double(nby2);
123 // h1->Sumw2();
124 // h1->Scale( 1.0 / ( n1 * dx1 * dy1 ) );
125  // scale histo 2 to scale of 1
126  h2->Sumw2();
127  h2->Scale( ( double(n1) * dx1 * dy1 ) / ( double(n2) * dx2 * dy2 ) );
128 
129 
130  if (global) {
131  // fill data structure for fit (coordinates + values + errors)
132  std::cout << "Do global fit" << std::endl;
133  // fit now all the function together
134 
135  // fill data structure for fit (coordinates + values + errors)
136  TAxis *xaxis1 = h1->GetXaxis();
137  TAxis *yaxis1 = h1->GetYaxis();
138  TAxis *xaxis2 = h2->GetXaxis();
139  TAxis *yaxis2 = h2->GetYaxis();
140 
141  int nbinX1 = h1->GetNbinsX();
142  int nbinY1 = h1->GetNbinsY();
143  int nbinX2 = h2->GetNbinsX();
144  int nbinY2 = h2->GetNbinsY();
145 
146  /// reset data structure
147  coords = std::vector<std::pair<double,double> >();
148  values = std::vector<double>();
149  errors = std::vector<double>();
150 
151 
152  for (int ix = 1; ix <= nbinX1; ++ix) {
153  for (int iy = 1; iy <= nbinY1; ++iy) {
154  if ( h1->GetBinContent(ix,iy) > 0 ) {
155  coords.push_back( std::make_pair(xaxis1->GetBinCenter(ix), yaxis1->GetBinCenter(iy) ) );
156  values.push_back( h1->GetBinContent(ix,iy) );
157  errors.push_back( h1->GetBinError(ix,iy) );
158  }
159  }
160  }
161  for (int ix = 1; ix <= nbinX2; ++ix) {
162  for (int iy = 1; iy <= nbinY2; ++iy) {
163  if ( h2->GetBinContent(ix,iy) > 0 ) {
164  coords.push_back( std::make_pair(xaxis2->GetBinCenter(ix), yaxis2->GetBinCenter(iy) ) );
165  values.push_back( h2->GetBinContent(ix,iy) );
166  errors.push_back( h2->GetBinError(ix,iy) );
167  }
168  }
169  }
170 
172  TVirtualFitter * minuit = TVirtualFitter::Fitter(0,10);
173  for (int i = 0; i < 10; ++i) {
174  minuit->SetParameter(i, func->GetParName(i), func->GetParameter(i), 0.01, 0,0);
175  }
176  minuit->SetFCN(myFcn);
177 
178  double arglist[100];
179  arglist[0] = 0;
180  // set print level
181  minuit->ExecuteCommand("SET PRINT",arglist,2);
182 
183  // minimize
184  arglist[0] = 5000; // number of function calls
185  arglist[1] = 0.01; // tolerance
186  minuit->ExecuteCommand("MIGRAD",arglist,2);
187 
188  //get result
189  double minParams[10];
190  double parErrors[10];
191  for (int i = 0; i < 10; ++i) {
192  minParams[i] = minuit->GetParameter(i);
193  parErrors[i] = minuit->GetParError(i);
194  }
195  double chi2, edm, errdef;
196  int nvpar, nparx;
197  minuit->GetStats(chi2,edm,errdef,nvpar,nparx);
198 
199  func->SetParameters(minParams);
200  func->SetParErrors(parErrors);
201  func->SetChisquare(chi2);
202  int ndf = coords.size()-nvpar;
203  func->SetNDF(ndf);
204 
205  std::cout << "Chi2 Fit = " << chi2 << " ndf = " << ndf << " " << func->GetNDF() << std::endl;
206 
207  // add to list of functions
208  h1->GetListOfFunctions()->Add(func);
209  h2->GetListOfFunctions()->Add(func);
210  }
211  else {
212  // fit independently
213  h1->Fit(func);
214  h2->Fit(func);
215  }
216 
217 
218 
219  // Create a new canvas.
220  TCanvas * c1 = new TCanvas("c1","Two HIstogram Fit example",100,10,900,800);
221  c1->Divide(2,2);
222  gStyle->SetOptFit();
223  gStyle->SetStatY(0.6);
224 
225  c1->cd(1);
226  h1->Draw();
227  func->SetRange(xlow1,ylow1,xup1,yup1);
228  func->DrawCopy("cont1 same");
229  c1->cd(2);
230  h1->Draw("lego");
231  func->DrawCopy("surf1 same");
232  c1->cd(3);
233  func->SetRange(xlow2,ylow2,xup2,yup2);
234  h2->Draw();
235  func->DrawCopy("cont1 same");
236  c1->cd(4);
237  h2->Draw("lego");
238  gPad->SetLogz();
239  func->Draw("surf1 same");
240 
241  return 0;
242 }
double par[1]
Definition: unuranDistr.cxx:38
virtual void Scale(Double_t c1=1, Option_t *option="")
Multiply this histogram by a constant c1.
Definition: TH1.cxx:6174
TList * GetListOfFunctions() const
Definition: TH1.h:244
virtual void SetParameters(const Double_t *params)
Definition: TF1.h:432
virtual void Rannor(Float_t &a, Float_t &b)
Return 2 numbers distributed following a gaussian with mean=0 and sigma=1.
Definition: TRandom.cxx:460
Random number generator class based on M.
Definition: TRandom3.h:29
virtual void Draw(Option_t *option="")
Draw this function with its current attributes.
Definition: TF2.cxx:216
TCanvas * c1
Definition: legend1.C:2
virtual Int_t SetParameter(Int_t ipar, const char *parname, Double_t value, Double_t verr, Double_t vlow, Double_t vhigh)=0
virtual const char * GetParName(Int_t ipar) const
Definition: TF1.h:363
R__EXTERN TStyle * gStyle
Definition: TStyle.h:423
TH1 * h
Definition: legend2.C:5
double gauss2D(double *x, double *par)
Definition: TwoHistoFit2D.C:16
std::vector< double > values
Definition: TwoHistoFit2D.C:32
static void SetDefaultFitter(const char *name="")
static: set name of default fitter
TVirtualPad * cd(Int_t subpadnumber=0)
Set current canvas & pad.
Definition: TCanvas.cxx:659
void myFcn(Int_t &, Double_t *, Double_t &fval, Double_t *p, Int_t)
Definition: TwoHistoFit2D.C:35
virtual Double_t GetParameter(Int_t ipar) const =0
int Int_t
Definition: RtypesCore.h:41
virtual void SetFCN(void *fcn)
To set the address of the minimization objective function.
void SetStatY(Float_t y=0)
Definition: TStyle.h:392
virtual Int_t GetNbinsX() const
Definition: TH1.h:296
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
Double_t x[n]
Definition: legend1.C:17
static double p2(double t, double a, double b, double c)
TH2D * h2
Definition: fit2dHist.C:45
virtual Double_t Rndm(Int_t i=0)
Machine independent random number generator.
Definition: TRandom3.cxx:94
TH1F * h1
Definition: legend1.C:5
virtual Double_t GetParError(Int_t ipar) const =0
virtual TF1 * DrawCopy(Option_t *option="") const
Draw a copy of this function with its current attributes-*.
Definition: TF2.cxx:243
void FillHisto(TH2D *h, int n, double *p)
Definition: TwoHistoFit2D.C:49
virtual Double_t GetBinContent(Int_t bin) const
Return content of bin number bin.
Definition: TH2.h:90
virtual void SetChisquare(Double_t chi2)
Definition: TF1.h:412
ROOT::R::TRInterface & r
Definition: Object.C:4
Class to manage histogram axis.
Definition: TAxis.h:36
virtual void Draw(Option_t *option="")
Draw this histogram with options.
Definition: TH1.cxx:2878
void SetOptFit(Int_t fit=1)
The type of information about fit parameters printed in the histogram statistics box can be selected ...
Definition: TStyle.cxx:1204
virtual Int_t ExecuteCommand(const char *command, Double_t *args, Int_t nargs)=0
TAxis * GetYaxis()
Definition: TH1.h:320
static double p1(double t, double a, double b)
A 2-Dim function with parameters.
Definition: TF2.h:33
virtual Int_t GetStats(Double_t &amin, Double_t &edm, Double_t &errdef, Int_t &nvpar, Int_t &nparx) const =0
TRandom3 rndm
Definition: TwoHistoFit2D.C:48
The Canvas class.
Definition: TCanvas.h:48
virtual void SetRange(Double_t xmin, Double_t xmax)
Initialize the upper and lower bounds to draw the function.
Definition: TF2.h:154
double Double_t
Definition: RtypesCore.h:55
Double_t y[n]
Definition: legend1.C:17
double func(double *x, double *p)
Definition: stressTF1.cxx:213
Abstract Base Class for Fitting.
virtual Double_t GetParameter(Int_t ipar) const
Definition: TF1.h:352
virtual Double_t GetBinCenter(Int_t bin) const
Return center of bin.
Definition: TAxis.cxx:449
virtual Int_t GetNbinsY() const
Definition: TH1.h:297
static TVirtualFitter * Fitter(TObject *obj, Int_t maxpar=25)
Static function returning a pointer to the current fitter.
virtual void Divide(Int_t nx=1, Int_t ny=1, Float_t xmargin=0.01, Float_t ymargin=0.01, Int_t color=0)
Automatic pad generation by division.
Definition: TPad.cxx:1073
std::vector< std::pair< double, double > > coords
Definition: TwoHistoFit2D.C:31
virtual void Add(TObject *obj)
Definition: TList.h:81
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:3169
virtual void Sumw2(Bool_t flag=kTRUE)
Create structure to store sum of squares of weights.
Definition: TH1.cxx:8350
#define gPad
Definition: TVirtualPad.h:288
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:3102
int TwoHistoFit2D(bool global=true)
Definition: TwoHistoFit2D.C:85
std::vector< double > errors
Definition: TwoHistoFit2D.C:33
double exp(double)
virtual Double_t GetBinError(Int_t bin) const
Return value of error associated to bin number bin.
Definition: TH1.cxx:8395
Int_t Fill(Double_t)
Invalid Fill method.
Definition: TH2.cxx:287
virtual TFitResultPtr Fit(const char *formula, Option_t *option="", Option_t *goption="", Double_t xmin=0, Double_t xmax=0)
Fit histogram with function fname.
Definition: TH1.cxx:3607
const Int_t n
Definition: legend1.C:16
TAxis * GetXaxis()
Definition: TH1.h:319
double my2Dfunc(double *x, double *par)
Definition: TwoHistoFit2D.C:21
2-D histogram with a double per channel (see TH1 documentation)}
Definition: TH2.h:297