Logo ROOT  
Reference Guide
 
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Modules Pages
Loading...
Searching...
No Matches
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///
6/// One unweighted histogram is compared with a weighted histogram.
7/// The normalized residuals are retrieved and plotted in a simple graph.
8/// The QQ plot of the normalized residual using the
9/// normal distribution is also plotted.
10///
11/// \macro_image
12/// \macro_output
13/// \macro_code
14///
15/// \author Nikolai Gagunashvili, Daniel Haertl, Lorenzo Moneta
16
17#include "TH1.h"
18#include "TH1D.h"
19#include "TF1.h"
20#include "TGraph.h"
21#include "TGraphQQ.h"
22#include "TCanvas.h"
23#include "TStyle.h"
24#include "TMath.h"
25
27{
28 // Note: The parameter w is used to produce the 2 pictures in
29 // the TH1::Chi2Test method. The 1st picture is produced with
30 // w=0 and the 2nd with w=17 (see TH1::Chi2Test() help).
31
32 // Define Histograms.
33 const Int_t n = 20;
34
35 TH1D *h1 = new TH1D("h1", "h1", n, 4, 16);
36 n, 4, 16);
37
38 h1->SetTitle("Unweighted Histogram");
39 h2->SetTitle("Weighted Histogram");
40
41 h1->SetBinContent(1, 0);
42 h1->SetBinContent(2, 1);
43 h1->SetBinContent(3, 0);
44 h1->SetBinContent(4, 1);
45 h1->SetBinContent(5, 1);
46 h1->SetBinContent(6, 6);
47 h1->SetBinContent(7, 7);
48 h1->SetBinContent(8, 2);
49 h1->SetBinContent(9, 22);
50 h1->SetBinContent(10, 30);
51 h1->SetBinContent(11, 27);
52 h1->SetBinContent(12, 20);
53 h1->SetBinContent(13, 13);
54 h1->SetBinContent(14, 9);
55 h1->SetBinContent(15, 9 + w);
56 h1->SetBinContent(16, 13);
57 h1->SetBinContent(17, 19);
58 h1->SetBinContent(18, 11);
59 h1->SetBinContent(19, 9);
60 h1->SetBinContent(20, 0);
61
62 h2->SetBinContent(1, 2.20173025 );
63 h2->SetBinContent(2, 3.30143857);
64 h2->SetBinContent(3, 2.5892849);
65 h2->SetBinContent(4, 2.99990201);
66 h2->SetBinContent(5, 4.92877054);
67 h2->SetBinContent(6, 8.33036995);
68 h2->SetBinContent(7, 6.95084763);
69 h2->SetBinContent(8, 15.206357);
70 h2->SetBinContent(9, 23.9236012);
71 h2->SetBinContent(10, 44.3848114);
72 h2->SetBinContent(11, 49.4465599);
73 h2->SetBinContent(12, 25.1868858);
74 h2->SetBinContent(13, 16.3129692);
75 h2->SetBinContent(14, 13.0289612);
76 h2->SetBinContent(15, 16.7857609);
77 h2->SetBinContent(16, 22.9914703);
78 h2->SetBinContent(17, 30.5279255);
79 h2->SetBinContent(18, 12.5252123);
80 h2->SetBinContent(19, 16.4104557);
81 h2->SetBinContent(20, 7.86067867);
82 h2->SetBinError(1, 0.38974303 );
83 h2->SetBinError(2, 0.536510944);
84 h2->SetBinError(3, 0.529702604);
85 h2->SetBinError(4, 0.642001867);
86 h2->SetBinError(5, 0.969341516);
87 h2->SetBinError(6, 1.47611344);
88 h2->SetBinError(7, 1.69797957);
89 h2->SetBinError(8, 3.28577447);
90 h2->SetBinError(9, 5.40784931);
91 h2->SetBinError(10, 9.10106468);
92 h2->SetBinError(11, 9.73541737);
93 h2->SetBinError(12, 5.55019951);
94 h2->SetBinError(13, 3.57914758);
95 h2->SetBinError(14, 2.77877331);
96 h2->SetBinError(15, 3.23697519);
97 h2->SetBinError(16, 4.3608489);
98 h2->SetBinError(17, 5.77172089);
99 h2->SetBinError(18, 3.38666105);
100 h2->SetBinError(19, 2.98861837);
101 h2->SetBinError(20, 1.58402085);
102
103 h1->SetEntries(217);
104 h2->SetEntries(500);
105
106 //apply the chi2 test and retrieve the residuals
107 Double_t res[n], x[20];
108 h1->Chi2Test(h2,"UW P",res);
109
110 //Graph for Residuals
111 for (Int_t i=0; i<n; i++) x[i]= 4.+i*12./20.+12./40.;
112 TGraph *resgr = new TGraph(n,x,res);
113 resgr->GetXaxis()->SetRangeUser(4,16);
114 resgr->GetYaxis()->SetRangeUser(-3.5,3.5);
115 resgr->GetYaxis()->SetTitle("Normalized Residuals");
116 resgr->SetMarkerStyle(21);
117 resgr->SetMarkerColor(2);
118 resgr->SetMarkerSize(.9);
119 resgr->SetTitle("Normalized Residuals");
120
121 //Quantile-Quantile plot
122 TF1 *f = new TF1("f","TMath::Gaus(x,0,1)",-10,10);
123 TGraphQQ *qqplot = new TGraphQQ(n,res,f);
124 qqplot->SetMarkerStyle(20);
125 qqplot->SetMarkerColor(2);
126 qqplot->SetMarkerSize(.9);
127 qqplot->SetTitle("Q-Q plot of Normalized Residuals");
128
129 //create Canvas
130 TCanvas *c1 = new TCanvas("c1","Chistat Plot",10,10,700,600);
131 c1->Divide(2,2);
132
133 // Draw Histogramms and Graphs
134 c1->cd(1);
135 h1->SetMarkerColor(4);
136 h1->SetMarkerStyle(20);
137
138 h1->Draw("E");
139
140 c1->cd(2);
141 h2->Draw("");
142 h2->SetMarkerColor(4);
143 h2->SetMarkerStyle(20);
144
145 c1->cd(3);
146 gPad->SetGridy();
147 resgr->Draw("APL");
148
149 c1->cd(4);
150 qqplot->Draw("AP");
151
152 c1->cd(0);
153
154 c1->Update();
155 return c1;
156}
#define f(i)
Definition RSha256.hxx:104
int Int_t
Definition RtypesCore.h:45
float Float_t
Definition RtypesCore.h:57
double Double_t
Definition RtypesCore.h:59
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
#define gPad
virtual void SetMarkerColor(Color_t mcolor=1)
Set the marker color.
Definition TAttMarker.h:39
virtual void SetMarkerStyle(Style_t mstyle=1)
Set the marker style.
Definition TAttMarker.h:41
The Canvas class.
Definition TCanvas.h:23
1-Dim function class
Definition TF1.h:234
This class allows to draw quantile-quantile plots.
Definition TGraphQQ.h:18
A TGraph is an object made of two arrays X and Y with npoints each.
Definition TGraph.h:41
1-D histogram with a double per channel (see TH1 documentation)
Definition TH1.h:698
void SetTitle(const char *title) override
Change/set the title.
Definition TH1.cxx:6721
void Draw(Option_t *option="") override
Draw this histogram with options.
Definition TH1.cxx:3038
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:9218
virtual Double_t Chi2Test(const TH1 *h2, Option_t *option="UU", Double_t *res=nullptr) const
test for comparing weighted and unweighted histograms.
Definition TH1.cxx:1980
virtual void SetEntries(Double_t n)
Definition TH1.h:411
return c1
Definition legend1.C:41
Double_t x[n]
Definition legend1.C:17
const Int_t n
Definition legend1.C:16
TH1F * h1
Definition legend1.C:5