Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
chi2test.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_math
3 /// \notebook
4 /// Example to use chi2 test for comparing two histograms
5 /// One unweighted histogram is compared with a weighted histogram.
6 /// The normalized residuals are retrieved and plotted in a simple graph.
7 /// The QQ plot of the normalized residual using the
8 /// normal distribution is also plotted.
9 ///
10 /// \macro_image
11 /// \macro_output
12 /// \macro_code
13 ///
14 /// \author Nikolai Gagunashvili, Daniel Haertl, Lorenzo Moneta
15 
16 #include "TH1.h"
17 #include "TH1D.h"
18 #include "TF1.h"
19 #include "TGraph.h"
20 #include "TGraphQQ.h"
21 #include "TCanvas.h"
22 #include "TStyle.h"
23 #include "TMath.h"
24 
25 TCanvas * chi2test(Float_t w=0)
26 {
27  // Note: The parameter w is used to produce the 2 pictures in
28  // the TH1::Chi2Test method. The 1st picture is produced with
29  // w=0 and the 2nd with w=17 (see TH1::Chi2Test() help).
30 
31  // Define Histograms.
32  const Int_t n = 20;
33 
34  TH1D *h1 = new TH1D("h1", "h1", n, 4, 16);
35  TH1D *h2 = new TH1D("h2", "h2", n, 4, 16);
36 
37  h1->SetTitle("Unweighted Histogram");
38  h2->SetTitle("Weighted Histogram");
39 
40  h1->SetBinContent(1, 0);
41  h1->SetBinContent(2, 1);
42  h1->SetBinContent(3, 0);
43  h1->SetBinContent(4, 1);
44  h1->SetBinContent(5, 1);
45  h1->SetBinContent(6, 6);
46  h1->SetBinContent(7, 7);
47  h1->SetBinContent(8, 2);
48  h1->SetBinContent(9, 22);
49  h1->SetBinContent(10, 30);
50  h1->SetBinContent(11, 27);
51  h1->SetBinContent(12, 20);
52  h1->SetBinContent(13, 13);
53  h1->SetBinContent(14, 9);
54  h1->SetBinContent(15, 9 + w);
55  h1->SetBinContent(16, 13);
56  h1->SetBinContent(17, 19);
57  h1->SetBinContent(18, 11);
58  h1->SetBinContent(19, 9);
59  h1->SetBinContent(20, 0);
60 
61  h2->SetBinContent(1, 2.20173025 );
62  h2->SetBinContent(2, 3.30143857);
63  h2->SetBinContent(3, 2.5892849);
64  h2->SetBinContent(4, 2.99990201);
65  h2->SetBinContent(5, 4.92877054);
66  h2->SetBinContent(6, 8.33036995);
67  h2->SetBinContent(7, 6.95084763);
68  h2->SetBinContent(8, 15.206357);
69  h2->SetBinContent(9, 23.9236012);
70  h2->SetBinContent(10, 44.3848114);
71  h2->SetBinContent(11, 49.4465599);
72  h2->SetBinContent(12, 25.1868858);
73  h2->SetBinContent(13, 16.3129692);
74  h2->SetBinContent(14, 13.0289612);
75  h2->SetBinContent(15, 16.7857609);
76  h2->SetBinContent(16, 22.9914703);
77  h2->SetBinContent(17, 30.5279255);
78  h2->SetBinContent(18, 12.5252123);
79  h2->SetBinContent(19, 16.4104557);
80  h2->SetBinContent(20, 7.86067867);
81  h2->SetBinError(1, 0.38974303 );
82  h2->SetBinError(2, 0.536510944);
83  h2->SetBinError(3, 0.529702604);
84  h2->SetBinError(4, 0.642001867);
85  h2->SetBinError(5, 0.969341516);
86  h2->SetBinError(6, 1.47611344);
87  h2->SetBinError(7, 1.69797957);
88  h2->SetBinError(8, 3.28577447);
89  h2->SetBinError(9, 5.40784931);
90  h2->SetBinError(10, 9.10106468);
91  h2->SetBinError(11, 9.73541737);
92  h2->SetBinError(12, 5.55019951);
93  h2->SetBinError(13, 3.57914758);
94  h2->SetBinError(14, 2.77877331);
95  h2->SetBinError(15, 3.23697519);
96  h2->SetBinError(16, 4.3608489);
97  h2->SetBinError(17, 5.77172089);
98  h2->SetBinError(18, 3.38666105);
99  h2->SetBinError(19, 2.98861837);
100  h2->SetBinError(20, 1.58402085);
101 
102  h1->SetEntries(217);
103  h2->SetEntries(500);
104 
105  //apply the chi2 test and retrieve the residuals
106  Double_t res[n], x[20];
107  h1->Chi2Test(h2,"UW P",res);
108 
109  //Graph for Residuals
110  for (Int_t i=0; i<n; i++) x[i]= 4.+i*12./20.+12./40.;
111  TGraph *resgr = new TGraph(n,x,res);
112  resgr->GetXaxis()->SetRangeUser(4,16);
113  resgr->GetYaxis()->SetRangeUser(-3.5,3.5);
114  resgr->GetYaxis()->SetTitle("Normalized Residuals");
115  resgr->SetMarkerStyle(21);
116  resgr->SetMarkerColor(2);
117  resgr->SetMarkerSize(.9);
118  resgr->SetTitle("Normalized Residuals");
119 
120  //Quantile-Quantile plot
121  TF1 *f = new TF1("f","TMath::Gaus(x,0,1)",-10,10);
122  TGraphQQ *qqplot = new TGraphQQ(n,res,f);
123  qqplot->SetMarkerStyle(20);
124  qqplot->SetMarkerColor(2);
125  qqplot->SetMarkerSize(.9);
126  qqplot->SetTitle("Q-Q plot of Normalized Residuals");
127 
128  //create Canvas
129  TCanvas *c1 = new TCanvas("c1","Chistat Plot",10,10,700,600);
130  c1->Divide(2,2);
131 
132  // Draw Histogramms and Graphs
133  c1->cd(1);
134  h1->SetMarkerColor(4);
135  h1->SetMarkerStyle(20);
136 
137  h1->Draw("E");
138 
139  c1->cd(2);
140  h2->Draw("");
141  h2->SetMarkerColor(4);
142  h2->SetMarkerStyle(20);
143 
144  c1->cd(3);
145  gPad->SetGridy();
146  resgr->Draw("APL");
147 
148  c1->cd(4);
149  qqplot->Draw("AP");
150 
151  c1->cd(0);
152 
153  c1->Update();
154  return c1;
155 }