Fork me on GitHub

source: svn/trunk/src/TriggerUtil.cc@ 857

Last change on this file since 857 was 537, checked in by Xavier Rouby, 15 years ago

adding #include<cstdlib> for compatibility with newer versions of g++

File size: 17.3 KB
Line 
1/***********************************************************************
2** **
3** /----------------------------------------------\ **
4** | Delphes, a framework for the fast simulation | **
5** | of a generic collider experiment | **
6** \------------- arXiv:0903.2225v1 ------------/ **
7** **
8** **
9** This package uses: **
10** ------------------ **
11** ROOT: Nucl. Inst. & Meth. in Phys. Res. A389 (1997) 81-86 **
12** FastJet algorithm: Phys. Lett. B641 (2006) [hep-ph/0512210] **
13** Hector: JINST 2:P09005 (2007) [physics.acc-ph:0707.1198v2] **
14** FROG: [hep-ex/0901.2718v1] **
15** HepMC: Comput. Phys. Commun.134 (2001) 41 **
16** **
17** ------------------------------------------------------------------ **
18** **
19** Main authors: **
20** ------------- **
21** **
22** Severine Ovyn Xavier Rouby **
23** severine.ovyn@uclouvain.be xavier.rouby@cern **
24** **
25** Center for Particle Physics and Phenomenology (CP3) **
26** Universite catholique de Louvain (UCL) **
27** Louvain-la-Neuve, Belgium **
28** **
29** Copyright (C) 2008-2009, **
30** All rights reserved. **
31** **
32***********************************************************************/
33
34
35// \brief Trigger class, and some generic definitions
36
37#include "TriggerUtil.h"
38
39#include <iostream>
40#include <iomanip>
41#include <sstream>
42#include <fstream>
43#include <string>
44#include <cstdlib>
45
46using namespace std;
47
48TriggerTable::TriggerTable() {
49 has_been_evaluated=false;
50}
51
52TriggerTable& TriggerTable::operator=(TriggerTable&tt) {
53 if (this == &tt) return *this;
54 has_been_evaluated = tt.has_been_evaluated;
55 GlobalResult = tt.GlobalResult;
56 for (unsigned int i=0; i<tt.list_of_trigger_bits.size(); i++)
57 list_of_trigger_bits.push_back(tt.list_of_trigger_bits[i]);
58 return *this;
59}
60
61TriggerTable::TriggerTable(const TriggerTable& tt) {
62 has_been_evaluated = tt.has_been_evaluated;
63 GlobalResult = tt.GlobalResult;
64 for (unsigned int i=0; i<tt.list_of_trigger_bits.size(); i++)
65 list_of_trigger_bits.push_back(tt.list_of_trigger_bits[i]);
66}
67
68
69void TriggerTable::TriggerCardReader(const string& filename){
70 ifstream trigger_file(filename.c_str());
71 if( !trigger_file.good() ) {
72 cerr << left << setw(30) <<"** ERROR: Trigger input file "<<""
73 << right << setw(15) <<filename<<""
74 << right << setw(24) <<"not found. Exiting... **"<<endl;
75 return;
76 }
77
78 TriggerBit mybit;
79 string temp_string, temp_string_no_comment, trigger_name, trigger_algorithm;
80 while ( getline(trigger_file,temp_string) ) {
81 string temp_string_no_comment = temp_string.substr(0,temp_string.find("#")); //remove comments
82 if(temp_string_no_comment.size()<2) continue; // avoid empty lines
83 istringstream temp_stream(temp_string_no_comment);
84
85 temp_stream >> trigger_name;
86 trigger_algorithm = temp_string_no_comment.substr( temp_string_no_comment.find(">>")+2) ;
87
88 mybit.TriggerBit::GetTrigCondition(trigger_algorithm);
89 list_of_trigger_bits.push_back(mybit);
90 }
91}
92
93
94void TriggerTable::PrintTriggerTable(const string& LogName) {
95
96 ofstream f_out(LogName.c_str(),ios::app);
97 f_out<<"* *"<<"\n";
98 f_out<<"#*********************************************** *"<<"\n";
99 f_out<<"# Trigger conditions defined in the trigger card *"<<"\n";
100 f_out<<"#*********************************************** *"<<"\n";
101 f_out<<"* *"<<"\n";
102
103 f_out.close();
104
105 for(unsigned int i=0; i < list_of_trigger_bits.size(); i++) {
106 list_of_trigger_bits[i].PrintTrigCondition(LogName,i+1);
107 }
108 f_out<<"* *"<<"\n";
109 f_out<<"#*********************************************************************"<<"\n";
110
111}
112
113bool TriggerTable::GetGlobalResult(TClonesArray *branchElecTrig, TClonesArray *branchMuonTrig,
114 TClonesArray *branchJetTrig, TClonesArray *branchTauJetTrig,
115 TClonesArray *branchPhotonTrig, TClonesArray *branchETmisTrig, ExRootTreeBranch *branchTrigger) {
116
117 TRootTrigger *elementTrigger;
118 GlobalResult=false;
119 for(unsigned int i=0; i < list_of_trigger_bits.size(); i++) {
120 if(list_of_trigger_bits[i].GetResult(branchElecTrig, branchMuonTrig, branchJetTrig, branchTauJetTrig, branchPhotonTrig, branchETmisTrig))
121 {
122 GlobalResult=true;
123 int flag=i+1;
124 elementTrigger = (TRootTrigger*) branchTrigger->NewEntry();
125 elementTrigger->Accepted=flag;
126 }
127 }
128 if(GlobalResult)
129 {
130 elementTrigger = (TRootTrigger*) branchTrigger->NewEntry();
131 elementTrigger->Accepted=0;
132 }
133 return GlobalResult;
134
135}
136
137//*************************************************************************
138TriggerBit::TriggerBit() {
139}
140
141TriggerBit::~TriggerBit() {
142 ElecValues.clear();
143 IElecValues.clear();
144 MuonValues.clear();
145 IMuonValues.clear();
146 JetValues.clear();
147 BjetValues.clear();
148 TauJetValues.clear();
149 EtmisValues.clear();
150 GammaValues.clear();
151}
152
153TriggerBit::TriggerBit(const TriggerBit& tb) {
154 for(unsigned int i=0; i<tb.ElecValues.size(); i++)
155 ElecValues.push_back(tb.ElecValues[i]);
156
157 for(unsigned int i=0; i<tb.IElecValues.size(); i++)
158 IElecValues.push_back(tb.IElecValues[i]);
159
160 for(unsigned int i=0; i<tb.MuonValues.size(); i++)
161 MuonValues.push_back(tb.MuonValues[i]);
162
163 for(unsigned int i=0; i<tb.IMuonValues.size(); i++)
164 IMuonValues.push_back(tb.IMuonValues[i]);
165
166 for(unsigned int i=0; i<tb.JetValues.size(); i++)
167 JetValues.push_back(tb.JetValues[i]);
168
169 for(unsigned int i=0; i<tb.BjetValues.size(); i++)
170 BjetValues.push_back(tb.BjetValues[i]);
171
172 for(unsigned int i=0; i<tb.TauJetValues.size(); i++)
173 TauJetValues.push_back(tb.TauJetValues[i]);
174
175 for(unsigned int i=0; i<tb.EtmisValues.size(); i++)
176 EtmisValues.push_back(tb.EtmisValues[i]);
177
178 for(unsigned int i=0; i<tb.GammaValues.size(); i++)
179 GammaValues.push_back(tb.GammaValues[i]);
180
181 Result = tb.Result;
182}
183
184TriggerBit& TriggerBit::operator=(const TriggerBit& tb) {
185 if (this==&tb) return *this;
186 for(unsigned int i=0; i<tb.ElecValues.size(); i++)
187 ElecValues.push_back(tb.ElecValues[i]);
188
189 for(unsigned int i=0; i<tb.IElecValues.size(); i++)
190 IElecValues.push_back(tb.IElecValues[i]);
191
192 for(unsigned int i=0; i<tb.MuonValues.size(); i++)
193 MuonValues.push_back(tb.MuonValues[i]);
194
195 for(unsigned int i=0; i<tb.IMuonValues.size(); i++)
196 IMuonValues.push_back(tb.IMuonValues[i]);
197
198 for(unsigned int i=0; i<tb.JetValues.size(); i++)
199 JetValues.push_back(tb.JetValues[i]);
200
201 for(unsigned int i=0; i<tb.BjetValues.size(); i++)
202 BjetValues.push_back(tb.BjetValues[i]);
203
204 for(unsigned int i=0; i<tb.TauJetValues.size(); i++)
205 TauJetValues.push_back(tb.TauJetValues[i]);
206
207 for(unsigned int i=0; i<tb.EtmisValues.size(); i++)
208 EtmisValues.push_back(tb.EtmisValues[i]);
209
210 for(unsigned int i=0; i<tb.GammaValues.size(); i++)
211 GammaValues.push_back(tb.GammaValues[i]);
212 Result = tb.Result;
213
214 return *this;
215}
216
217
218
219bool TriggerBit::GetResult(TClonesArray *branchElecTrig, TClonesArray *branchMuonTrig,
220 TClonesArray *branchJetTrig, TClonesArray *branchTauJetTrig,
221 TClonesArray *branchPhotonTrig, TClonesArray *branchETmisTrig)
222{
223 TSimpleArray<TRootJet> bjets=SubArrayBjets(branchJetTrig);
224 TSimpleArray<TRootElectron> Ielectron=SubArrayIElec(branchElecTrig);
225 TSimpleArray<TRootMuon> Imuon=SubArrayIMuon(branchMuonTrig);
226
227 int elec_size = ElecValues.size();
228 int Ielec_size = IElecValues.size();
229 int muon_size = MuonValues.size();
230 int Imuon_size = IMuonValues.size();
231 int jet_size = JetValues.size();
232 int bjet_size = BjetValues.size();
233 int taujet_size = TauJetValues.size();
234 int gamma_size = GammaValues.size();
235 int etmis_size = EtmisValues.size();
236
237 Result=false;
238 if(
239 (branchElecTrig->GetEntries() >= elec_size )
240 && (Ielectron.GetEntries() >= Ielec_size )
241 && (branchMuonTrig->GetEntries() >= muon_size )
242 && (Imuon.GetEntries() >= Imuon_size )
243 && (branchJetTrig->GetEntries() >= jet_size )
244 && (bjets.GetEntries() >= bjet_size )
245 && (branchTauJetTrig->GetEntries() >= taujet_size )
246 && (branchPhotonTrig->GetEntries() >= gamma_size )
247 )
248 {
249 Result=true;
250 if(elec_size!=0){
251 TRootElectron *electron;
252 for(int i=0;i<elec_size;i++){
253 electron = (TRootElectron*)branchElecTrig->At(i);
254 if(electron->PT < ElecValues[i])Result=false;}}
255
256 if(Ielec_size!=0)
257 {
258 TRootElectron *i_elec; Ielectron.Reset();
259 int i=0;
260 while((i_elec = Ielectron.Next()))
261 {
262 if(i_elec->PT < IElecValues[i])Result=false;
263 i++;
264 }
265 }
266
267 if(muon_size!=0){
268 TRootMuon *muon;
269 for(int i=0;i<muon_size;i++){
270 muon = (TRootMuon*)branchMuonTrig->At(i);
271 if(muon->PT < MuonValues[i])Result=false;}}
272
273 if(Imuon_size!=0)
274 {
275 TRootMuon *i_muon; Imuon.Reset();
276 int i=0;
277 while((i_muon = Imuon.Next()))
278 {
279 if(i_muon->PT < IMuonValues[i])Result=false;
280 i++;
281 }
282 }
283
284 if(jet_size!=0){
285 TRootJet *jet;
286 for(int i=0;i<jet_size;i++){
287 jet = (TRootJet*)branchJetTrig->At(i);
288 if(jet->PT < JetValues[i])Result=false;}}
289
290 if(bjet_size!=0)
291 {
292 TRootJet *i_bjets; bjets.Reset();
293 int i=0;
294 while((i_bjets = bjets.Next()))
295 {
296 if(i_bjets->PT < BjetValues[i])Result=false;
297 i++;
298 }
299 }
300
301 if(taujet_size!=0){
302 TRootTauJet *taujet;
303 for(int i=0;i<taujet_size;i++){
304 taujet = (TRootTauJet*)branchTauJetTrig->At(i);
305 if(taujet->PT < TauJetValues[i])Result=false;}}
306
307 if(gamma_size!=0){
308 TRootPhoton *gamma;
309 for(int i=0;i<gamma_size;i++){
310 gamma = (TRootPhoton*)branchPhotonTrig->At(i);
311 if(gamma->PT < GammaValues[i])Result=false;}}
312
313
314 if(etmis_size!=0){
315 TRootETmis *etmis = (TRootETmis*)branchETmisTrig->At(0);
316 if(etmis->ET < EtmisValues[0])Result=false;}
317
318 }
319 return Result;
320}
321
322void TriggerBit::GetTrigCondition(const string& trigger_algorithm) {
323
324 vector<string> ElecSequences;
325 vector<string> IElecSequences;
326 vector<string> MuonSequences;
327 vector<string> IMuonSequences;
328 vector<string> JetSequences;
329 vector<string> TauJetSequences;
330 vector<string> EtmisSequences;
331 vector<string> GammaSequences;
332 vector<string> BjetSequences;
333
334 char * result = new char[256];
335 result = strtok( (char*) trigger_algorithm.c_str(),"&");
336 while( result != NULL )
337 {
338 if(strstr (result,"ELEC"))ElecSequences.push_back(result);
339 if(strstr (result,"IElec"))IElecSequences.push_back(result);
340 if(strstr (result,"MUON"))MuonSequences.push_back(result);
341 if(strstr (result,"IMuon"))IMuonSequences.push_back(result);
342 if(strstr (result,"JET"))JetSequences.push_back(result);
343 if(strstr (result,"TAU"))TauJetSequences.push_back(result);
344 if(strstr (result,"ETMIS"))EtmisSequences.push_back(result);
345 if(strstr (result,"GAMMA"))GammaSequences.push_back(result);
346 if(strstr (result,"Bjet"))BjetSequences.push_back(result);
347 result = strtok( NULL,"&");
348 }
349 delete [] result;
350
351 ElecValues = GetCuts(ElecSequences);
352 IElecValues = GetCuts(IElecSequences);
353 MuonValues = GetCuts(MuonSequences);
354 IMuonValues = GetCuts(IMuonSequences);
355 JetValues = GetCuts(JetSequences);
356 TauJetValues = GetCuts(TauJetSequences);
357 EtmisValues = GetCuts(EtmisSequences);
358 GammaValues = GetCuts(GammaSequences);
359 BjetValues = GetCuts(BjetSequences);
360
361}
362
363
364void TriggerBit::PrintTrigCondition(const string& LogName,const int i)
365{
366 int elec_size = TriggerBit::ElecValues.size();
367 int Ielec_size = TriggerBit::IElecValues.size();
368 int muon_size = TriggerBit::MuonValues.size();
369 int Imuon_size = TriggerBit::IMuonValues.size();
370 int jet_size = TriggerBit::JetValues.size();
371 int taujet_size = TriggerBit::TauJetValues.size();
372 int gamma_size = TriggerBit::GammaValues.size();
373 int etmis_size = TriggerBit::EtmisValues.size();
374 int bjets_size = TriggerBit::BjetValues.size();
375
376 ofstream f_out(LogName.c_str(),ios::app);
377
378 f_out <<"* *"<<"\n";
379 f_out << left << setw(25) <<"# Trigger bit number "<<""
380 << left << setw(35) <<i <<""<< right << setw(10)<<"*"<<"\n";
381 f_out <<"#*************************"<<"\n";
382 f_out << left << setw(45) <<"* # Number of reconstructed objects required:"<<""<<right << setw(23)<<"*"<<"\n";
383
384 if(elec_size!=0){
385 f_out << left << setw(5) <<"** -"<<""<< left << setw(2) << elec_size <<""<< left << setw(20)<<"Electron(s) with p_T: "<<"";
386 for(int i=0;i<elec_size;i++){f_out << left << setw(5) << TriggerBit::ElecValues[i]<<"";}f_out <<"\n";}
387
388 if(Ielec_size!=0){
389 f_out << left << setw(5) <<"** -"<<""<< left << setw(2) << Ielec_size <<""<< left << setw(20)<<"Isolated electron(s) with p_T: "<<"";
390 for(int i=0;i<Ielec_size;i++){f_out << left << setw(5) << TriggerBit::IElecValues[i]<<"";}f_out <<"\n";}
391
392 if(muon_size!=0){
393 f_out << left << setw(5) <<"** -"<<""<< left << setw(2) << muon_size <<""<< left << setw(20)<<"Muon(s) with p_T: "<<"";
394 for(int i=0;i<muon_size;i++){f_out << left << setw(5) << TriggerBit::MuonValues[i]<<"";}f_out <<"\n";}
395
396 if(Imuon_size!=0){
397 f_out << left << setw(5) <<"** -"<<""<< left << setw(2) << Imuon_size <<""<< left << setw(20)<<"Isolated muon(s) with p_T: "<<"";
398 for(int i=0;i<Imuon_size;i++){f_out << left << setw(5) << TriggerBit::IMuonValues[i]<<"";}f_out <<"\n";}
399
400 if(jet_size!=0){
401 f_out << left << setw(5) <<"** -"<<""<< left << setw(2) << jet_size <<""<< left << setw(20)<<"Jet(s) with p_T: "<<"";
402 for(int i=0;i<jet_size;i++){f_out << left << setw(5) << TriggerBit::JetValues[i]<<"";}f_out <<"\n";}
403
404 if(taujet_size!=0){
405 f_out << left << setw(5) <<"** -"<<""<< left << setw(2) << taujet_size <<""<< left << setw(20)<<"Tau-jet(s) with p_T: "<<"";
406 for(int i=0;i<taujet_size;i++){f_out << left << setw(5) << TriggerBit::TauJetValues[i]<<"";}f_out <<"\n";}
407
408 if(bjets_size!=0){
409 f_out << left << setw(5) <<"** -"<<""<< left << setw(2) << bjets_size <<""<< left << setw(20)<<"B-jet(s) with p_T: "<<"";
410 for(int i=0;i<bjets_size;i++){f_out << left << setw(5) << TriggerBit::BjetValues[i]<<"";}f_out <<"\n";}
411
412 if(gamma_size!=0){
413 f_out << left << setw(5) <<"** -"<<""<< left << setw(2) << gamma_size <<""<< left << setw(20)<<"Photon(s) with p_T: "<<"";
414 for(int i=0;i<gamma_size;i++){f_out << left << setw(5) << TriggerBit::GammaValues[i]<<"";}f_out <<"\n";}
415
416 if(etmis_size!=0){
417 f_out << left << setw(5) <<"** -"<<""<< left << setw(20)<<"Transverse missing energy with p_T: "<<"";
418 for(int i=0;i<etmis_size;i++){f_out << left << setw(5) << TriggerBit::EtmisValues[i]<<"";}f_out <<"\n";}
419 f_out <<"* *"<<"\n";
420 f_out <<"**********************************************************************"<<"\n";
421
422}
423
424vector<float> TriggerBit::GetCuts(const vector<string> &Sequences)
425{
426 vector<float> OrderedValue;
427 string ptVal;
428 Int_t PtVal=0;
429 for(unsigned int i=0; i < Sequences.size(); i++) {
430 ptVal = Sequences[i].substr( Sequences[i].find("'")+1) ;
431 PtVal = atoi( ptVal.c_str() );
432 OrderedValue.push_back(PtVal);
433 }
434 return OrderedValue;
435
436}
437
438TSimpleArray<TRootJet> TriggerBit::SubArrayBjets(TClonesArray *JET)
439{
440 TIter itJet((TCollection*)JET);
441 TRootJet *jet;
442 itJet.Reset();
443 TSimpleArray<TRootJet> array;
444 while( (jet = (TRootJet*) itJet.Next()) )
445 {
446 if(jet->Btag==1)array.Add(jet);
447 }
448 return array;
449}
450
451TSimpleArray<TRootElectron> TriggerBit::SubArrayIElec(TClonesArray *ELEC)
452{
453 TIter itElec((TCollection*)ELEC);
454 TRootElectron *electron;
455 itElec.Reset();
456 TSimpleArray<TRootElectron> array;
457 while( (electron = (TRootElectron*) itElec.Next()) )
458 {
459 if(electron->IsolFlag==1)array.Add(electron);
460 }
461 return array;
462}
463
464TSimpleArray<TRootMuon> TriggerBit::SubArrayIMuon(TClonesArray *MUON)
465{
466 TIter itMuon((TCollection*)MUON);
467 TRootMuon *muon;
468 itMuon.Reset();
469 TSimpleArray<TRootMuon> array;
470 while( (muon = (TRootMuon*) itMuon.Next()) )
471 {
472 if(muon->IsolFlag==1)array.Add(muon);
473 }
474 return array;
475}
476
Note: See TracBrowser for help on using the repository browser.