ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
TestBinomial.C
Go to the documentation of this file.
1 // Perform a fit to a set of data with binomial errors
2 // like those derived from the division of two histograms.
3 // Three different fits are performed and compared:
4 //
5 // - simple least square fit to the divided histogram obtained
6 // from TH1::Divide with option b
7 // - least square fit to the TGraphAsymmErrors obtained from
8 // TGraphAsymmErrors::BayesDivide
9 // - likelihood fit performed on the dividing histograms using
10 // binomial statistics with the TBinomialEfficiency class
11 //
12 // The first two methods are biased while the last one is statistical correct.
13 // Running the script passing an integer value n larger than 1, n fits are
14 // performed and the bias are also shown.
15 // To run the script :
16 //
17 // to show the bias performing 100 fits for 1000 events per "experiment"
18 // root[0]: .x TestBinomial.C+
19 //
20 // to show the bias performing 100 fits for 1000 events per "experiment"
21 // .x TestBinomial.C+(100, 1000)
22 //
23 //
25 #include "TVirtualFitter.h"
26 #include "TH1.h"
27 #include "TRandom3.h"
28 #include "TF1.h"
29 #include "TFitResult.h"
30 #include "TStyle.h"
31 #include "TCanvas.h"
32 #include "TLegend.h"
33 #include "TPaveStats.h"
34 #include "Math/IntegratorOptions.h"
35 #include <cassert>
36 #include <iostream>
37 
38 void TestBinomial(int nloop = 100, int nevts = 100, bool plot=false, bool debug = false, int seed = 111)
39 {
41  gStyle->SetLineWidth(2.0);
42  gStyle->SetOptStat(11);
43 
44  TObjArray hbiasNorm;
45  hbiasNorm.Add(new TH1D("h0Norm", "Bias Histogram fit",100,-5,5));
46  hbiasNorm.Add(new TH1D("h1Norm","Bias Binomial fit",100,-5,5));
47  TObjArray hbiasThreshold;
48  hbiasThreshold.Add(new TH1D("h0Threshold", "Bias Histogram fit",100,-5,5));
49  hbiasThreshold.Add(new TH1D("h1Threshold","Bias Binomial fit",100,-5,5));
50  TObjArray hbiasWidth;
51  hbiasWidth.Add(new TH1D("h0Width", "Bias Histogram fit",100,-5,5));
52  hbiasWidth.Add(new TH1D("h1Width","Bias Binomial fit",100,-5,5));
53  TH1D* hChisquared = new TH1D("hChisquared",
54  "#chi^{2} probability (Baker-Cousins)", 200, 0.0, 1.0);
55 
58 
59  // Note: in order to be able to use TH1::FillRandom() to generate
60  // pseudo-experiments, we use a trick: generate "selected"
61  // and "non-selected" samples independently. These are
62  // statistically independent and therefore can be safely
63  // added to yield the "before selection" sample.
64 
65 
66  // Define (arbitrarily?) a distribution of input events.
67  // Here: assume a x^(-2) distribution. Boundaries: [10, 100].
68 
69  Double_t xmin =10, xmax = 100;
70  TH1D* hM2D = new TH1D("hM2D", "x^(-2) denominator distribution",
71  45, xmin, xmax);
72  TH1D* hM2N = new TH1D("hM2N", "x^(-2) numerator distribution",
73  45, xmin, xmax);
74  TH1D* hM2E = new TH1D("hM2E", "x^(-2) efficiency",
75  45, xmin, xmax);
76 
77  TF1* fM2D = new TF1("fM2D", "(1-[0]/(1+exp(([1]-x)/[2])))/(x*x)",
78  xmin, xmax);
79  TF1* fM2N = new TF1("fM2N", "[0]/(1+exp(([1]-x)/[2]))/(x*x)",
80  xmin, xmax);
81  TF1* fM2Fit = new TF1("fM2Fit", "[0]/(1+exp(([1]-x)/[2]))",
82  xmin, xmax);
83  TF1* fM2Fit2 = 0;
84 
85  TRandom3 rb(seed);
86 
87  // First try: use a single set of parameters.
88  // For each try, we need to find the overall normalization
89 
90  Double_t norm = 0.80;
91  Double_t threshold = 25.0;
92  Double_t width = 5.0;
93 
94  fM2D->SetParameter(0, norm);
95  fM2D->SetParameter(1, threshold);
96  fM2D->SetParameter(2, width);
97  fM2N->SetParameter(0, norm);
98  fM2N->SetParameter(1, threshold);
99  fM2N->SetParameter(2, width);
100  Double_t integralN = fM2N->Integral(xmin, xmax);
101  Double_t integralD = fM2D->Integral(xmin, xmax);
102  Double_t fracN = integralN/(integralN+integralD);
103  Int_t nevtsN = rb.Binomial(nevts, fracN);
104  Int_t nevtsD = nevts - nevtsN;
105 
106  std::cout << nevtsN << " " << nevtsD << std::endl;
107 
108  gStyle->SetOptFit(1111);
109 
110  // generate many times to see the bias
111  for (int iloop = 0; iloop < nloop; ++iloop) {
112 
113  // generate pseudo-experiments
114  hM2D->Reset();
115  hM2N->Reset();
116  hM2D->FillRandom(fM2D->GetName(), nevtsD);
117  hM2N->FillRandom(fM2N->GetName(), nevtsN);
118  hM2D->Add(hM2N);
119 
120  // construct the "efficiency" histogram
121  hM2N->Sumw2();
122  hM2E->Divide(hM2N, hM2D, 1, 1, "b");
123 
124  // Fit twice, using the same fit function.
125  // In the first (standard) fit, initialize to (arbitrary) values.
126  // In the second fit, use the results from the first fit (this
127  // makes it easier for the fit -- but the purpose here is not to
128  // show how easy or difficult it is to obtain results, but whether
129  // the CORRECT results are obtained or not!).
130 
131  fM2Fit->SetParameter(0, 0.5);
132  fM2Fit->SetParameter(1, 15.0);
133  fM2Fit->SetParameter(2, 2.0);
134  fM2Fit->SetParError(0, 0.1);
135  fM2Fit->SetParError(1, 1.0);
136  fM2Fit->SetParError(2, 0.2);
137  TH1 * hf = fM2Fit->GetHistogram();
138  // std::cout << "Function values " << std::endl;
139  // for (int i = 1; i <= hf->GetNbinsX(); ++i)
140  // std::cout << hf->GetBinContent(i) << " ";
141  // std::cout << std::endl;
142 
143  TCanvas* cEvt;
144  if (plot) {
145  cEvt = new TCanvas(Form("cEnv%d",iloop),
146  Form("plots for experiment %d", iloop),
147  1000, 600);
148  cEvt->Divide(1,2);
149  cEvt->cd(1);
150  hM2D->DrawCopy("HIST");
151  hM2N->SetLineColor(kRed);
152  hM2N->DrawCopy("HIST SAME");
153  cEvt->cd(2);
154  }
155  for (int fit = 0; fit < 2; ++fit) {
156  Int_t status = 0;
157  switch (fit) {
158  case 0:
159  {
160  // TVirtualPad * pad = gPad;
161  // new TCanvas();
162  // fM2Fit->Draw();
163  // gPad = pad;
164  TString optFit = "RN";
165  if (debug) optFit += TString("SV");
166  TFitResultPtr res = hM2E->Fit(fM2Fit, optFit);
167  if (plot) {
168  hM2E->DrawCopy("E");
169  fM2Fit->DrawCopy("SAME");
170  }
171  if (debug) res->Print();
172  status = res;
173  break;
174  }
175  case 1:
176  {
177  // if (fM2Fit2) delete fM2Fit2;
178  // fM2Fit2 = dynamic_cast<TF1*>(fM2Fit->Clone("fM2Fit2"));
179  fM2Fit2 = fM2Fit; // do not clone/copy the function
180  if (fM2Fit2->GetParameter(0) >= 1.0)
181  fM2Fit2->SetParameter(0, 0.95);
182  fM2Fit2->SetParLimits(0, 0.0, 1.0);
183 
184  // TVirtualPad * pad = gPad;
185  // new TCanvas();
186  // fM2Fit2->Draw();
187  // gPad = pad;
188 
189  TBinomialEfficiencyFitter bef(hM2N, hM2D);
190  TString optFit = "RI";
191  if (debug) optFit += TString("SV");
192  TFitResultPtr res = bef.Fit(fM2Fit2,optFit);
193  status = res;
194  if (status !=0) {
195  std::cerr << "Error performing binomial efficiency fit, result = "
196  << status << std::endl;
197  res->Print();
198  continue;
199  }
200  if (plot) {
201  fM2Fit2->SetLineColor(kRed);
202  fM2Fit2->DrawCopy("SAME");
203  }
204  if (debug) {
205  res->Print();
206  }
207  }
208  }
209 
210  if (status != 0) break;
211 
212  Double_t fnorm = fM2Fit->GetParameter(0);
213  Double_t enorm = fM2Fit->GetParError(0);
214  Double_t fthreshold = fM2Fit->GetParameter(1);
215  Double_t ethreshold = fM2Fit->GetParError(1);
216  Double_t fwidth = fM2Fit->GetParameter(2);
217  Double_t ewidth = fM2Fit->GetParError(2);
218  if (fit == 1) {
219  fnorm = fM2Fit2->GetParameter(0);
220  enorm = fM2Fit2->GetParError(0);
221  fthreshold = fM2Fit2->GetParameter(1);
222  ethreshold = fM2Fit2->GetParError(1);
223  fwidth = fM2Fit2->GetParameter(2);
224  ewidth = fM2Fit2->GetParError(2);
225  hChisquared->Fill(fM2Fit2->GetProb());
226  }
227 
228  TH1D* h = dynamic_cast<TH1D*>(hbiasNorm[fit]);
229  h->Fill((fnorm-norm)/enorm);
230  h = dynamic_cast<TH1D*>(hbiasThreshold[fit]);
231  h->Fill((fthreshold-threshold)/ethreshold);
232  h = dynamic_cast<TH1D*>(hbiasWidth[fit]);
233  h->Fill((fwidth-width)/ewidth);
234  }
235  }
236 
237 
238  TCanvas* c1 = new TCanvas("c1",
239  "Efficiency fit biases",10,10,1000,800);
240  c1->Divide(2,2);
241 
242  TH1D *h0, *h1;
243  c1->cd(1);
244  h0 = dynamic_cast<TH1D*>(hbiasNorm[0]);
245  h0->Draw("HIST");
246  h1 = dynamic_cast<TH1D*>(hbiasNorm[1]);
247  h1->SetLineColor(kRed);
248  h1->Draw("HIST SAMES");
249  TLegend* l1 = new TLegend(0.1, 0.75, 0.5, 0.9,
250  "plateau parameter", "ndc");
251  l1->AddEntry(h0, Form("histogram: mean = %4.2f RMS = \
252  %4.2f", h0->GetMean(), h0->GetRMS()), "l");
253  l1->AddEntry(h1, Form("binomial : mean = %4.2f RMS = \
254  %4.2f", h1->GetMean(), h1->GetRMS()), "l");
255  l1->Draw();
256 
257  c1->cd(2);
258  h0 = dynamic_cast<TH1D*>(hbiasThreshold[0]);
259  h0->Draw("HIST");
260  h1 = dynamic_cast<TH1D*>(hbiasThreshold[1]);
261  h1->SetLineColor(kRed);
262  h1->Draw("HIST SAMES");
263  TLegend* l2 = new TLegend(0.1, 0.75, 0.5, 0.9,
264  "threshold parameter", "ndc");
265  l2->AddEntry(h0, Form("histogram: mean = %4.2f RMS = \
266  %4.2f", h0->GetMean(), h0->GetRMS()), "l");
267  l2->AddEntry(h1, Form("binomial : mean = %4.2f RMS = \
268  %4.2f", h1->GetMean(), h1->GetRMS()), "l");
269  l2->Draw();
270 
271  c1->cd(3);
272  h0 = dynamic_cast<TH1D*>(hbiasWidth[0]);
273  h0->Draw("HIST");
274  h1 = dynamic_cast<TH1D*>(hbiasWidth[1]);
275  h1->SetLineColor(kRed);
276  h1->Draw("HIST SAMES");
277  TLegend* l3 = new TLegend(0.1, 0.75, 0.5, 0.9, "width parameter", "ndc");
278  l3->AddEntry(h0, Form("histogram: mean = %4.2f RMS = \
279  %4.2f", h0->GetMean(), h0->GetRMS()), "l");
280  l3->AddEntry(h1, Form("binomial : mean = %4.2f RMS = \
281  %4.2f", h1->GetMean(), h1->GetRMS()), "l");
282  l3->Draw();
283 
284  c1->cd(4);
285  hChisquared->Draw("HIST");
286 }
287 
288 int main() {
289  TestBinomial();
290 }
virtual void SetLineWidth(Width_t lwidth)
Definition: TAttLine.h:57
virtual Int_t Fill(Double_t x)
Increment bin with abscissa X by 1.
Definition: TH1.cxx:3159
An array of TObjects.
Definition: TObjArray.h:39
float xmin
Definition: THbookFile.cxx:93
Random number generator class based on M.
Definition: TRandom3.h:29
This class displays a legend box (TPaveText) containing several legend entries.
Definition: TLegend.h:35
virtual Int_t Binomial(Int_t ntot, Double_t prob)
Generates a random integer N according to the binomial law.
Definition: TRandom.cxx:172
TCanvas * c1
Definition: legend1.C:2
Definition: Rtypes.h:61
Binomial fitter for the division of two histograms.
R__EXTERN TStyle * gStyle
Definition: TStyle.h:423
virtual void Draw(Option_t *option="")
Draw this legend with its current attributes.
Definition: TLegend.cxx:373
TH1 * h
Definition: legend2.C:5
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
static void SetDefaultIntegrator(const char *name)
Basic string class.
Definition: TString.h:137
TAlienJobStatus * status
Definition: TAlienJob.cxx:51
int Int_t
Definition: RtypesCore.h:41
virtual Double_t GetParError(Int_t ipar) const
Return value of parameter number ipar.
Definition: TF1.cxx:1608
virtual Double_t Integral(Double_t a, Double_t b, Double_t epsrel=1.e-12)
IntegralOneDim or analytical integral.
Definition: TF1.cxx:2241
virtual TH1 * GetHistogram() const
Return a pointer to the histogram used to vusualize the function.
Definition: TF1.cxx:1274
void plot(TString fname="data.root", TString var0="var0", TString var1="var1")
Definition: createData.C:17
TLine l1(2.5, 4.5, 15.5, 4.5)
virtual TH1 * DrawCopy(Option_t *option="", const char *name_postfix="_copy") const
Copy this histogram and Draw in the current pad.
Definition: TH1.cxx:2925
TH1F * h1
Definition: legend1.C:5
virtual void Reset(Option_t *option="")
Reset.
Definition: TH1.cxx:9698
virtual Bool_t Divide(TF1 *f1, Double_t c1=1)
Performs the operation: this = this/(c1*f1) if errors are defined (see TH1::Sumw2), errors are also recalculated.
Definition: TH1.cxx:2630
virtual void SetLineColor(Color_t lcolor)
Definition: TAttLine.h:54
virtual void SetParLimits(Int_t ipar, Double_t parmin, Double_t parmax)
Set limits for parameter ipar.
Definition: TF1.cxx:3183
virtual Double_t GetMean(Int_t axis=1) const
For axis = 1,2 or 3 returns the mean value of the histogram along X,Y or Z axis.
Definition: TH1.cxx:7014
virtual void FillRandom(const char *fname, Int_t ntimes=5000)
Fill histogram following distribution in function fname.
Definition: TH1.cxx:3330
virtual void Draw(Option_t *option="")
Draw this histogram with options.
Definition: TH1.cxx:2878
virtual void SetParError(Int_t ipar, Double_t error)
Set error for parameter number ipar.
Definition: TF1.cxx:3158
Provides an indirection to the TFitResult class and with a semantics identical to a TFitResult pointe...
Definition: TFitResultPtr.h:33
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
char * Form(const char *fmt,...)
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:51
virtual void SetMarkerStyle(Style_t mstyle=1)
Definition: TAttMarker.h:53
float xmax
Definition: THbookFile.cxx:93
Double_t GetRMS(Int_t axis=1) const
Definition: TH1.h:315
virtual Double_t GetProb() const
Return the fit probability.
Definition: TF1.cxx:1633
1-D histogram with a double per channel (see TH1 documentation)}
Definition: TH1.h:613
The Canvas class.
Definition: TCanvas.h:48
double Double_t
Definition: RtypesCore.h:55
TLegendEntry * AddEntry(const TObject *obj, const char *label="", Option_t *option="lpf")
Add a new entry to this legend.
Definition: TLegend.cxx:280
The TH1 histogram class.
Definition: TH1.h:80
virtual Bool_t Add(TF1 *h1, Double_t c1=1, Option_t *option="")
Performs the operation: this = this + c1*f1 if errors are defined (see TH1::Sumw2), errors are also recalculated.
Definition: TH1.cxx:780
virtual Double_t GetParameter(Int_t ipar) const
Definition: TF1.h:352
virtual void Print(Option_t *option="") const
Print result of the fit, by default chi2, parameter values and errors.
Definition: TFitResult.cxx:42
TFitResultPtr Fit(TF1 *f1, Option_t *option="")
Carry out the fit of the given function to the given histograms.
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
1-Dim function class
Definition: TF1.h:149
virtual void Sumw2(Bool_t flag=kTRUE)
Create structure to store sum of squares of weights.
Definition: TH1.cxx:8350
void SetOptStat(Int_t stat=1)
The type of information printed in the histogram statistics box can be selected via the parameter mod...
Definition: TStyle.cxx:1252
bool debug
void TestBinomial(int nloop=100, int nevts=100, bool plot=false, bool debug=false, int seed=111)
Definition: TestBinomial.C:38
void Add(TObject *obj)
Definition: TObjArray.h:75
virtual void SetParameter(Int_t param, Double_t value)
Definition: TF1.h:424
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
double norm(double *x, double *p)
Definition: unuranDistr.cxx:40
int main()
Definition: TestBinomial.C:288
virtual TF1 * DrawCopy(Option_t *option="") const
Draw a copy of this function with its current attributes.
Definition: TF1.cxx:1080