Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
testUnfold5a.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_unfold
3/// \notebook
4/// Test program for the classes TUnfoldDensity and TUnfoldBinning
5///
6/// A toy test of the TUnfold package
7///
8/// This is an example of unfolding a two-dimensional distribution
9/// also using an auxiliary measurement to constrain some background
10///
11/// The example comprises several macros
12/// - testUnfold5a.C create root files with TTree objects for
13/// signal, background and data
14/// - write files testUnfold5_signal.root
15/// testUnfold5_background.root
16/// testUnfold5_data.root
17///
18/// - testUnfold5b.C create a root file with the TUnfoldBinning objects
19/// - write file testUnfold5_binning.root
20///
21/// - testUnfold5c.C loop over trees and fill histograms based on the
22/// TUnfoldBinning objects
23/// - read testUnfold5_binning.root
24/// testUnfold5_signal.root
25/// testUnfold5_background.root
26/// testUnfold5_data.root
27///
28/// - write testUnfold5_histograms.root
29///
30/// - testUnfold5d.C run the unfolding
31/// - read testUnfold5_histograms.root
32/// - write testUnfold5_result.root
33/// testUnfold5_result.ps
34///
35/// \macro_output
36/// \macro_code
37///
38/// **Version 17.6, in parallel to changes in TUnfold**
39///
40/// #### History:
41/// - Version 17.5, in parallel to changes in TUnfold
42/// - Version 17.4, in parallel to changes in TUnfold
43/// - Version 17.3, in parallel to changes in TUnfold
44/// - Version 17.2, in parallel to changes in TUnfold
45/// - Version 17.1, in parallel to changes in TUnfold
46/// - Version 17.0 example for multi-dimensional unfolding
47///
48/// This file is part of TUnfold.
49///
50/// TUnfold is free software: you can redistribute it and/or modify
51/// it under the terms of the GNU General Public License as published by
52/// the Free Software Foundation, either version 3 of the License, or
53/// (at your option) any later version.
54///
55/// TUnfold is distributed in the hope that it will be useful,
56/// but WITHOUT ANY WARRANTY; without even the implied warranty of
57/// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
58/// GNU General Public License for more details.
59///
60/// You should have received a copy of the GNU General Public License
61/// along with TUnfold. If not, see <http://www.gnu.org/licenses/>.
62///
63/// \author Stefan Schmitt DESY, 14.10.2008
64
65#include <iostream>
66#include <cmath>
67#include <TMath.h>
68#include <TRandom3.h>
69#include <TFile.h>
70#include <TTree.h>
71
72using std::cout;
73
74TRandom *g_rnd=nullptr;
75
76class ToyEvent {
77public:
78 void GenerateDataEvent(TRandom *rnd);
79 void GenerateSignalEvent(TRandom *rnd);
80 void GenerateBgrEvent(TRandom *rnd);
81 // reconstructed quantities
82 inline Double_t GetPtRec(void) const { return fPtRec; }
83 inline Double_t GetEtaRec(void) const { return fEtaRec; }
84 inline Double_t GetDiscriminator(void) const {return fDiscriminator; }
85 inline Bool_t IsTriggered(void) const { return fIsTriggered; }
86
87 // generator level quantities
88 inline Double_t GetPtGen(void) const {
89 if(IsSignal()) return fPtGen;
90 else return -1.0;
91 }
92 inline Double_t GetEtaGen(void) const {
93 if(IsSignal()) return fEtaGen;
94 else return 999.0;
95 }
96 inline Bool_t IsSignal(void) const { return fIsSignal; }
97protected:
98
99 void GenerateSignalKinematics(TRandom *rnd,Bool_t isData);
100 void GenerateBgrKinematics(TRandom *rnd,Bool_t isData);
101 void GenerateReco(TRandom *rnd);
102
103 // reconstructed quantities
104 Double_t fPtRec;
105 Double_t fEtaRec;
106 Double_t fDiscriminator;
107 Bool_t fIsTriggered;
108 // generated quantities
109 Double_t fPtGen;
110 Double_t fEtaGen;
111 Bool_t fIsSignal;
112
113 static Double_t kDataSignalFraction;
114
115};
116
117void testUnfold5a()
118{
119 // random generator
120 g_rnd=new TRandom3();
121
122 // data and MC number of events
123 const Int_t neventData = 20000;
124 Double_t const neventSignalMC =2000000;
125 Double_t const neventBgrMC =2000000;
126
127 Float_t etaRec,ptRec,discr,etaGen,ptGen;
128 Int_t istriggered,issignal;
129
130 //==================================================================
131 // Step 1: generate data TTree
132
133 TFile *dataFile=new TFile("testUnfold5_data.root","recreate");
134 TTree *dataTree=new TTree("data","event");
135
136 dataTree->Branch("etarec",&etaRec,"etarec/F");
137 dataTree->Branch("ptrec",&ptRec,"ptrec/F");
138 dataTree->Branch("discr",&discr,"discr/F");
139
140 // for real data, only the triggered events are available
141 dataTree->Branch("istriggered",&istriggered,"istriggered/I");
142 // data truth parameters
143 dataTree->Branch("etagen",&etaGen,"etagen/F");
144 dataTree->Branch("ptgen",&ptGen,"ptgen/F");
145 dataTree->Branch("issignal",&issignal,"issignal/I");
146
147 cout<<"fill data tree\n";
148
149 Int_t nEvent=0,nTriggered=0;
150 while(nTriggered<neventData) {
151 ToyEvent event;
152 event.GenerateDataEvent(g_rnd);
153
154 etaRec=event.GetEtaRec();
155 ptRec=event.GetPtRec();
156 discr=event.GetDiscriminator();
157 istriggered=event.IsTriggered() ? 1 : 0;
158 etaGen=event.GetEtaGen();
159 ptGen=event.GetPtGen();
160 issignal=event.IsSignal() ? 1 : 0;
161
162 dataTree->Fill();
163
164 if(!(nEvent%100000)) cout<<" data event "<<nEvent<<"\n";
165
166 if(istriggered) nTriggered++;
167 nEvent++;
168
169 }
170
171 dataTree->Write();
172 delete dataTree;
173 delete dataFile;
174
175 //==================================================================
176 // Step 2: generate signal TTree
177
178 TFile *signalFile=new TFile("testUnfold5_signal.root","recreate");
179 TTree *signalTree=new TTree("signal","event");
180
181 signalTree->Branch("etarec",&etaRec,"etarec/F");
182 signalTree->Branch("ptrec",&ptRec,"ptrec/F");
183 signalTree->Branch("discr",&discr,"discr/F");
184 signalTree->Branch("istriggered",&istriggered,"istriggered/I");
185 signalTree->Branch("etagen",&etaGen,"etagen/F");
186 signalTree->Branch("ptgen",&ptGen,"ptgen/F");
187
188 cout<<"fill signal tree\n";
189
190 for(int ievent=0;ievent<neventSignalMC;ievent++) {
191 ToyEvent event;
192 event.GenerateSignalEvent(g_rnd);
193
194 etaRec=event.GetEtaRec();
195 ptRec=event.GetPtRec();
196 discr=event.GetDiscriminator();
197 istriggered=event.IsTriggered() ? 1 : 0;
198 etaGen=event.GetEtaGen();
199 ptGen=event.GetPtGen();
200
201 if(!(ievent%100000)) cout<<" signal event "<<ievent<<"\n";
202
203 signalTree->Fill();
204 }
205
206 signalTree->Write();
207 delete signalTree;
208 delete signalFile;
209
210 // ==============================================================
211 // Step 3: generate background MC TTree
212
213 TFile *bgrFile=new TFile("testUnfold5_background.root","recreate");
214 TTree *bgrTree=new TTree("background","event");
215
216 bgrTree->Branch("etarec",&etaRec,"etarec/F");
217 bgrTree->Branch("ptrec",&ptRec,"ptrec/F");
218 bgrTree->Branch("discr",&discr,"discr/F");
219 bgrTree->Branch("istriggered",&istriggered,"istriggered/I");
220
221 cout<<"fill background tree\n";
222
223 for(int ievent=0;ievent<neventBgrMC;ievent++) {
224 ToyEvent event;
225 event.GenerateBgrEvent(g_rnd);
226 etaRec=event.GetEtaRec();
227 ptRec=event.GetPtRec();
228 discr=event.GetDiscriminator();
229 istriggered=event.IsTriggered() ? 1 : 0;
230
231 if(!(ievent%100000)) cout<<" background event "<<ievent<<"\n";
232
233 bgrTree->Fill();
234 }
235
236 bgrTree->Write();
237 delete bgrTree;
238 delete bgrFile;
239}
240
241Double_t ToyEvent::kDataSignalFraction=0.8;
242
243void ToyEvent::GenerateDataEvent(TRandom *rnd) {
244 fIsSignal=rnd->Uniform()<kDataSignalFraction;
245 if(IsSignal()) {
246 GenerateSignalKinematics(rnd,kTRUE);
247 } else {
248 GenerateBgrKinematics(rnd,kTRUE);
249 }
250 GenerateReco(rnd);
251}
252
253void ToyEvent::GenerateSignalEvent(TRandom *rnd) {
254 fIsSignal=true;
255 GenerateSignalKinematics(rnd,kFALSE);
256 GenerateReco(rnd);
257}
258
259void ToyEvent::GenerateBgrEvent(TRandom *rnd) {
260 fIsSignal=false;
261 GenerateBgrKinematics(rnd,kFALSE);
262 GenerateReco(rnd);
263}
264
265void ToyEvent::GenerateSignalKinematics(TRandom *rnd,Bool_t isData) {
266 Double_t e_T0=0.5;
267 Double_t e_T0_eta=0.0;
268 Double_t e_n=2.0;
269 Double_t e_n_eta=0.0;
270 Double_t eta_p2=0.0;
271 Double_t etaMax=4.0;
272 if(isData) {
273 e_T0=0.6;
274 e_n=2.5;
275 e_T0_eta=0.05;
276 e_n_eta=-0.05;
277 eta_p2=1.5;
278 }
279 if(eta_p2>0.0) {
280 fEtaGen=TMath::Power(rnd->Uniform(),eta_p2)*etaMax;
281 if(rnd->Uniform()>=0.5) fEtaGen= -fEtaGen;
282 } else {
283 fEtaGen=rnd->Uniform(-etaMax,etaMax);
284 }
285 Double_t n=e_n + e_n_eta*fEtaGen;
286 Double_t T0=e_T0 + e_T0_eta*fEtaGen;
287 fPtGen=(TMath::Power(rnd->Uniform(),1./(1.-n))-1.)*T0;
288 /* static int print=100;
289 if(print) {
290 cout<<fEtaGen
291 <<" "<<fPtGen
292 <<"\n";
293 print--;
294 } */
295}
296
297void ToyEvent::GenerateBgrKinematics(TRandom *rnd,Bool_t isData) {
298 fPtGen=0.0;
299 fEtaGen=0.0;
300 fPtRec=rnd->Exp(isData ? 2.5 : 2.5);
301 fEtaRec=rnd->Uniform(-3.,3.);
302}
303
304void ToyEvent::GenerateReco(TRandom *rnd) {
305 if(fIsSignal) {
306 Double_t expEta=TMath::Exp(fEtaGen);
307 Double_t eGen=fPtGen*(expEta+1./expEta);
308 Double_t sigmaE=0.1*TMath::Sqrt(eGen)+(0.01+0.002*TMath::Abs(fEtaGen))
309 *eGen;
310 Double_t eRec;
311 do {
312 eRec=rnd->Gaus(eGen,sigmaE);
313 } while(eRec<=0.0);
314 Double_t sigmaEta=0.1+0.02*fEtaGen;
315 fEtaRec=rnd->Gaus(fEtaGen,sigmaEta);
316 fPtRec=eRec/(expEta+1./expEta);
317 do {
318 Double_t tauDiscr=0.08-0.04/(1.+fPtRec/10.0);
319 Double_t sigmaDiscr=0.01;
320 fDiscriminator=1.0-rnd->Exp(tauDiscr)+rnd->Gaus(0.,sigmaDiscr);
321 } while((fDiscriminator<=0.)||(fDiscriminator>=1.));
322 /* static int print=100;
323 if(print) {
324 cout<<fEtaGen<<" "<<fPtGen
325 <<" -> "<<fEtaRec<<" "<<fPtRec
326 <<"\n";
327 print--;
328 } */
329 } else {
330 do {
331 Double_t tauDiscr=0.15-0.05/(1.+fPtRec/5.0)+0.1*fEtaRec;
332 Double_t sigmaDiscr=0.02+0.01*fEtaRec;
333 fDiscriminator=rnd->Exp(tauDiscr)+rnd->Gaus(0.,sigmaDiscr);
334 } while((fDiscriminator<=0.)||(fDiscriminator>=1.));
335 }
336 fIsTriggered=(rnd->Uniform()<1./(TMath::Exp(-fPtRec+3.5)+1.));
337}
bool Bool_t
Definition RtypesCore.h:63
int Int_t
Definition RtypesCore.h:45
float Float_t
Definition RtypesCore.h:57
constexpr Bool_t kFALSE
Definition RtypesCore.h:101
double Double_t
Definition RtypesCore.h:59
constexpr Bool_t kTRUE
Definition RtypesCore.h:100
A ROOT file is an on-disk file, usually with extension .root, that stores objects in a file-system-li...
Definition TFile.h:53
Random number generator class based on M.
Definition TRandom3.h:27
This is the base class for the ROOT Random number generators.
Definition TRandom.h:27
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:275
virtual Double_t Exp(Double_t tau)
Returns an exponential deviate.
Definition TRandom.cxx:252
virtual Double_t Uniform(Double_t x1=1)
Returns a uniform deviate on the interval (0, x1).
Definition TRandom.cxx:682
A TTree represents a columnar dataset.
Definition TTree.h:79
virtual Int_t Fill()
Fill all branches.
Definition TTree.cxx:4603
TBranch * Branch(const char *name, T *obj, Int_t bufsize=32000, Int_t splitlevel=99)
Add a new branch, and infer the data type from the type of obj being passed.
Definition TTree.h:353
Int_t Write(const char *name=nullptr, Int_t option=0, Int_t bufsize=0) override
Write this object to the current directory.
Definition TTree.cxx:9753
const Int_t n
Definition legend1.C:16
T etaMax()
Function providing the maximum possible value of pseudorapidity for a non-zero rho,...
Definition etaMax.h:51
Double_t Exp(Double_t x)
Returns the base-e exponential function of x, which is e raised to the power x.
Definition TMath.h:709
Double_t Sqrt(Double_t x)
Returns the square root of x.
Definition TMath.h:662
LongDouble_t Power(LongDouble_t x, LongDouble_t y)
Returns x raised to the power y.
Definition TMath.h:721
Short_t Abs(Short_t d)
Returns the absolute value of parameter Short_t d.
Definition TMathBase.h:123