Logo ROOT   6.07/09
Reference Guide
testUnfold4.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_unfold
3 /// \notebook
4 /// Test program for the class TUnfoldSys
5 ///
6 /// Simple toy tests of the TUnfold package
7 ///
8 /// Pseudo data (5000 events) are unfolded into three components
9 /// The unfolding is performed once without and once with area constraint
10 ///
11 /// Ideally, the pulls may show that the result is biased if no constraint
12 /// is applied. This is expected because the true data errors are not known,
13 /// and instead the sqrt(data) errors are used.
14 ///
15 /// History:
16 /// - Version 16.1, parallel to changes in TUnfold
17 /// - Version 16.0, parallel to changes in TUnfold
18 /// - Version 15, use L-curve scan to scan the average correlation
19 ///
20 /// \macro_image
21 /// \macro_output
22 /// \macro_code
23 ///
24 /// \author Stefan Schmitt, DESY
25 
26 #include <TMath.h>
27 #include <TCanvas.h>
28 #include <TRandom3.h>
29 #include <TFitter.h>
30 #include <TF1.h>
31 #include <TStyle.h>
32 #include <TVector.h>
33 #include <TGraph.h>
34 
35 #include "TUnfoldDensity.h"
36 
37 using namespace std;
38 
39 TRandom *rnd=0;
40 
41 Int_t GenerateGenEvent(Int_t nmax,const Double_t *probability) {
42  // choose an integer random number in the range [0,nmax]
43  // (the generator level bin)
44  // depending on the probabilities
45  // probability[0],probability[1],...probability[nmax-1]
46  Double_t f=rnd->Rndm();
47  Int_t r=0;
48  while((r<nmax)&&(f>=probability[r])) {
49  f -= probability[r];
50  r++;
51  }
52  return r;
53 }
54 
55 Double_t GenerateRecEvent(const Double_t *shapeParm) {
56  // return a coordinate (the reconstructed variable)
57  // depending on shapeParm[]
58  // shapeParm[0]: fraction of events with Gaussian distribution
59  // shapeParm[1]: mean of Gaussian
60  // shapeParm[2]: width of Gaussian
61  // (1-shapeParm[0]): fraction of events with flat distribution
62  // shapeParm[3]: minimum of flat component
63  // shapeParm[4]: maximum of flat component
64  Double_t f=rnd->Rndm();
65  Double_t r;
66  if(f<shapeParm[0]) {
67  r=rnd->Gaus(shapeParm[1],shapeParm[2]);
68  } else {
69  r=rnd->Rndm()*(shapeParm[4]-shapeParm[3])+shapeParm[3];
70  }
71  return r;
72 }
73 
74 void testUnfold4()
75 {
76  // switch on histogram errors
78 
79  // random generator
80  rnd=new TRandom3();
81 
82  // data and MC number of events
83  Double_t const nData0= 500.0;
84  Double_t const nMC0 = 50000.0;
85 
86  // Binning
87  // reconstructed variable (0-10)
88  Int_t const nDet=15;
89  Double_t const xminDet=0.0;
90  Double_t const xmaxDet=15.0;
91 
92  // signal binning (three shapes: 0,1,2)
93  Int_t const nGen=3;
94  Double_t const xminGen=-0.5;
95  Double_t const xmaxGen= 2.5;
96 
97  // parameters
98  // fraction of events per signal shape
99  static const Double_t genFrac[]={0.3,0.6,0.1};
100 
101  // signal shapes
102  static const Double_t genShape[][5]=
103  {{1.0,2.0,1.5,0.,15.},
104  {1.0,7.0,2.5,0.,15.},
105  {0.0,0.0,0.0,0.,15.}};
106 
107  // define DATA histograms
108  // observed data distribution
109  TH1D *histDetDATA=new TH1D("Yrec",";DATA(Yrec)",nDet,xminDet,xmaxDet);
110 
111  // define MC histograms
112  // matrix of migrations
113  TH2D *histGenDetMC=new TH2D("Yrec%Xgen","MC(Xgen,Yrec)",
114  nGen,xminGen,xmaxGen,nDet,xminDet,xmaxDet);
115 
116  TH1D *histUnfold=new TH1D("Xgen",";DATA(Xgen)",nGen,xminGen,xmaxGen);
117 
118  TH1D **histPullNC=new TH1D* [nGen];
119  TH1D **histPullArea=new TH1D* [nGen];
120  for(int i=0;i<nGen;i++) {
121  histPullNC[i]=new TH1D(TString::Format("PullNC%d",i),"pull",15,-3.,3.);
122  histPullArea[i]=new TH1D(TString::Format("PullArea%d",i),"pull",15,-3.,3.);
123  }
124 
125  // this method is new in version 16 of TUnfold
126  cout<<"TUnfold version is "<<TUnfold::GetTUnfoldVersion()<<"\n";
127 
128  for(int itoy=0;itoy<1000;itoy++) {
129  if(!(itoy %10)) cout<<"toy iteration: "<<itoy<<"\n";
130  histDetDATA->Reset();
131  histGenDetMC->Reset();
132 
133  Int_t nData=rnd->Poisson(nData0);
134  for(Int_t i=0;i<nData;i++) {
135  Int_t iGen=GenerateGenEvent(nGen,genFrac);
136  Double_t yObs=GenerateRecEvent(genShape[iGen]);
137  histDetDATA->Fill(yObs);
138  }
139 
140  Int_t nMC=rnd->Poisson(nMC0);
141  for(Int_t i=0;i<nMC;i++) {
142  Int_t iGen=GenerateGenEvent(nGen,genFrac);
143  Double_t yObs=GenerateRecEvent(genShape[iGen]);
144  histGenDetMC->Fill(iGen,yObs);
145  }
146  /* for(Int_t ix=0;ix<=histGenDetMC->GetNbinsX()+1;ix++) {
147  for(Int_t iy=0;iy<=histGenDetMC->GetNbinsY()+1;iy++) {
148  cout<<ix<<iy<<" : "<<histGenDetMC->GetBinContent(ix,iy)<<"\n";
149  }
150  } */
151  //========================
152  // unfolding
153 
154  // switch off info messages
155 #ifndef DEBUG
156  gErrorIgnoreLevel = 1001;
157 #endif
158 
159  TUnfoldSys unfold(histGenDetMC,TUnfold::kHistMapOutputHoriz,
161  // define the input vector (the measured data distribution)
162  unfold.SetInput(histDetDATA,0.0,1.0);
163 
164  // run the unfolding
165  unfold.ScanLcurve(50,0.,0.,0,0,0);
166 
167  // fill pull distributions without constraint
168  unfold.GetOutput(histUnfold);
169 
170  for(int i=0;i<nGen;i++) {
171  histPullNC[i]->Fill((histUnfold->GetBinContent(i+1)-genFrac[i]*nData0)/
172  histUnfold->GetBinError(i+1));
173  }
174 
175  // repeat unfolding on the same data, now with Area constraint
176  unfold.SetConstraint(TUnfold::kEConstraintArea);
177 
178  // run the unfolding
179  unfold.ScanLcurve(50,0.,0.,0,0,0);
180 
181  // fill pull distributions with constraint
182  unfold.GetOutput(histUnfold);
183 
184  for(int i=0;i<nGen;i++) {
185  histPullArea[i]->Fill((histUnfold->GetBinContent(i+1)-genFrac[i]*nData0)/
186  histUnfold->GetBinError(i+1));
187  }
188 
189  }
190  TCanvas *output = new TCanvas();
191  output->Divide(3,2);
192 
193  gStyle->SetOptFit(1111);
194 
195  for(int i=0;i<nGen;i++) {
196  output->cd(i+1);
197  histPullNC[i]->Fit("gaus");
198  histPullNC[i]->Draw();
199  }
200  for(int i=0;i<nGen;i++) {
201  output->cd(i+4);
202  histPullArea[i]->Fit("gaus");
203  histPullArea[i]->Draw();
204  }
205 }
virtual Int_t Fill(Double_t x)
Increment bin with abscissa X by 1.
Definition: TH1.cxx:3127
Random number generator class based on M.
Definition: TRandom3.h:29
virtual Double_t GetBinContent(Int_t bin) const
Return content of bin number bin.
Definition: TH1.cxx:4640
R__EXTERN Int_t gErrorIgnoreLevel
Definition: TError.h:107
virtual Double_t Gaus(Double_t mean=0, Double_t sigma=1)
Samples a random number from the standard Normal (Gaussian) Distribution with the given mean and sigm...
Definition: TRandom.cxx:235
R__EXTERN TStyle * gStyle
Definition: TStyle.h:418
TVirtualPad * cd(Int_t subpadnumber=0)
Set current canvas & pad.
Definition: TCanvas.cxx:659
TUnfold is used to decompose a measurement y into several sources x given the measurement uncertainti...
Definition: TUnfoldSys.h:51
int Int_t
Definition: RtypesCore.h:41
STL namespace.
static void SetDefaultSumw2(Bool_t sumw2=kTRUE)
When this static function is called with sumw2=kTRUE, all new histograms will automatically activate ...
Definition: TH1.cxx:5969
static TString Format(const char *fmt,...)
Static method which formats a string using a printf style format descriptor and return a TString...
Definition: TString.cxx:2335
This is the base class for the ROOT Random number generators.
Definition: TRandom.h:31
virtual void Reset(Option_t *option="")
Reset.
Definition: TH1.cxx:9343
virtual Double_t Rndm()
Machine independent random number generator.
Definition: TRandom.cxx:512
TRandom2 r(17)
virtual void Draw(Option_t *option="")
Draw this histogram with options.
Definition: TH1.cxx:2853
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:1209
static const char * GetTUnfoldVersion(void)
Definition: TUnfold.cxx:289
tomato 1-D histogram with a double per channel (see TH1 documentation)}
Definition: TH1.h:618
The Canvas class.
Definition: TCanvas.h:41
double f(double x)
double Double_t
Definition: RtypesCore.h:55
THist< 2, double, THistStatContent, THistStatUncertainty > TH2D
Definition: THist.hxx:307
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:1089
THist< 1, double, THistStatContent, THistStatUncertainty > TH1D
Definition: THist.hxx:301
virtual void Reset(Option_t *option="")
Reset this histogram: contents, errors, etc.
Definition: TH2.cxx:3811
virtual Int_t Poisson(Double_t mean)
Generates a random integer N according to a Poisson law.
Definition: TRandom.cxx:362
static void output(int code)
Definition: gifencode.c:226
virtual Double_t GetBinError(Int_t bin) const
Return value of error associated to bin number bin.
Definition: TH1.cxx:8130
Int_t Fill(Double_t)
Invalid Fill method.
Definition: TH2.cxx:292
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:3565
tomato 2-D histogram with a double per channel (see TH1 documentation)}
Definition: TH2.h:296