CMS 3D CMS Logo

DTTrigTest.cc
Go to the documentation of this file.
1 //-------------------------------------------------
2 //
13 //
14 //--------------------------------------------------
15 
16 // This class's header
18 
19 // Framework related classes
22 
25 
26 // Trigger and DataFormats headers
31 
32 // ROOT headers
33 #include "TROOT.h"
34 #include "TTree.h"
35 #include "TFile.h"
36 
37 // Collaborating classes
39 #include <CLHEP/Vector/LorentzVector.h>
40 
41 // C++ headers
42 #include <iostream>
43 #include <cmath>
44 #include<ctime>
45 
46 using namespace std;
47 using namespace edm;
48 
49 const double DTTrigTest::my_TtoTDC = 32./25.;
50 
52 
53  my_debug= pset.getUntrackedParameter<bool>("debug");
54  string outputfile = pset.getUntrackedParameter<string>("outputFileName");
55  if (my_debug)
56  cout << "[DTTrigTest] Creating rootfile " << outputfile <<endl;
57  my_rootfile = new TFile(outputfile.c_str(),"RECREATE");
58  my_tree = new TTree("h1","GMT",0);
59  my_params = pset;
60  if (my_debug) cout << "[DTTrigTest] Constructor executed!!!" << endl;
61 
62 }
63 
65 
66  if (my_trig != nullptr) delete my_trig;
67  delete my_rootfile;
68  if (my_debug)
69  cout << "[DTTrigTest] Destructor executed!!!" << endl;
70 
71 }
72 
74 
75  if (my_debug)
76  cout << "[DTTrigTest] Writing Tree and Closing File" << endl;
77  my_tree->Write();
78  delete my_tree;
79  my_rootfile->Close();
80 
81 }
82 
83 //void DTTrigTest::beginJob(const EventSetup & iEventSetup){
85  // get DTConfigManager
86  // ESHandle< DTConfigManager > confManager ;
87  // iEventSetup.get< DTConfigManagerRcd >().get( confManager ) ;
88 
89  //for testing purpose....
90  //DTBtiId btiid(1,1,1,1,1);
91  //confManager->getDTConfigBti(btiid)->print();
92 
93 // my_trig = new DTTrig(my_params);
94 
95 // my_trig->createTUs(iEventSetup);
96 // if (my_debug)
97 // cout << "[DTTrigTest] TU's Created" << endl;
98 
99  // BOOKING of the tree's varables
100  // GENERAL block branches
101  my_tree->Branch("Run",&runn,"Run/I");
102  my_tree->Branch("Event",&eventn,"Event/I");
103  my_tree->Branch("Weight",&weight,"Weight/F");
104  // GEANT block branches
105  my_tree->Branch("Ngen",&ngen,"Ngen/I");
106  my_tree->Branch("Pxgen",pxgen,"Pxgen[Ngen]/F");
107  my_tree->Branch("Pygen",pygen,"Pygen[Ngen]/F");
108  my_tree->Branch("Pzgen",pzgen,"Pzgen[Ngen]/F");
109  my_tree->Branch("Ptgen",ptgen,"Ptgen[Ngen]/F");
110  my_tree->Branch("Etagen",etagen,"Etagen[Ngen]/F");
111  my_tree->Branch("Phigen",phigen,"Phigen[Ngen]/F");
112  my_tree->Branch("Chagen",chagen,"Chagen[Ngen]/I");
113  my_tree->Branch("Vxgen",vxgen,"Vxgen[Ngen]/F");
114  my_tree->Branch("Vygen",vygen,"Vygen[Ngen]/F");
115  my_tree->Branch("Vzgen",vzgen,"Vzgen[Ngen]/F");
116  // L1MuDTBtiChipS block
117  my_tree->Branch("Nbti",&nbti,"Nbti/I");
118  my_tree->Branch("bwh",bwh,"bwh[Nbti]/I");
119  my_tree->Branch("bstat",bstat,"bstat[Nbti]/I");
120  my_tree->Branch("bsect",bsect,"bsect[Nbti]/I");
121  my_tree->Branch("bsl",bsl,"bsl[Nbti]/I");
122  my_tree->Branch("bnum",bnum,"bnum[Nbti]/I");
123  my_tree->Branch("bbx",bbx,"bbx[Nbti]/I");
124  my_tree->Branch("bcod",bcod,"bcod[Nbti]/I");
125  my_tree->Branch("bk",bk,"bk[Nbti]/I");
126  my_tree->Branch("bx",bx,"bx[Nbti]/I");
127  my_tree->Branch("bposx",bposx,"bposx[Nbti]/F");
128  my_tree->Branch("bposy",bposy,"bposy[Nbti]/F");
129  my_tree->Branch("bposz",bposz,"bposz[Nbti]/F");
130  my_tree->Branch("bdirx",bdirx,"bdirx[Nbti]/F");
131  my_tree->Branch("bdiry",bdiry,"bdiry[Nbti]/F");
132  my_tree->Branch("bdirz",bdirz,"bdirz[Nbti]/F");
133  // L1MuDTTracoChipS block
134  my_tree->Branch("Ntraco",&ntraco,"Ntraco/I");
135  my_tree->Branch("twh",twh,"twh[Ntraco]/I");
136  my_tree->Branch("tstat",tstat,"tstat[Ntraco]/I");
137  my_tree->Branch("tsect",tsect,"tsect[Ntraco]/I");
138  my_tree->Branch("tnum",tnum,"tnum[Ntraco]/I");
139  my_tree->Branch("tbx",tbx,"tbx[Ntraco]/I");
140  my_tree->Branch("tcod",tcod,"tcod[Ntraco]/I");
141  my_tree->Branch("tk",tk,"tk[Ntraco]/I");
142  my_tree->Branch("tx",tx,"tx[Ntraco]/I");
143  my_tree->Branch("tposx",tposx,"tposx[Ntraco]/F");
144  my_tree->Branch("tposy",tposy,"tposy[Ntraco]/F");
145  my_tree->Branch("tposz",tposz,"tposz[Ntraco]/F");
146  my_tree->Branch("tdirx",tdirx,"tdirx[Ntraco]/F");
147  my_tree->Branch("tdiry",tdiry,"tdiry[Ntraco]/F");
148  my_tree->Branch("tdirz",tdirz,"tdirz[Ntraco]/F");
149  // TSPHI block
150  my_tree->Branch("Ntsphi",&ntsphi,"Ntsphi/I");
151  my_tree->Branch("swh",swh,"swh[Ntsphi]/I");
152  my_tree->Branch("sstat",sstat,"sstat[Ntsphi]/I");
153  my_tree->Branch("ssect",ssect,"ssect[Ntsphi]/I");
154  my_tree->Branch("sbx",sbx,"sbx[Ntsphi]/I");
155  my_tree->Branch("scod",scod,"scod[Ntsphi]/I");
156  my_tree->Branch("sphi",sphi,"sphi[Ntsphi]/I");
157  my_tree->Branch("sphib",sphib,"sphib[Ntsphi]/I");
158  my_tree->Branch("sposx",sposx,"sposx[Ntsphi]/F");
159  my_tree->Branch("sposy",sposy,"sposy[Ntsphi]/F");
160  my_tree->Branch("sposz",sposz,"sposz[Ntsphi]/F");
161  my_tree->Branch("sdirx",sdirx,"sdirx[Ntsphi]/F");
162  my_tree->Branch("sdiry",sdiry,"sdiry[Ntsphi]/F");
163  my_tree->Branch("sdirz",sdirz,"sdirz[Ntsphi]/F");
164  // TSTHETA block
165  my_tree->Branch("Ntstheta",&ntstheta,"Ntstheta/I");
166  my_tree->Branch("thwh",thwh,"thwh[Ntstheta]/I");
167  my_tree->Branch("thstat",thstat,"thstat[Ntstheta]/I");
168  my_tree->Branch("thsect",thsect,"thsect[Ntstheta]/I");
169  my_tree->Branch("thbx",thbx,"thbx[Ntstheta]/I");
170  my_tree->Branch("thcode",thcode,"thcode[Ntstheta][7]/I");
171  my_tree->Branch("thpos",thpos,"thpos[Ntstheta][7]/I");
172  my_tree->Branch("thqual",thqual,"thqual[Ntstheta][7]/I");
173  // SC PHI block
174  my_tree->Branch("Nscphi",&nscphi,"Nscphi/I");
175  my_tree->Branch("scphwh",scphwh,"scphwh[Nscphi]/I");
176  my_tree->Branch("scphstat",scphstat,"scphstat[Nscphi]/I");
177  my_tree->Branch("scphsect",scphsect,"scphsect[Nscphi]/I");
178  my_tree->Branch("scphbx",scphbx,"scphbx[Nscphi]/I");
179  my_tree->Branch("scphcod",scphcod,"scphcod[Nscphi]/I");
180  my_tree->Branch("scphphi",scphphi,"scphphi[Nscphi]/I");
181  my_tree->Branch("scphphib",scphphib,"scphphib[Nscphi]/I");
182  my_tree->Branch("scphposx",scphposx,"scphposx[Nscphi]/F");
183  my_tree->Branch("scphposy",scphposy,"scphposy[Nscphi]/F");
184  my_tree->Branch("scphposz",scphposz,"scphposz[Nscphi]/F");
185  my_tree->Branch("scphdirx",scphdirx,"scphdirx[Nscphi]/F");
186  my_tree->Branch("scphdiry",scphdiry,"scphdiry[Nscphi]/F");
187  my_tree->Branch("scphdirz",scphdirz,"scphdirz[Nscphi]/F");
188  // SC THETA block
189  my_tree->Branch("Nsctheta",&nsctheta,"Nsctheta/I");
190  my_tree->Branch("scthwh",scthwh,"scthwh[Nsctheta]/I");
191  my_tree->Branch("scthstat",scthstat,"scthstat[Nsctheta]/I");
192  my_tree->Branch("scthsect",scthsect,"scthsect[Nsctheta]/I");
193  my_tree->Branch("scthbx",scthbx,"scthbx[Nsctheta]/I");
194  my_tree->Branch("scthcode",scthcode,"scthcode[Nsctheta][7]/I");
195  my_tree->Branch("scthpos",scthpos,"scthpos[Nsctheta][7]/I");
196  my_tree->Branch("scthqual",scthqual,"scthqual[Nsctheta][7]/I");
197 
198 }
199 
200 void DTTrigTest::beginRun(const edm::Run& iRun, const edm::EventSetup& iEventSetup) {
201 
202  if (!my_trig) {
204  my_trig->createTUs(iEventSetup);
205  if (my_debug)
206  cout << "[DTTrigTest] TU's Created" << endl;
207  }
208 
209 }
210 
211 
212 
213 void DTTrigTest::analyze(const Event & iEvent, const EventSetup& iEventSetup){
214 
215  const int MAXGEN = 10;
216  const float ptcut = 1.0;
217  const float etacut = 2.4;
218  my_trig->triggerReco(iEvent,iEventSetup);
219  if (my_debug)
220  cout << "[DTTrigTest] Trigger algorithm executed for run " << iEvent.id().run()
221  <<" event " << iEvent.id().event() << endl;
222 
223  // GENERAL Block
224  runn = iEvent.id().run();
225  eventn = iEvent.id().event();
226  weight = 1; // FIXME what to do with this varable?
227 
228  // GEANT Block
229  Handle<vector<SimTrack> > MyTracks;
230  Handle<vector<SimVertex> > MyVertexes;
231  iEvent.getByLabel("g4SimHits",MyTracks);
232  iEvent.getByLabel("g4SimHits",MyVertexes);
233  vector<SimTrack>::const_iterator itrack;
234  ngen=0;
235  if (my_debug)
236  cout << "[DTTrigTest] Tracks found in the detector (not only muons) " << MyTracks->size() <<endl;
237 
238  for (itrack=MyTracks->begin(); itrack!=MyTracks->end(); itrack++){
239  if ( abs(itrack->type())==13){
240  math::XYZTLorentzVectorD momentum = itrack->momentum();
241  float pt = momentum.Pt();
242  float eta = momentum.eta();
243  if ( pt>ptcut && fabs(eta)<etacut ){
244  float phi = momentum.phi();
245  int charge = static_cast<int> (-itrack->type()/13); //static_cast<int> (itrack->charge());
246  if ( phi<0 ) phi = 2*M_PI + phi;
247  int vtxindex = itrack->vertIndex();
248  float gvx=0,gvy=0,gvz=0;
249  if (vtxindex >-1){
250  gvx=MyVertexes->at(vtxindex).position().x();
251  gvy=MyVertexes->at(vtxindex).position().y();
252  gvz=MyVertexes->at(vtxindex).position().z();
253  }
254  if ( ngen < MAXGEN ) {
255  pxgen[ngen]=momentum.x();
256  pygen[ngen]=momentum.y();
257  pzgen[ngen]=momentum.z();
258  ptgen[ngen]=pt;
259  etagen[ngen]=eta;
260  phigen[ngen]=phi;
261  chagen[ngen]=charge;
262  vxgen[ngen]=gvx;
263  vygen[ngen]=gvy;
264  vzgen[ngen]=gvz;
265  ngen++;
266  }
267  }
268  }
269  }
270 
271  // L1 Local Trigger Block
272  // BTI
273  vector<DTBtiTrigData> btitrigs = my_trig->BtiTrigs();
274  vector<DTBtiTrigData>::const_iterator pbti;
275  int ibti = 0;
276  if (my_debug)
277  cout << "[DTTrigTest] " << btitrigs.size() << " BTI triggers found" << endl;
278 
279  for ( pbti = btitrigs.begin(); pbti != btitrigs.end(); pbti++ ) {
280  if ( ibti < 100 ) {
281  bwh[ibti]=pbti->wheel();
282  bstat[ibti]=pbti->station();
283  bsect[ibti]=pbti->sector();
284  bsl[ibti]=pbti->btiSL();
285  bnum[ibti]=pbti->btiNumber();
286  bbx[ibti]=pbti->step();
287  bcod[ibti]=pbti->code();
288  bk[ibti]=pbti->K();
289  bx[ibti]=pbti->X();
290  GlobalPoint pos = my_trig->CMSPosition(&(*pbti));
291  GlobalVector dir = my_trig->CMSDirection(&(*pbti));
292  bposx[ibti] = pos.x();
293  bposy[ibti] = pos.y();
294  bposz[ibti] = pos.z();
295  bdirx[ibti] = dir.x();
296  bdiry[ibti] = dir.y();
297  bdirz[ibti] = dir.z();
298  ibti++;
299  }
300  }
301  nbti = ibti;
302 
303  //TRACO
304  vector<DTTracoTrigData> tracotrigs = my_trig->TracoTrigs();
305  vector<DTTracoTrigData>::const_iterator ptc;
306  int itraco = 0;
307  if (my_debug)
308  cout << "[DTTrigTest] " << tracotrigs.size() << " TRACO triggers found" << endl;
309 
310  for (ptc=tracotrigs.begin(); ptc!=tracotrigs.end(); ptc++) {
311  if (itraco<80) {
312  twh[itraco]=ptc->wheel();
313  tstat[itraco]=ptc->station();
314  tsect[itraco]=ptc->sector();
315  tnum[itraco]=ptc->tracoNumber();
316  tbx[itraco]=ptc->step();
317  tcod[itraco]=ptc->code();
318  tk[itraco]=ptc->K();
319  tx[itraco]=ptc->X();
320  GlobalPoint pos = my_trig->CMSPosition(&(*ptc));
321  GlobalVector dir = my_trig->CMSDirection(&(*ptc));
322  tposx[itraco] = pos.x();
323  tposy[itraco] = pos.y();
324  tposz[itraco] = pos.z();
325  tdirx[itraco] = dir.x();
326  tdiry[itraco] = dir.y();
327  tdirz[itraco] = dir.z();
328  itraco++;
329  }
330  }
331  ntraco = itraco;
332 
333  //TSPHI
334  vector<DTChambPhSegm> tsphtrigs = my_trig->TSPhTrigs();
335  vector<DTChambPhSegm>::const_iterator ptsph;
336  int itsphi = 0;
337  if (my_debug )
338  cout << "[DTTrigTest] " << tsphtrigs.size() << " TSPhi triggers found" << endl;
339 
340  for (ptsph=tsphtrigs.begin(); ptsph!=tsphtrigs.end(); ptsph++) {
341  if (itsphi<40 ) {
342  swh[itsphi] = ptsph->wheel();
343  sstat[itsphi] = ptsph->station();
344  ssect[itsphi] = ptsph->sector();
345  sbx[itsphi] = ptsph->step();
346  scod[itsphi] = ptsph->oldCode();
347  sphi[itsphi] = ptsph->phi();
348  sphib[itsphi] = ptsph->phiB();
349  GlobalPoint pos = my_trig->CMSPosition(&(*ptsph));
350  GlobalVector dir = my_trig->CMSDirection(&(*ptsph));
351  sposx[itsphi] = pos.x();
352  sposy[itsphi] = pos.y();
353  sposz[itsphi] = pos.z();
354  sdirx[itsphi] = dir.x();
355  sdiry[itsphi] = dir.y();
356  sdirz[itsphi] = dir.z();
357  itsphi++;
358  }
359  }
360  ntsphi = itsphi;
361 
362  //TSTHETA
363  vector<DTChambThSegm> tsthtrigs = my_trig->TSThTrigs();
364  vector<DTChambThSegm>::const_iterator ptsth;
365  int itstheta = 0;
366  if (my_debug)
367  cout << "[DTTrigTest] " << tsthtrigs.size() << " TSTheta triggers found" << endl;
368 
369  for (ptsth=tsthtrigs.begin(); ptsth!=tsthtrigs.end(); ptsth++) {
370  if (itstheta<40 ) {
371  thwh[itstheta] = ptsth->ChamberId().wheel();
372  thstat[itstheta] = ptsth->ChamberId().station();
373  thsect[itstheta] = ptsth->ChamberId().sector();
374  thbx[itstheta] = ptsth->step();
375  for(int i=0;i<7;i++) {
376  thcode[itstheta][i] = ptsth->code(i);
377  thpos[itstheta][i] = ptsth->position(i);
378  thqual[itstheta][i] = ptsth->quality(i);
379  }
380  itstheta++;
381  }
382  }
383  ntstheta = itstheta;
384 
385  //SCPHI
386  vector<DTSectCollPhSegm> scphtrigs = my_trig->SCPhTrigs();
387  vector<DTSectCollPhSegm>::const_iterator pscph;
388  int iscphi = 0;
389  if (my_debug )
390  cout << "[DTTrigTest] " << scphtrigs.size() << " SectCollPhi triggers found" << endl;
391 
392  for (pscph=scphtrigs.begin(); pscph!=scphtrigs.end(); pscph++) {
393  if (iscphi<40 ) {
394  const DTChambPhSegm *seg = (*pscph).tsPhiTrig();
395  scphwh[iscphi] = pscph->wheel();
396  scphstat[iscphi] = pscph->station();
397  scphsect[iscphi] = pscph->sector();
398  scphbx[iscphi] = pscph->step();
399  scphcod[iscphi] = pscph->oldCode();
400  scphphi[iscphi] = pscph->phi();
401  scphphib[iscphi] = pscph->phiB();
404  scphposx[iscphi] = pos.x();
405  scphposy[iscphi] = pos.y();
406  scphposz[iscphi] = pos.z();
407  scphdirx[iscphi] = dir.x();
408  scphdiry[iscphi] = dir.y();
409  scphdirz[iscphi] = dir.z();
410  iscphi++;
411  }
412  }
413  nscphi = iscphi;
414 
415  //SCTHETA
416  vector<DTSectCollThSegm> scthtrigs = my_trig->SCThTrigs();
417  vector<DTSectCollThSegm>::const_iterator pscth;
418  int isctheta = 0;
419  if (my_debug)
420  cout << "[DTTrigTest] " << scthtrigs.size() << " SectCollTheta triggers found" << endl;
421 
422  for (pscth=scthtrigs.begin(); pscth!=scthtrigs.end(); pscth++) {
423  if (isctheta<40 ) {
424  scthwh[isctheta] = pscth->ChamberId().wheel();
425  scthstat[isctheta] = pscth->ChamberId().station();
426  scthsect[isctheta] = pscth->ChamberId().sector();
427  scthbx[isctheta] = pscth->step();
428  for(int i=0;i<7;i++) {
429  scthcode[isctheta][i] = pscth->code(i);
430  scthpos[isctheta][i] = pscth->position(i);
431  scthqual[isctheta][i] = pscth->quality(i);
432  }
433  isctheta++;
434  }
435  }
436  nsctheta = isctheta;
437 
438  //Fill the tree
439  my_tree->Fill();
440 
441 }
RunNumber_t run() const
Definition: EventID.h:39
float scphdirx[100]
Definition: DTTrigTest.h:165
int thpos[40][7]
Definition: DTTrigTest.h:150
EventNumber_t event() const
Definition: EventID.h:41
T getUntrackedParameter(std::string const &, T const &) const
int tstat[80]
Definition: DTTrigTest.h:113
int sbx[40]
Definition: DTTrigTest.h:132
int scphbx[40]
Definition: DTTrigTest.h:158
int scthsect[40]
Definition: DTTrigTest.h:173
edm::ParameterSet my_params
Definition: DTTrigTest.h:66
float scphposz[100]
Definition: DTTrigTest.h:164
float etagen[10]
Definition: DTTrigTest.h:85
float tdirz[100]
Definition: DTTrigTest.h:125
ROOT::Math::LorentzVector< ROOT::Math::PxPyPzE4D< double > > XYZTLorentzVectorD
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:14
int bstat[100]
Definition: DTTrigTest.h:95
int bsl[100]
Definition: DTTrigTest.h:97
int bnum[100]
Definition: DTTrigTest.h:98
int thstat[40]
Definition: DTTrigTest.h:146
int twh[80]
Definition: DTTrigTest.h:112
float vygen[10]
Definition: DTTrigTest.h:89
int scod[40]
Definition: DTTrigTest.h:133
float tdiry[100]
Definition: DTTrigTest.h:124
int scthbx[40]
Definition: DTTrigTest.h:174
int ssect[40]
Definition: DTTrigTest.h:131
int sphib[40]
Definition: DTTrigTest.h:135
float bdirz[100]
Definition: DTTrigTest.h:108
int tcod[80]
Definition: DTTrigTest.h:117
float vxgen[10]
Definition: DTTrigTest.h:88
int scthpos[40][7]
Definition: DTTrigTest.h:176
int scphphib[40]
Definition: DTTrigTest.h:161
float pxgen[10]
Definition: DTTrigTest.h:81
float scphdirz[100]
Definition: DTTrigTest.h:167
int scphcod[40]
Definition: DTTrigTest.h:159
T y() const
Definition: PV3DBase.h:63
void beginRun(const edm::Run &iRun, const edm::EventSetup &iEventSetup) override
Create DTTrig instance and TUs.
Definition: DTTrigTest.cc:200
float scphposx[100]
Definition: DTTrigTest.h:162
int bx[100]
Definition: DTTrigTest.h:102
float scphposy[100]
Definition: DTTrigTest.h:163
Definition: weight.py:1
int scphwh[40]
Definition: DTTrigTest.h:155
int bwh[100]
Definition: DTTrigTest.h:94
float sposx[100]
Definition: DTTrigTest.h:136
#define nullptr
int sstat[40]
Definition: DTTrigTest.h:130
float ptgen[10]
Definition: DTTrigTest.h:84
void endJob() override
Close Tree and write File.
Definition: DTTrigTest.cc:73
int bcod[100]
Definition: DTTrigTest.h:100
int bsect[100]
Definition: DTTrigTest.h:96
float sposy[100]
Definition: DTTrigTest.h:137
int tnum[80]
Definition: DTTrigTest.h:115
std::vector< DTSectCollPhSegm > SCPhTrigs()
Return a copy of all the Sector Collector (Phi) triggers.
Definition: DTTrig.cc:553
int thbx[40]
Definition: DTTrigTest.h:148
int iEvent
Definition: GenABIO.cc:230
float sposz[100]
Definition: DTTrigTest.h:138
float tposz[100]
Definition: DTTrigTest.h:122
void beginJob() override
Create tree and Branches.
Definition: DTTrigTest.cc:84
int thwh[40]
Definition: DTTrigTest.h:145
ConsumesCollector consumesCollector()
Use a ConsumesCollector to gather consumes information from helper functions.
int scphphi[40]
Definition: DTTrigTest.h:160
T z() const
Definition: PV3DBase.h:64
float scphdiry[100]
Definition: DTTrigTest.h:166
int scthwh[40]
Definition: DTTrigTest.h:171
int scthqual[40][7]
Definition: DTTrigTest.h:177
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
static const double my_TtoTDC
Definition: DTTrigTest.h:57
int tk[80]
Definition: DTTrigTest.h:118
int scphsect[40]
Definition: DTTrigTest.h:157
int swh[40]
Definition: DTTrigTest.h:129
TTree * my_tree
Definition: DTTrigTest.h:69
int tx[80]
Definition: DTTrigTest.h:119
int bbx[100]
Definition: DTTrigTest.h:99
float sdirz[100]
Definition: DTTrigTest.h:141
void triggerReco(const edm::Event &iEvent, const edm::EventSetup &iSetup)
Run the whole trigger reconstruction chain.
Definition: DTTrig.cc:152
int nsctheta
Definition: DTTrigTest.h:170
const int MAXGEN
int thcode[40][7]
Definition: DTTrigTest.h:149
GlobalPoint CMSPosition(const DTTrigData *trig) const
Coordinate of a trigger-data object in CMS frame.
Definition: DTTrig.h:211
int ntstheta
Definition: DTTrigTest.h:144
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:535
GlobalVector CMSDirection(const DTTrigData *trig) const
Direction of a trigger-data object in CMS frame.
Definition: DTTrig.h:221
#define M_PI
float bdirx[100]
Definition: DTTrigTest.h:106
float tdirx[100]
Definition: DTTrigTest.h:123
int scphstat[40]
Definition: DTTrigTest.h:156
float tposx[100]
Definition: DTTrigTest.h:120
std::vector< DTChambThSegm > TSThTrigs()
Return a copy of all the Trigger Server (Theta) triggers.
Definition: DTTrig.cc:537
Definition: DTTrig.h:53
float sdiry[100]
Definition: DTTrigTest.h:140
float pygen[10]
Definition: DTTrigTest.h:82
std::vector< DTTracoTrigData > TracoTrigs()
Return a copy of all the TRACO triggers.
Definition: DTTrig.cc:505
float bposx[100]
Definition: DTTrigTest.h:103
~DTTrigTest() override
Destructor.
Definition: DTTrigTest.cc:64
TFile * my_rootfile
Definition: DTTrigTest.h:71
int tsect[80]
Definition: DTTrigTest.h:114
float phigen[10]
Definition: DTTrigTest.h:86
void analyze(const edm::Event &iEvent, const edm::EventSetup &iEventSetup) override
Analyze function executed on all the events.
Definition: DTTrigTest.cc:213
int scthstat[40]
Definition: DTTrigTest.h:172
DTTrigTest(const edm::ParameterSet &pset)
Constructor.
Definition: DTTrigTest.cc:51
int sphi[40]
Definition: DTTrigTest.h:134
edm::EventID id() const
Definition: EventBase.h:60
float bdiry[100]
Definition: DTTrigTest.h:107
HLT enums.
int thsect[40]
Definition: DTTrigTest.h:147
float sdirx[100]
Definition: DTTrigTest.h:139
int scthcode[40][7]
Definition: DTTrigTest.h:175
float pzgen[10]
Definition: DTTrigTest.h:83
int eventn
Definition: DTTrigTest.h:76
float bposz[100]
Definition: DTTrigTest.h:105
dbl *** dir
Definition: mlp_gen.cc:35
float tposy[100]
Definition: DTTrigTest.h:121
int chagen[10]
Definition: DTTrigTest.h:87
float bposy[100]
Definition: DTTrigTest.h:104
T x() const
Definition: PV3DBase.h:62
void createTUs(const edm::EventSetup &iSetup)
Create the trigger units and store them in the cache.
Definition: DTTrig.cc:76
float vzgen[10]
Definition: DTTrigTest.h:90
bool my_debug
Definition: DTTrigTest.h:63
Definition: Run.h:44
DTTrig * my_trig
Definition: DTTrigTest.h:60
std::vector< DTBtiTrigData > BtiTrigs()
Return a copy of all the BTI triggers.
Definition: DTTrig.cc:489
std::vector< DTSectCollThSegm > SCThTrigs()
Return a copy of all the Sector Collector (Theta) triggers.
Definition: DTTrig.cc:582
int thqual[40][7]
Definition: DTTrigTest.h:151
std::vector< DTChambPhSegm > TSPhTrigs()
Return a copy of all the Trigger Server (Phi) triggers.
Definition: DTTrig.cc:521
int tbx[80]
Definition: DTTrigTest.h:116