78 void GenerateDataEvent(
TRandom *rnd);
79 void GenerateSignalEvent(
TRandom *rnd);
80 void GenerateBgrEvent(
TRandom *rnd);
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; }
88 inline Double_t GetPtGen(
void)
const {
89 if(IsSignal())
return fPtGen;
92 inline Double_t GetEtaGen(
void)
const {
93 if(IsSignal())
return fEtaGen;
96 inline Bool_t IsSignal(
void)
const {
return fIsSignal; }
101 void GenerateReco(
TRandom *rnd);
113 static Double_t kDataSignalFraction;
123 const Int_t neventData = 20000;
124 Double_t const neventSignalMC =2000000;
125 Double_t const neventBgrMC =2000000;
127 Float_t etaRec,ptRec,discr,etaGen,ptGen;
128 Int_t istriggered,issignal;
133 TFile *dataFile=
new TFile(
"testUnfold5_data.root",
"recreate");
136 dataTree->
Branch(
"etarec",&etaRec,
"etarec/F");
137 dataTree->
Branch(
"ptrec",&ptRec,
"ptrec/F");
138 dataTree->
Branch(
"discr",&discr,
"discr/F");
141 dataTree->
Branch(
"istriggered",&istriggered,
"istriggered/I");
143 dataTree->
Branch(
"etagen",&etaGen,
"etagen/F");
144 dataTree->
Branch(
"ptgen",&ptGen,
"ptgen/F");
145 dataTree->
Branch(
"issignal",&issignal,
"issignal/I");
147 cout<<
"fill data tree\n";
149 Int_t nEvent=0,nTriggered=0;
150 while(nTriggered<neventData) {
152 event.GenerateDataEvent(g_rnd);
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;
164 if(!(nEvent%100000)) cout<<
" data event "<<nEvent<<
"\n";
166 if(istriggered) nTriggered++;
178 TFile *signalFile=
new TFile(
"testUnfold5_signal.root",
"recreate");
179 TTree *signalTree=
new TTree(
"signal",
"event");
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");
188 cout<<
"fill signal tree\n";
190 for(
int ievent=0;ievent<neventSignalMC;ievent++) {
192 event.GenerateSignalEvent(g_rnd);
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();
201 if(!(ievent%100000)) cout<<
" signal event "<<ievent<<
"\n";
213 TFile *bgrFile=
new TFile(
"testUnfold5_background.root",
"recreate");
214 TTree *bgrTree=
new TTree(
"background",
"event");
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");
221 cout<<
"fill background tree\n";
223 for(
int ievent=0;ievent<neventBgrMC;ievent++) {
225 event.GenerateBgrEvent(g_rnd);
226 etaRec=
event.GetEtaRec();
227 ptRec=
event.GetPtRec();
228 discr=
event.GetDiscriminator();
229 istriggered=
event.IsTriggered() ? 1 : 0;
231 if(!(ievent%100000)) cout<<
" background event "<<ievent<<
"\n";
241Double_t ToyEvent::kDataSignalFraction=0.8;
243void ToyEvent::GenerateDataEvent(
TRandom *rnd) {
244 fIsSignal=rnd->
Uniform()<kDataSignalFraction;
246 GenerateSignalKinematics(rnd,
kTRUE);
248 GenerateBgrKinematics(rnd,
kTRUE);
253void ToyEvent::GenerateSignalEvent(
TRandom *rnd) {
255 GenerateSignalKinematics(rnd,
kFALSE);
259void ToyEvent::GenerateBgrEvent(
TRandom *rnd) {
261 GenerateBgrKinematics(rnd,
kFALSE);
265void ToyEvent::GenerateSignalKinematics(
TRandom *rnd,
Bool_t isData) {
281 if(rnd->
Uniform()>=0.5) fEtaGen= -fEtaGen;
283 fEtaGen=rnd->
Uniform(-etaMax,etaMax);
286 Double_t T0=e_T0 + e_T0_eta*fEtaGen;
297void ToyEvent::GenerateBgrKinematics(
TRandom *rnd,
Bool_t isData) {
300 fPtRec=rnd->
Exp(isData ? 2.5 : 2.5);
304void ToyEvent::GenerateReco(
TRandom *rnd) {
307 Double_t eGen=fPtGen*(expEta+1./expEta);
312 eRec=rnd->
Gaus(eGen,sigmaE);
315 fEtaRec=rnd->
Gaus(fEtaGen,sigmaEta);
316 fPtRec=eRec/(expEta+1./expEta);
318 Double_t tauDiscr=0.08-0.04/(1.+fPtRec/10.0);
320 fDiscriminator=1.0-rnd->
Exp(tauDiscr)+rnd->
Gaus(0.,sigmaDiscr);
321 }
while((fDiscriminator<=0.)||(fDiscriminator>=1.));
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.));
A ROOT file is an on-disk file, usually with extension .root, that stores objects in a file-system-li...
Random number generator class based on M.
This is the base class for the ROOT Random number generators.
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...
virtual Double_t Exp(Double_t tau)
Returns an exponential deviate.
virtual Double_t Uniform(Double_t x1=1)
Returns a uniform deviate on the interval (0, x1).
A TTree represents a columnar dataset.
virtual Int_t Fill()
Fill all branches.
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.
Int_t Write(const char *name=nullptr, Int_t option=0, Int_t bufsize=0) override
Write this object to the current directory.
T etaMax()
Function providing the maximum possible value of pseudorapidity for a non-zero rho,...
Double_t Exp(Double_t x)
Returns the base-e exponential function of x, which is e raised to the power x.
Double_t Sqrt(Double_t x)
Returns the square root of x.
LongDouble_t Power(LongDouble_t x, LongDouble_t y)
Returns x raised to the power y.
Short_t Abs(Short_t d)
Returns the absolute value of parameter Short_t d.