Fork me on GitHub

source: svn/trunk/Examples/src/Analysis_Ex.cc@ 286

Last change on this file since 286 was 270, checked in by severine ovyn, 16 years ago

remove double GEN pfff

File size: 16.0 KB
RevLine 
[260]1/***********************************************************************
2** **
3** /----------------------------------------------\ **
4** | Delphes, a framework for the fast simulation | **
5** | of a generic collider experiment | **
6** \----------------------------------------------/ **
7** **
8** **
9** This package uses: **
10** ------------------ **
11** FastJet algorithm: Phys. Lett. B641 (2006) [hep-ph/0512210] **
12** Hector: JINST 2:P09005 (2007) [physics.acc-ph:0707.1198v2] **
13** FROG: [hep-ex/0901.2718v1] **
14** **
15** ------------------------------------------------------------------ **
16** **
17** Main authors: **
18** ------------- **
19** **
20** Severine Ovyn Xavier Rouby **
21** severine.ovyn@uclouvain.be xavier.rouby@cern **
22** **
23** Center for Particle Physics and Phenomenology (CP3) **
24** Universite catholique de Louvain (UCL) **
25** Louvain-la-Neuve, Belgium **
26** **
27** Copyright (C) 2008-2009, **
28** All rights reserved. **
29** **
30***********************************************************************/
31
32
[81]33#include "Examples/interface/Analysis_Ex.h"
[227]34#include <iostream>
35#include <sstream>
36#include <fstream>
37#include <iomanip>
[81]38
39using namespace std;
40
41//******************************Debut de l'analyse****************************************
42//*****************************************************************************************
43
44Analysis_Ex::Analysis_Ex(string CardWithCuts,string LogName)
45{
46 string temp_string;
47 istringstream curstring;
48
49 ifstream fichier_a_lire(CardWithCuts.c_str());
50 if(!fichier_a_lire.good()) {
51 cout << "DataCardname " << CardWithCuts << " not found, use default values" << endl;
52 return;
53 }
54
55 while (getline(fichier_a_lire,temp_string)) {
56 curstring.clear(); // needed when using several times istringstream::str(string)
57 curstring.str(temp_string);
58 string varname;
59 float value;
60
61 if(strstr(temp_string.c_str(),"#")) { }//remove comments
62 else if(strstr(temp_string.c_str(),"PT_ELEC")){curstring >> varname >> value; PT_ELEC = value;}
63 else if(strstr(temp_string.c_str(),"PT_MUON")){curstring >> varname >> value; PT_MUON = value;}
64 else if(strstr(temp_string.c_str(),"INV_MASS_LL")){curstring >> varname >> value; INV_MASS_LL = value;}
65 }
66
67 ofstream f_out(LogName.c_str(),ofstream::app);
68
69 f_out<<"*******************************************************************"<<endl;
70 f_out << left << setw(30) <<"Cut values used in the analysis: "<<""
71 << right << setw(37) <<"-------------------------------------"<<"\n";
72 f_out << left <<setw(50) << "Invariant mass of the leptons: "<<""
73 << right <<setw(17) << INV_MASS_LL <<"\n";
74 f_out<<"*******************************************************************"<<endl;
75
76}
77
78Analysis_Ex::~Analysis_Ex()
79{
80}
81
[84]82void Analysis_Ex::Run(ExRootTreeReader *treeReaderGen, ExRootTreeReader *treeReaderRec, ExRootTreeReader *treeReaderTrig, ExRootTreeWriter *treeWriter)
[81]83{
84 total=0;//initialisation of total number of events
85 cut_trig=0;
86 cut_1=0;//initialisation of counter for cut 1
87 cut_2=0;
88 //access the branches************************
89 //to get the generator level information
90 const TClonesArray *GEN = treeReaderGen->UseBranch("Particle");
91
92 //to get the reconstructed level information
93 const TClonesArray *JET = treeReaderRec->UseBranch("Jet");
94 const TClonesArray *TAUJET = treeReaderRec->UseBranch("TauJet");
95 const TClonesArray *PHOTO = treeReaderRec->UseBranch("Photon");
96 const TClonesArray *ELEC = treeReaderRec->UseBranch("Electron");
97 const TClonesArray *MUON = treeReaderRec->UseBranch("Muon");
98 const TClonesArray *TRACKS = treeReaderRec->UseBranch("Tracks");
99 const TClonesArray *CALO = treeReaderRec->UseBranch("CaloTower");
100
101 //to get the VFD reconstructed level information
102 const TClonesArray *ZDC = treeReaderRec->UseBranch("ZDChits");
103 const TClonesArray *RP220 = treeReaderRec->UseBranch("RP220hits");
104 const TClonesArray *FP420 = treeReaderRec->UseBranch("FP420hits");
105
106 //to get the trigger information
107 const TClonesArray *TRIGGER = treeReaderTrig->UseBranch("TrigResult");
[84]108
109 //Define the branches that will be filled during the analysis
110 ExRootTreeBranch *INVMASS = treeWriter->NewBranch("INVMass", TRootInvm::Class());
111 TRootInvm *inv_mass;
[81]112 //*******************************************
113
114 //run on the events
115 Long64_t entry, allEntries = treeReaderRec->GetEntries();
116 cout << "** Chain contains " << allEntries << " events" << endl;
117 total=allEntries;
118
119 //general information
120 float E,Px,Py,Pz;
121 float PT,Eta,Phi;
122
123 //lepton information
124 bool IsolFlag;
125
126 //bjet information
127 bool Btag;
128
129 //Particle level information
130 int PID, Status, M1,M2,D1,D2;
131 float Charge, T, X, Y, Z, M;
132
133 //VFD information
134 float S,q2,Tx,Ty;
135 int side;
136
137 for(entry = 0; entry < allEntries; ++entry)
138 {
139 treeReaderGen->ReadEntry(entry);//access information of generated information
140 treeReaderRec->ReadEntry(entry);//access information of reconstructed information
141 treeReaderTrig->ReadEntry(entry);//access information of Trigger information
142
143 //*****************************************************
144 //Example how to run on the generator level information
145 //*****************************************************
146 TIter itGen((TCollection*)GEN);
[270]147 GenParticle *gen;
[81]148 itGen.Reset();
[270]149 while( (gen = (GenParticle*) itGen.Next()) )
[81]150 {
151 PID = gen->PID; // particle HEP ID number
152 Status = gen->Status; // particle status
153 M1 = gen->M1; // particle 1st mother
154 M2 = gen->M2; // particle 2nd mother
155 D1 = gen->D1; // particle 1st daughter
156 D2 = gen->D2; // particle 2nd daughter
[270]157 Charge = gen->Charge; // electrical charge
[81]158
159 T = gen->T; // particle vertex position (t component)
160 X = gen->X; // particle vertex position (x component)
161 Y = gen->Y; // particle vertex position (y component)
162 Z = gen->Z; // particle vertex position (z component)
163 M = gen->M; // particle mass
164 }
165
166
167 //***********************************************
168 //Example how to run on the reconstructed objects
169 //***********************************************
170
171 //access the Electron branch
172 TIter itElec((TCollection*)ELEC);
173 TRootElectron *elec;
174 itElec.Reset();
175 while( (elec = (TRootElectron*) itElec.Next()) )
176 {
177 E = elec->E; // particle energy in GeV
178 Px = elec->Px; // particle momentum vector (x component) in GeV
179 Py = elec->Py; // particle momentum vector (y component) in GeV
180 Pz = elec->Pz; // particle momentum vector (z component) in GeV
181
182 PT = elec->PT; // particle transverse momentum in GeV
183 Eta = elec->Eta; // particle pseudorapidity
184 Phi = elec->Phi; // particle azimuthal angle in rad
185 IsolFlag = elec->IsolFlag; // is the particule isolated?
186 }
187 //Running on the muon branch is identical:
188 TIter itMuon((TCollection*)MUON);
189 TRootMuon *muon;
190 itMuon.Reset();
191 while( (muon = (TRootMuon*) itMuon.Next()) ){}
192
193 //access the Photon branch
194 TIter itGam((TCollection*)PHOTO);
195 TRootPhoton *gam;
196 itGam.Reset();
197 while( (gam = (TRootPhoton*) itGam.Next()) )
198 {
199 E = gam->E; // particle energy in GeV
200 Px = gam->Px; // particle momentum vector (x component) in GeV
201 Py = gam->Py; // particle momentum vector (y component) in GeV
202 Pz = gam->Pz; // particle momentum vector (z component) in GeV
203
204 PT = gam->PT; // particle transverse momentum in GeV
205 Eta = gam->Eta; // particle pseudorapidity
206 Phi = gam->Phi; // particle azimuthal angle in rad
207 }
208
209 //access the jet branch
210 TIter itJet((TCollection*)JET);
211 TRootJet *jet;
212 itJet.Reset();
213 while( (jet = (TRootJet*) itJet.Next()) )
214 {
215 E = jet->E; // particle energy in GeV
216 Px = jet->Px; // particle momentum vector (x component) in GeV
217 Py = jet->Py; // particle momentum vector (y component) in GeV
218 Pz = jet->Pz; // particle momentum vector (z component) in GeV
219
220 PT = jet->PT; // particle transverse momentum in GeV
221 Eta = jet->Eta; // particle pseudorapidity
222 Phi = jet->Phi; // particle azimuthal angle in rad
223 Btag = jet->Btag; // is the jet BTagged
224 }
225 //Running on the tau-jet branch is identical:
226 TIter itTaujet((TCollection*)TAUJET);
227 TRootTauJet *taujet;
228 itTaujet.Reset();
229 while( (taujet = (TRootTauJet*) itTaujet.Next()) ){}
230
231 //access the track branch
232 TIter itTrack((TCollection*)TRACKS);
233 TRootTracks *tracks;
234 itTrack.Reset();
235 while( (tracks = (TRootTracks*) itTrack.Next()) )
236 {
237 E = tracks->E; // particle energy in GeV
238 Px = tracks->Px; // particle momentum vector (x component) in GeV
239 Py = tracks->Py; // particle momentum vector (y component) in GeV
240 Pz = tracks->Pz; // particle momentum vector (z component) in GeV
241
242 PT = tracks->PT; // particle transverse momentum in GeV
243 Eta = tracks->Eta; // particle pseudorapidity
244 Phi = tracks->Phi; // particle azimuthal angle in rad
245 }
246
247 //Running on the calo branch is identical:
248 TIter itCalo((TCollection*)CALO);
249 TRootCalo *calo;
250 itCalo.Reset();
251 while( (calo = (TRootCalo*) itCalo.Next()) ){}
252
253 //***************************************************
254 //Example how to run on the VFD reconstructed objects
255 //***************************************************
256
257 //access the ZDC branch
258 TIter itZdc((TCollection*)ZDC);
259 TRootZdcHits *zdc;
260 itZdc.Reset();
261 while( (zdc = (TRootZdcHits*) itZdc.Next()) )
262 {
263 E = zdc->E; // particle energy in GeV
264 Px = zdc->Px; // particle momentum vector (x component) in GeV
265 Py = zdc->Py; // particle momentum vector (y component) in GeV
266 Pz = zdc->Pz; // particle momentum vector (z component) in GeV
267
268 PT = zdc->PT; // particle transverse momentum in GeV
269 Eta = zdc->Eta; // particle pseudorapidity
270 Phi = zdc->Phi; // particle azimuthal angle in rad
271
272 T = zdc->T; // time of flight [s]
273 side = zdc->side; // -1 or +1
274 }
275
276 //access the RP220 branch
277 TIter itRp220((TCollection*)RP220);
278 TRootRomanPotHits *rp220;
279 itRp220.Reset();
280 //TRootRomanPotHits.MakeClass("test_xav");
281 while( (rp220 = (TRootRomanPotHits*) itRp220.Next()) )
282 {
283 T = rp220->T; // time of flight to the detector [s]
284 S = rp220->S; // distance to the IP [m]
285 E = rp220->E; // reconstructed energy [GeV]
286 q2 = rp220->q2; // reconstructed squared momentum transfer [GeV^2]
287
288 X = rp220->X; // horizontal distance to the beam [um]
289 Y = rp220->Y; // vertical distance to the beam [um]
290
291 Tx = rp220->Tx; // angle of the momentum in the horizontal (x,z) plane [urad]
292 Ty = rp220->Ty; // angle of the momentum in the verical (y,z) plane [urad]
293
294 side = rp220->side; // -1 or 1
295 }
296 //running on FP420 branch is identical
297 TIter itFp420((TCollection*)FP420);
298 TRootRomanPotHits *fp420;
299 itFp420.Reset();
300 while( (fp420 = (TRootRomanPotHits*) itFp420.Next()) ){}
301
302 //*********************************************
303 //Example how to run on the trigger information
304 //*********************************************
305
306 TRootTrigger *trig;
307 int NumTrigBit = TRIGGER->GetEntries();
308 //get the global response of the trigger
309
310 bool GlobalResponse=false;
311 if(NumTrigBit!=0)GlobalResponse=true;
312 cout<<"GlobalResponse "<<GlobalResponse<<endl;
313 for(int i=0; i < NumTrigBit-1; i++){
314 trig = (TRootTrigger*)TRIGGER->At(i);
315 cout<<"The event has been accepted by the trigger number: "<<trig->Accepted<<endl;
316 }
317
318 //********************************
319 //Example of a very small analysis
320 //********************************
321
322 TLorentzVector Lept[2];
323
324 if(NumTrigBit==0)continue; //event not accepted by the trigger
325 cut_trig++;//event accepted
326
327 TSimpleArray<TRootElectron> el=SubArrayEl(ELEC,PT_ELEC);//the central isolated electrons, pt > PT_ELEC GeV
328 TSimpleArray<TRootMuon> mu=SubArrayMu(MUON,PT_MUON);//the central isolated electrons, pt > PT_MUON GeV
329
330 Int_t numElec=el.GetEntries();
331
332 if(el.GetEntries()+mu.GetEntries()!=2)continue;//Exactly 2 isolated leptons are needed
333 cut_1++;//event accepted
334 for(Int_t i=0;i < numElec; i++)Lept[i].SetPxPyPzE(el[i]->Px,el[i]->Py,el[i]->Pz,el[i]->E);
335 for(Int_t k = numElec; k < (numElec+mu.GetEntries()); k++)Lept[k].SetPxPyPzE(mu[k-numElec]->Px,mu[k-numElec]->Py,mu[k-numElec]->Pz,mu[k-numElec]->E);
[264]336 cout<<"normalement il y a quelque chose... "<<endl;
[84]337 //Example how to white a branch in the output file
338 inv_mass=(TRootInvm*) INVMASS->NewEntry();
339 inv_mass->M=(Lept[0]+Lept[1]).M();
340
[81]341 if((Lept[0]+Lept[1]).M() > INV_MASS_LL )continue;// the invariant mass should be < INV_MASS_LL
342 cut_2++;//event accepted
[84]343
344 treeWriter->Fill();
[81]345 }
[84]346 treeWriter->Write();
[81]347
348}
349
350void Analysis_Ex::WriteOutput(string LogName)
351{
352 ofstream f_out(LogName.c_str(),ofstream::app);
353
354 f_out<<"*******************************************************************"<<endl;
355 f_out << left << setw(20) << "Numer of Events "<<""
356 << right << setw(15) << total <<"\n";
357 f_out << left << setw(17) << " Accepted by the trigger "<<""
358 << right << setw(20) << cut_trig <<"\n";
359 f_out << left << setw(17) <<" 2 leptons "<< ""
360 << right << setw(20) << cut_1 << ""
361 << right << setw(15) << cut_1/total << "\n";
362 f_out << left << setw(17) <<" Invariant mass "<< ""
363 << right << setw(20) << cut_2 << ""
364 << right << setw(15) << cut_2/total << "\n";
365 f_out<<"*******************************************************************"<<endl;
366 f_out<<" "<<endl;
367
368}
369
370TSimpleArray<TRootElectron> Analysis_Ex::SubArrayEl(const TClonesArray *ELEC,float pt)
371{
372 TIter itElec((TCollection*)ELEC);
373 TRootElectron *elec;
374 itElec.Reset();
375 TSimpleArray<TRootElectron> array;
376 while( (elec = (TRootElectron*) itElec.Next()) )
377 {
378 if(elec->PT<pt)continue;
379 array.Add(elec);
380 }
381 return array;
382}
383
384TSimpleArray<TRootMuon> Analysis_Ex::SubArrayMu(const TClonesArray *MUON,float pt)
385{
386 TIter itMuon((TCollection*)MUON);
387 TRootMuon *muon;
388 itMuon.Reset();
389 TSimpleArray<TRootMuon> array;
390 while( (muon = (TRootMuon*) itMuon.Next()) )
391 {
392 if(muon->PT<pt)continue;
393 array.Add(muon);
394 }
395 return array;
396}
397
Note: See TracBrowser for help on using the repository browser.