ROOT logo
// Example to use chi2 test for comparing two histograms 
// One unweighted histogram is compared with a weighted histogram. 
// The normalized residuals are retrieved and plotted in a simple graph. 
// The QQ plot of the normalized residual using the 
// normal distribution is also plotted. 
//
//Authors: Nikolai Gagunashvili, Daniel Haertl, Lorenzo Moneta


#include "TH1.h"
#include "TH1D.h"
#include "TF1.h"
#include "TGraph.h"
#include "TGraphQQ.h"
#include "TCanvas.h"
#include "TStyle.h"
#include "TMath.h"


TCanvas * chi2test(Float_t w=0)
{
  // Note: The parameter w is used to produce the 2 pictures in
  // the TH1::Chi2Test method. The 1st picture is produced with 
  // w=0 and the 2nd with w=17 (see TH1::Chi2Test() help).

  // Define Histograms.
  const Int_t n = 20;

  TH1D *h1 = new TH1D("h1", "h1", n, 4, 16);
  TH1D *h2 = new TH1D("h2", "h2", n, 4, 16);

  h1->SetTitle("Unweighted Histogram");
  h2->SetTitle("Weighted Histogram");

  h1->SetBinContent(1, 0);
  h1->SetBinContent(2, 1);
  h1->SetBinContent(3, 0);
  h1->SetBinContent(4, 1);
  h1->SetBinContent(5, 1);
  h1->SetBinContent(6, 6);
  h1->SetBinContent(7, 7);
  h1->SetBinContent(8, 2);
  h1->SetBinContent(9, 22);
  h1->SetBinContent(10, 30);
  h1->SetBinContent(11, 27);
  h1->SetBinContent(12, 20);
  h1->SetBinContent(13, 13);
  h1->SetBinContent(14, 9);
  h1->SetBinContent(15, 9 + w);
  h1->SetBinContent(16, 13);
  h1->SetBinContent(17, 19);
  h1->SetBinContent(18, 11);
  h1->SetBinContent(19, 9);
  h1->SetBinContent(20, 0);

  // error is automatically the sqrt of the bin content
//   h1->SetBinError(1, 0);
//   h1->SetBinError(2, 1);
//   h1->SetBinError(3, 0);
//   h1->SetBinError(4, 1);
//   h1->SetBinError(5, 1);
//   h1->SetBinError(6, 2.44948983);
//   h1->SetBinError(7, 2.64575124 );
//   h1->SetBinError(8, 1.41421354);
//   h1->SetBinError(9, 4.69041586);
//   h1->SetBinError(10, 5.47722578);
//   h1->SetBinError(11, 5.19615221);
//   h1->SetBinError(12, 4.47213602);
//   h1->SetBinError(13, 3.60555124);
//   h1->SetBinError(14, 3.);
//   h1->SetBinError(15, 5.38516474);
//   h1->SetBinError(16, 3.60555124);
//   h1->SetBinError(17, 4.35889912);
//   h1->SetBinError(18, 3.31662488);
//   h1->SetBinError(19, 3.);
//   h1->SetBinError(20, 0);

  h2->SetBinContent(1, 2.20173025 );
  h2->SetBinContent(2, 3.30143857);
  h2->SetBinContent(3, 2.5892849);
  h2->SetBinContent(4, 2.99990201);
  h2->SetBinContent(5, 4.92877054);
  h2->SetBinContent(6, 8.33036995);
  h2->SetBinContent(7, 6.95084763);
  h2->SetBinContent(8, 15.206357);
  h2->SetBinContent(9, 23.9236012);
  h2->SetBinContent(10, 44.3848114);
  h2->SetBinContent(11, 49.4465599);
  h2->SetBinContent(12, 25.1868858);
  h2->SetBinContent(13, 16.3129692);
  h2->SetBinContent(14, 13.0289612);
  h2->SetBinContent(15, 16.7857609);
  h2->SetBinContent(16, 22.9914703);
  h2->SetBinContent(17, 30.5279255);
  h2->SetBinContent(18, 12.5252123);
  h2->SetBinContent(19, 16.4104557);
  h2->SetBinContent(20, 7.86067867);
  h2->SetBinError(1, 0.38974303 );
  h2->SetBinError(2, 0.536510944);
  h2->SetBinError(3, 0.529702604);
  h2->SetBinError(4, 0.642001867);
  h2->SetBinError(5, 0.969341516);
  h2->SetBinError(6, 1.47611344);
  h2->SetBinError(7, 1.69797957);
  h2->SetBinError(8, 3.28577447);
  h2->SetBinError(9, 5.40784931);
  h2->SetBinError(10, 9.10106468);
  h2->SetBinError(11, 9.73541737);
  h2->SetBinError(12, 5.55019951);
  h2->SetBinError(13, 3.57914758);
  h2->SetBinError(14, 2.77877331);
  h2->SetBinError(15, 3.23697519);
  h2->SetBinError(16, 4.3608489);
  h2->SetBinError(17, 5.77172089);
  h2->SetBinError(18, 3.38666105);
  h2->SetBinError(19, 2.98861837);
  h2->SetBinError(20, 1.58402085);

  h1->SetEntries(217);
  h2->SetEntries(500);

//apply the chi2 test and retrieve the residuals
  Double_t res[n], x[20];
  h1->Chi2Test(h2,"UW P",res);

//Graph for Residuals
  for (Int_t i=0; i<n; i++) x[i]= 4.+i*12./20.+12./40.;
  TGraph *resgr = new TGraph(n,x,res);
  resgr->GetXaxis()->SetRangeUser(4,16);
  resgr->GetYaxis()->SetRangeUser(-3.5,3.5);
  resgr->GetYaxis()->SetTitle("Normalized Residuals");
  resgr->SetMarkerStyle(21);
  resgr->SetMarkerColor(2);
  resgr->SetMarkerSize(.9);
  resgr->SetTitle("Normalized Residuals");

//Quantile-Quantile plot
  TF1 *f = new TF1("f","TMath::Gaus(x,0,1)",-10,10);
  TGraphQQ *qqplot = new TGraphQQ(n,res,f);
  qqplot->SetMarkerStyle(20);
  qqplot->SetMarkerColor(2);
  qqplot->SetMarkerSize(.9);
  qqplot->SetTitle("Q-Q plot of Normalized Residuals");

//create Canvas
  TCanvas *c1 = new TCanvas("c1","Chistat Plot",10,10,700,600);
  c1->SetFillColor(33);
  c1->Divide(2,2);

// Draw Histogramms and Graphs
  c1->cd(1);
  gPad->SetFrameFillColor(21);
  h1->SetMarkerColor(4);
  h1->SetMarkerStyle(20);

  h1->Draw("E");

  c1->cd(2);
  gPad->SetFrameFillColor(21);
  h2->Draw("");
  h2->SetMarkerColor(4);
  h2->SetMarkerStyle(20);

  c1->cd(3);
  gPad->SetFrameFillColor(21);
  gPad->SetGridy();
  resgr->Draw("APL");

  c1->cd(4);
  gPad->SetFrameFillColor(21);
  qqplot->Draw("AP");

  c1->cd(0);
  
  c1->Update();
   return c1;
}
 chi2test.C:1
 chi2test.C:2
 chi2test.C:3
 chi2test.C:4
 chi2test.C:5
 chi2test.C:6
 chi2test.C:7
 chi2test.C:8
 chi2test.C:9
 chi2test.C:10
 chi2test.C:11
 chi2test.C:12
 chi2test.C:13
 chi2test.C:14
 chi2test.C:15
 chi2test.C:16
 chi2test.C:17
 chi2test.C:18
 chi2test.C:19
 chi2test.C:20
 chi2test.C:21
 chi2test.C:22
 chi2test.C:23
 chi2test.C:24
 chi2test.C:25
 chi2test.C:26
 chi2test.C:27
 chi2test.C:28
 chi2test.C:29
 chi2test.C:30
 chi2test.C:31
 chi2test.C:32
 chi2test.C:33
 chi2test.C:34
 chi2test.C:35
 chi2test.C:36
 chi2test.C:37
 chi2test.C:38
 chi2test.C:39
 chi2test.C:40
 chi2test.C:41
 chi2test.C:42
 chi2test.C:43
 chi2test.C:44
 chi2test.C:45
 chi2test.C:46
 chi2test.C:47
 chi2test.C:48
 chi2test.C:49
 chi2test.C:50
 chi2test.C:51
 chi2test.C:52
 chi2test.C:53
 chi2test.C:54
 chi2test.C:55
 chi2test.C:56
 chi2test.C:57
 chi2test.C:58
 chi2test.C:59
 chi2test.C:60
 chi2test.C:61
 chi2test.C:62
 chi2test.C:63
 chi2test.C:64
 chi2test.C:65
 chi2test.C:66
 chi2test.C:67
 chi2test.C:68
 chi2test.C:69
 chi2test.C:70
 chi2test.C:71
 chi2test.C:72
 chi2test.C:73
 chi2test.C:74
 chi2test.C:75
 chi2test.C:76
 chi2test.C:77
 chi2test.C:78
 chi2test.C:79
 chi2test.C:80
 chi2test.C:81
 chi2test.C:82
 chi2test.C:83
 chi2test.C:84
 chi2test.C:85
 chi2test.C:86
 chi2test.C:87
 chi2test.C:88
 chi2test.C:89
 chi2test.C:90
 chi2test.C:91
 chi2test.C:92
 chi2test.C:93
 chi2test.C:94
 chi2test.C:95
 chi2test.C:96
 chi2test.C:97
 chi2test.C:98
 chi2test.C:99
 chi2test.C:100
 chi2test.C:101
 chi2test.C:102
 chi2test.C:103
 chi2test.C:104
 chi2test.C:105
 chi2test.C:106
 chi2test.C:107
 chi2test.C:108
 chi2test.C:109
 chi2test.C:110
 chi2test.C:111
 chi2test.C:112
 chi2test.C:113
 chi2test.C:114
 chi2test.C:115
 chi2test.C:116
 chi2test.C:117
 chi2test.C:118
 chi2test.C:119
 chi2test.C:120
 chi2test.C:121
 chi2test.C:122
 chi2test.C:123
 chi2test.C:124
 chi2test.C:125
 chi2test.C:126
 chi2test.C:127
 chi2test.C:128
 chi2test.C:129
 chi2test.C:130
 chi2test.C:131
 chi2test.C:132
 chi2test.C:133
 chi2test.C:134
 chi2test.C:135
 chi2test.C:136
 chi2test.C:137
 chi2test.C:138
 chi2test.C:139
 chi2test.C:140
 chi2test.C:141
 chi2test.C:142
 chi2test.C:143
 chi2test.C:144
 chi2test.C:145
 chi2test.C:146
 chi2test.C:147
 chi2test.C:148
 chi2test.C:149
 chi2test.C:150
 chi2test.C:151
 chi2test.C:152
 chi2test.C:153
 chi2test.C:154
 chi2test.C:155
 chi2test.C:156
 chi2test.C:157
 chi2test.C:158
 chi2test.C:159
 chi2test.C:160
 chi2test.C:161
 chi2test.C:162
 chi2test.C:163
 chi2test.C:164
 chi2test.C:165
 chi2test.C:166
 chi2test.C:167
 chi2test.C:168
 chi2test.C:169
 chi2test.C:170
 chi2test.C:171
 chi2test.C:172
 chi2test.C:173
 chi2test.C:174
 chi2test.C:175
 chi2test.C:176
 chi2test.C:177
 chi2test.C:178
thumb