CMS 3D CMS Logo

DTTrig.cc
Go to the documentation of this file.
1 //-------------------------------------------------
2 //
3 // Class: DTTrig
4 //
5 // Description: Steering routine for L1 trigger simulation in
6 // a muon barrel station
7 //
8 //
9 // Author List:
10 // C. Grandi
11 // Modifications:
12 // S Vanini, S. Marcellini, D. Bonacorsi, C.Battilana
13 //
14 // 07/03/30 : configuration now through DTConfigManager SV
15 //-------------------------------------------------------
16 
17 //-----------------------
18 // This Class's Header --
19 //-----------------------
21 
27 
30 
31 //---------------
32 // C++ Headers --
33 //---------------
34 #include <iostream>
35 #include <string>
36 
37 //-------------------------------
38 // Collaborating Class Headers --
39 //-------------------------------
41 
42 //----------------
43 // Constructors --
44 //----------------
45 
47  : _inputexist(true), _configid(0), _geomid(0) {
48  // Set configuration parameters
49  _debug = params.getUntrackedParameter<bool>("debug");
50 
51  if (_debug) {
52  std::cout << std::endl;
53  std::cout << "**** Initialization of DTTrigger ****" << std::endl;
54  }
55 
56  _digitag = params.getParameter<edm::InputTag>("digiTag");
57  iC.consumes<DTDigiCollection>(_digitag);
58  dtGeomToken_ = iC.esConsumes<DTGeometry, MuonGeometryRecord>();
59  confToken_ = iC.esConsumes<DTConfigManager, DTConfigManagerRcd>();
61 }
62 
63 void DTTrig::createTUs(const edm::EventSetup& iSetup) {
64  // build up Sector Collectors and then
65  // build the trrigger units (one for each chamber)
66  for (int iwh = -2; iwh <= 2; iwh++) {
67  for (int ise = 1; ise <= 12; ise++) {
68  if (_debug) {
69  std::cout << "calling sectcollid wh sc " << iwh << " " << ise << std::endl;
70  }
71  DTSectCollId scid(iwh, ise);
72  {
73  SC_iterator it = _cache1.find(scid);
74  if (it != _cache1.end()) {
75  if (_debug) {
76  std::cout << "DTTrig::createTUs: Sector Collector unit already exists" << std::endl;
77  }
78  continue;
79  }
80  }
81  {
82  auto element = _cache1.emplace(scid, scid);
83  if (_debug) {
84  std::cout << " DTTrig::createTUs new SC sc = " << &(element.first->second) << " at scid.sector() "
85  << scid.sector() << " at scid.wheel() " << scid.wheel() << std::endl;
86  }
87  }
88  }
89  }
90 
92  for (std::vector<const DTChamber*>::const_iterator ich = dtGeom->chambers().begin(); ich != dtGeom->chambers().end();
93  ich++) {
94  const DTChamber* chamb = (*ich);
95  DTChamberId chid = chamb->id();
96  TU_iterator it = _cache.find(chid);
97  if (it != _cache.end()) {
98  if (_debug)
99  std::cout << "DTTrig::init: Trigger unit already exists" << std::endl;
100  continue;
101  }
102 
103  auto info = _cache.emplace(chid, chamb);
104  auto tru = &(info.first->second);
105 
106  //----------- add TU to corresponding SC
107  // returning correspondent SC id
108  DTSectCollId scid;
109  if (chid.sector() == 13) {
110  scid = DTSectCollId(chid.wheel(), 4);
111  } else if (chid.sector() == 14) {
112  scid = DTSectCollId(chid.wheel(), 10);
113  } else {
114  scid = DTSectCollId(chid.wheel(), chid.sector());
115  }
116 
117  SC_iterator it1 = _cache1.find(scid);
118 
119  if (it1 != _cache1.end()) {
120  auto& sc = (*it1).second;
121  if (_debug) {
122  std::cout << "DTTrig::init: adding TU in SC << "
123  << " sector = " << scid.sector() << " wheel = " << scid.wheel() << std::endl;
124  }
125  sc.addTU(tru);
126  } else {
127  std::cout << "DTTrig::createTUs: Trigger Unit not in the map: ";
128  }
129  }
130 }
131 
133  updateES(iSetup);
134  if (!_inputexist)
135  return;
136 
137  DTDigiMap digiMap;
138  //Sort digis by chamber so they can be used by BTIs
140  iEvent.getByLabel(_digitag, dtDigis);
141 
142  if (!dtDigis.isValid()) {
143  LogDebug("DTTrig") << "DTTrig::triggerReco DTDigiCollection with input tag " << _digitag
144  << "requested in configuration, but not found in the event." << std::endl;
145  _inputexist = false;
146  return;
147  }
148 
150 
151  for (detUnitIt = dtDigis->begin(); detUnitIt != dtDigis->end(); ++detUnitIt) {
152  const DTLayerId& layId = (*detUnitIt).first;
153  const DTChamberId chambId = layId.superlayerId().chamberId();
154  const DTDigiCollection::Range& range = (*detUnitIt).second;
155  digiMap[chambId].put(range, layId);
156  }
157 
158  //Run reconstruct for single trigger subsystem (Bti, Traco TS)
159  for (TU_iterator it = _cache.begin(); it != _cache.end(); it++) {
160  DTSCTrigUnit& thisTU = (*it).second;
161  if (thisTU.BtiTrigs()->size() > 0) {
162  thisTU.BtiTrigs()->clearCache();
163  thisTU.TSThTrigs()->clearCache();
164  thisTU.TracoTrigs()->clearCache();
165  thisTU.TSPhTrigs()->clearCache();
166  }
167  DTChamberId chid = thisTU.statId();
168  DTDigiMap_iterator dmit = digiMap.find(chid);
169  if (dmit != digiMap.end()) {
170  thisTU.BtiTrigs()->reconstruct((*dmit).second);
171  if (thisTU.BtiTrigs()->size() > 0) {
172  thisTU.TSThTrigs()->reconstruct();
173  thisTU.TracoTrigs()->reconstruct();
174  if (thisTU.TracoTrigs()->size() > 0)
175  thisTU.TSPhTrigs()->reconstruct();
176  }
177  }
178  }
179  //Run reconstruct for Sector Collector
180  for (SC_iterator it = _cache1.begin(); it != _cache1.end(); it++) {
181  DTSectColl& sectcoll = (*it).second;
182  DTSectCollId scid = (*it).first;
183  if (sectcoll.sizePh() > 0 || sectcoll.sizeTh() > 0)
184  sectcoll.clearCache();
185  bool mustreco = false;
186  for (int i = 1; i < 5; i++) {
187  if (sectcoll.getTSPhi(i)->size() > 0) {
188  mustreco = true;
189  break;
190  }
191  }
192  for (int i = 1; i < 4; i++) {
193  if (sectcoll.getTSTheta(i)->size() > 0) {
194  mustreco = true;
195  break;
196  }
197  }
198  if (scid.sector() == 4 || scid.sector() == 10) {
199  if (sectcoll.getTSPhi(5)->size() > 0)
200  mustreco = true;
201  }
202  if (mustreco)
203  sectcoll.reconstruct();
204  }
205 }
206 
207 void DTTrig::updateES(const edm::EventSetup& iSetup) {
208  // Check for updatets in config
210  edm::ESHandle<DTGeometry> geomHandle;
211 
212  if (iSetup.get<DTConfigManagerRcd>().cacheIdentifier() != _configid) {
213  if (_debug)
214  std::cout << "DTTrig::updateES updating DTTPG configuration" << std::endl;
215 
216  _configid = iSetup.get<DTConfigManagerRcd>().cacheIdentifier();
217  confHandle = iSetup.getHandle(confToken_);
218  _conf_manager = confHandle.product();
219  for (TU_iterator it = _cache.begin(); it != _cache.end(); it++) {
220  (*it).second.setConfig(_conf_manager);
221  }
222  for (SC_iterator it = _cache1.begin(); it != _cache1.end(); it++) {
223  (*it).second.setConfig(_conf_manager);
224  }
225  }
226 
227  if (iSetup.get<MuonGeometryRecord>().cacheIdentifier() != _configid) {
228  if (_debug)
229  std::cout << "DTTrig::updateES updating muon geometry" << std::endl;
230 
231  _geomid = iSetup.get<MuonGeometryRecord>().cacheIdentifier();
232  geomHandle = iSetup.getHandle(dtGeomToken_);
233  for (TU_iterator it = _cache.begin(); it != _cache.end(); it++) {
234  (*it).second.setGeom(geomHandle->chamber((*it).second.statId()));
235  }
236  }
237 }
238 
240  // Delete the map
241  _cache.clear();
242  _cache1.clear();
243 }
244 
245 DTSCTrigUnit* DTTrig::trigUnit(DTChamberId chid) { /*check();*/ return const_cast<DTSCTrigUnit*>(constTrigUnit(chid)); }
246 
248  // std::cout << " SC: running DTTrig::constTrigUnit(DTChamberId chid)" << std::endl;
249  TU_const_iterator it = _cache.find(chid);
250  if (it == _cache.end()) {
251  std::cout << "DTTrig::trigUnit: Trigger Unit not in the map: ";
252  std::cout << " wheel=" << chid.wheel();
253  std::cout << ", station=" << chid.station();
254  std::cout << ", sector=" << chid.sector();
255  std::cout << std::endl;
256  return nullptr;
257  }
258 
259  return &(*it).second;
260 }
261 
263  SC_const_iterator it = _cache1.find(scid);
264  if (it == _cache1.end()) {
265  std::cout << "DTTrig::SCUnit: Trigger Unit not in the map: ";
266  std::cout << " wheel=" << scid.wheel();
267  std::cout << ", sector=" << scid.sector();
268  std::cout << std::endl;
269  return nullptr;
270  }
271 
272  return &(*it).second;
273 }
274 
275 DTSCTrigUnit* DTTrig::trigUnit(int wheel, int stat, int sect) {
276  return const_cast<DTSCTrigUnit*>(constTrigUnit(wheel, stat, sect));
277 }
278 
279 DTSectColl const* DTTrig::SCUnit(int wheel, int sect) const {
280  sect++;
281  return SCUnit(DTSectCollId(wheel, sect));
282 }
283 
284 DTSCTrigUnit const* DTTrig::constTrigUnit(int wheel, int stat, int sect) const {
285  sect++; // offset 1 for sector number ([0,11] --> [1,12])
286  return constTrigUnit(DTChamberId(wheel, stat, sect));
287 }
288 
290  if (unit == nullptr)
291  return nullptr;
292  if (unit->nPhiSegm(step) < 1)
293  return nullptr;
294  return const_cast<DTChambPhSegm*>(unit->phiSegment(step, 1));
295 }
296 
298  if (unit == nullptr)
299  return nullptr;
300  if (unit->nPhiSegm(step) < 2)
301  return nullptr;
302  return const_cast<DTChambPhSegm*>(unit->phiSegment(step, 2));
303 }
304 
306  if (unit == nullptr)
307  return nullptr;
308  if (unit->nThetaSegm(step) < 1)
309  return nullptr;
310  return const_cast<DTChambThSegm*>(unit->thetaSegment(step, 1));
311 }
312 
314 
316 
318  if (sid.station() == 4)
319  return nullptr;
320  return chThetaSegm(trigUnit(sid), step);
321 }
322 
323 DTChambPhSegm* DTTrig::chPhiSegm1(int wheel, int stat, int sect, int step) {
324  return chPhiSegm1(trigUnit(wheel, stat, sect), step);
325  // to make it transparent to the outside world
326  // return chSectCollSegm1(wheel,stat,sect,step);
327 }
328 
329 DTChambPhSegm* DTTrig::chPhiSegm2(int wheel, int stat, int sect, int step) {
330  // if(stat==4&&(sect==3||sect==9)) {
331  // if hrizontal chambers of MB4 get first track of twin chamber (flag=1)
332  // return chPhiSegm1(trigUnit(wheel,stat,sect,1),step);
333  // } else {
334  return chPhiSegm2(trigUnit(wheel, stat, sect), step);
335  // to make it transparent to the outside world
336  // return chSectCollSegm2(wheel,stat,sect,step);
337  //}
338 }
339 
340 DTChambThSegm* DTTrig::chThetaSegm(int wheel, int stat, int sect, int step) {
341  if (stat == 4)
342  return nullptr;
343  return chThetaSegm(trigUnit(wheel, stat, sect), step);
344 }
345 
346 // SM sector collector section
348  if (unit == nullptr)
349  return nullptr;
350  if (unit->nSegmPh(step) < 1)
351  return nullptr;
352  return const_cast<DTSectCollPhSegm*>(unit->SectCollPhSegment(step, 1));
353 }
354 
356  if (unit == nullptr)
357  return nullptr;
358  if (unit->nSegmPh(step) < 2)
359  return nullptr;
360  return const_cast<DTSectCollPhSegm*>(unit->SectCollPhSegment(step, 2));
361 }
362 
364  return chSectCollPhSegm1(const_cast<DTSectColl*>(SCUnit(wheel, sect)), step);
365 }
366 
368  // if(stat==4&&(sect==3||sect==9)) {
369  // if hrizontal chambers of MB4 get first track of twin chamber (flag=1)
370  //return chSectCollSegm1(trigUnit(wheel,stat,sect,1),step);
371  //} else {
372  return chSectCollPhSegm2(const_cast<DTSectColl*>(SCUnit(wheel, sect)), step);
373  //}
374 }
375 
377  if (unit == nullptr)
378  return nullptr;
379  if (unit->nSegmTh(step) < 1)
380  return nullptr;
381  return const_cast<DTSectCollThSegm*>(unit->SectCollThSegment(step));
382 }
383 
385  return chSectCollThSegm(const_cast<DTSectColl*>(SCUnit(wheel, sect)), step);
386 }
387 
388 // end SM
389 
390 void DTTrig::dumpGeom() const {
391  /*check();*/
392  for (TU_const_iterator it = _cache.begin(); it != _cache.end(); it++) {
393  ((*it).second).dumpGeom();
394  }
395 }
396 
397 void DTTrig::dumpLuts(short int lut_btic, const DTConfigManager* conf) const {
398  for (TU_const_iterator it = _cache.begin(); it != _cache.end(); it++) {
399  const DTSCTrigUnit& thisTU = (*it).second;
400 
401  // dump lut command file from geometry
402  thisTU.dumpLUT(lut_btic);
403 
404  // dump lut command file from parameters (DB or CMSSW)
405  DTChamberId chid = thisTU.statId();
406  conf->dumpLUTParam(chid);
407  }
408 
409  return;
410 }
411 
412 std::vector<DTBtiTrigData> DTTrig::BtiTrigs() const {
413  /*check();*/
414  std::vector<DTBtiTrigData> trigs;
415  for (auto ptu = _cache.begin(); ptu != _cache.end(); ptu++) {
416  const DTSCTrigUnit& tu = (*ptu).second;
417  auto peb = tu.BtiTrigs()->end();
418  for (auto p = tu.BtiTrigs()->begin(); p != peb; p++) {
419  trigs.push_back(*p);
420  }
421  }
422  return trigs;
423 }
424 
425 std::vector<DTTracoTrigData> DTTrig::TracoTrigs() const {
426  std::vector<DTTracoTrigData> trigs;
427  /*check();*/
428  for (auto ptu = _cache.begin(); ptu != _cache.end(); ptu++) {
429  const DTSCTrigUnit& tu = (*ptu).second;
430  auto peb = tu.TracoTrigs()->end();
431  for (auto p = tu.TracoTrigs()->begin(); p != peb; p++) {
432  trigs.push_back(*p);
433  }
434  }
435  return trigs;
436 }
437 
438 std::vector<DTChambPhSegm> DTTrig::TSPhTrigs() const {
439  /*check();*/
440  std::vector<DTChambPhSegm> trigs;
441  for (auto ptu = _cache.begin(); ptu != _cache.end(); ptu++) {
442  const DTSCTrigUnit& tu = (*ptu).second;
443  auto peb = tu.TSPhTrigs()->end();
444  for (auto p = tu.TSPhTrigs()->begin(); p != peb; p++) {
445  trigs.push_back(*p);
446  }
447  }
448  return trigs;
449 }
450 
451 std::vector<DTChambThSegm> DTTrig::TSThTrigs() const {
452  /*check();*/
453  std::vector<DTChambThSegm> trigs;
454  for (auto ptu = _cache.begin(); ptu != _cache.end(); ptu++) {
455  const DTSCTrigUnit& tu = (*ptu).second;
456  auto peb = tu.TSThTrigs()->end();
457  for (auto p = tu.TSThTrigs()->begin(); p != peb; p++) {
458  trigs.push_back(*p);
459  }
460  }
461  return trigs;
462 }
463 
464 std::vector<DTSectCollPhSegm> DTTrig::SCPhTrigs() const {
465  /*check();*/
466  std::vector<DTSectCollPhSegm> trigs;
467  for (auto psc = _cache1.begin(); psc != _cache1.end(); psc++) {
468  // DTSCTrigUnit* tu = (*ptu).second;
469  //
470  // old SMDB:
471  // DTSectColl* tu = (*ptu).second;
472  // std::vector<DTChambPhSegm>::const_iterator p=0;
473  // std::vector<DTChambPhSegm>::const_iterator peb=tu->SCTrigs()->end();
474  // for(p=tu->SCTrigs()->begin();p!=peb;p++){
475  // trigs.push_back(*p);
476  // }
477 
478  const DTSectColl& sc = (*psc).second;
479  auto peb = sc.endPh();
480  for (auto p = sc.beginPh(); p != peb; p++) {
481  trigs.push_back(*p);
482  }
483  }
484  return trigs;
485 }
486 
487 std::vector<DTSectCollThSegm> DTTrig::SCThTrigs() const {
488  /*check();*/
489  std::vector<DTSectCollThSegm> trigs;
490  for (auto psc = _cache1.begin(); psc != _cache1.end(); psc++) {
491  const DTSectColl& sc = (*psc).second;
492  auto peb = sc.endTh();
493  for (auto p = sc.beginTh(); p != peb; p++) {
494  trigs.push_back(*p);
495  }
496  }
497  return trigs;
498 }
void reconstruct() override
Load BTIs triggers and run TSTheta algoritm.
Definition: DTTSTheta.h:81
virtual void reconstruct(const DTDigiCollection dtDigis)
Definition: DTBtiCard.h:101
std::vector< DTBtiTrigData > BtiTrigs() const
Return a copy of all the BTI triggers.
Definition: DTTrig.cc:412
int station() const
Return the station number.
Definition: DTChamberId.h:45
void clearCache()
Clear all traco stuff (cache & map)
Definition: DTTracoCard.cc:61
SCcontainer::iterator SC_iterator
Definition: DTTrig.h:65
void reconstruct() override
Load Trigger Units triggers and run Sector Collector algorithm.
Definition: DTSectColl.h:178
static const TGPicture * info(bool iBackgroundIsBlack)
void reconstruct() override
Load TRACO triggers and run TSPhi algorithm.
Definition: DTTSPhi.h:80
void dumpLUT(short int btic) const
Dump the Lut file.
Definition: DTSCTrigUnit.h:117
SCcontainer _cache1
Definition: DTTrig.h:242
std::map< DTChamberId, DTDigiCollection, std::less< DTChamberId > > DTDigiMap
Definition: DTTrig.h:69
DTChamberId id() const
Return the DTChamberId of this chamber.
Definition: DTChamber.cc:32
DTBtiCard * BtiTrigs() const
Return container of BTI triggers.
Definition: DTSCTrigUnit.h:84
bool _debug
Definition: DTTrig.h:248
void clearCache()
Clear all BTI stuff (map & cache)
Definition: DTBtiCard.cc:85
DTChambPhSegm * chPhiSegm2(DTChamberId sid, int step)
Return the second phi track segment in req. chamber/step.
Definition: DTTrig.cc:315
std::vector< DTChambThSegm > TSThTrigs() const
Return a copy of all the Trigger Server (Theta) triggers.
Definition: DTTrig.cc:451
DTTrig(const edm::ParameterSet &params, edm::ConsumesCollector &&ix)
Constructors.
Definition: DTTrig.cc:46
TUcontainer _cache
Definition: DTTrig.h:241
void reconstruct() override
Load BTIs triggers and run TRACOs algorithm.
Definition: DTTracoCard.h:101
int sizeTh() const
Return Theta cache size.
Definition: DTSectColl.h:160
TUcontainer::iterator TU_iterator
Definition: DTTrig.h:62
unsigned long long cacheIdentifier() const
DTSectCollPhSegm * chSectCollPhSegm1(DTSectColl *unit, int step)
Return the first phi track segment in req. chamber/step [SC step].
Definition: DTTrig.cc:347
DTChambPhSegm * chPhiSegm1(DTChamberId sid, int step)
Return the first phi track segment in req. chamber/step.
Definition: DTTrig.cc:313
unsigned long long _geomid
Definition: DTTrig.h:252
int wheel() const
Definition: DTSectCollId.h:30
int sizePh() const
Return Phi cache size.
Definition: DTSectColl.h:151
DTSCTrigUnit const * constTrigUnit(DTChamberId sid) const
Return a trigger unit - Muon numbering - const version.
Definition: DTTrig.cc:247
DTChamberId statId() const
Identifier of the associated chamber.
Definition: DTSCTrigUnit.h:72
int size() const
Get cache vector&#39;s size.
Definition: DTCache.h:44
int iEvent
Definition: GenABIO.cc:224
T const * product() const
Definition: ESHandle.h:86
std::vector< DTSectCollPhSegm >::const_iterator endPh() const
Return iterator to the end of Phi cache.
Definition: DTSectColl.h:154
std::vector< DTSectCollPhSegm > SCPhTrigs() const
Return a copy of all the Sector Collector (Phi) triggers.
Definition: DTTrig.cc:464
DTChamberId chamberId() const
Return the corresponding ChamberId.
std::vector< DTTracoTrigData > TracoTrigs() const
Return a copy of all the TRACO triggers.
Definition: DTTrig.cc:425
DTTSPhi * TSPhTrigs() const
Return the chamber Trigger Server (Phi)
Definition: DTSCTrigUnit.h:90
std::vector< DTSectCollThSegm > SCThTrigs() const
Return a copy of all the Sector Collector (Theta) triggers.
Definition: DTTrig.cc:487
void updateES(const edm::EventSetup &iSetup)
update the eventsetup info
Definition: DTTrig.cc:207
void dumpLUTParam(DTChamberId &chambid) const
Dump luts string commands from configuration parameters.
std::vector< DTSectCollThSegm >::const_iterator beginTh() const
Return iterator to the begni of Theta cache.
Definition: DTSectColl.h:157
edm::InputTag _digitag
Definition: DTTrig.h:244
std::vector< DTSectCollThSegm >::const_iterator endTh() const
Return iterator to the end of Theta cache.
Definition: DTSectColl.h:163
T get() const
Definition: EventSetup.h:79
const_iterator end() const
Get last cache element.
Definition: DTCache.h:41
void triggerReco(const edm::Event &iEvent, const edm::EventSetup &iSetup)
Run the whole trigger reconstruction chain.
Definition: DTTrig.cc:132
std::vector< DTChambPhSegm > TSPhTrigs() const
Return a copy of all the Trigger Server (Phi) triggers.
Definition: DTTrig.cc:438
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:130
Basic3DVector unit() const
DTSectCollPhSegm * chSectCollPhSegm2(DTSectColl *unit, int step)
Return the second phi track segment in req. chamber/step [SC step].
Definition: DTTrig.cc:355
DTDigiMap::iterator DTDigiMap_iterator
Definition: DTTrig.h:70
int sector() const
Definition: DTSectCollId.h:31
DTTSPhi * getTSPhi(int istat) const
Return TSPhi.
Definition: DTSectColl.h:69
TUcontainer::const_iterator TU_const_iterator
Definition: DTTrig.h:63
DTSCTrigUnit * trigUnit(DTChamberId sid)
Return a trigger unit - Muon numbering.
Definition: DTTrig.cc:245
edm::ESGetToken< DTConfigManager, DTConfigManagerRcd > confToken_
Definition: DTTrig.h:247
std::pair< const_iterator, const_iterator > Range
DTTracoCard * TracoTrigs() const
Return container of TRACO triggers.
Definition: DTSCTrigUnit.h:87
bool isValid() const
Definition: HandleBase.h:70
edm::ESGetToken< DTGeometry, MuonGeometryRecord > dtGeomBeginRunToken_
Definition: DTTrig.h:246
DTChambThSegm * chThetaSegm(DTChamberId sid, int step)
Return the theta candidates in req. chamber/step.
Definition: DTTrig.cc:317
void clearCache()
Local position in chamber of a trigger-data object.
Definition: DTSectColl.h:172
int wheel() const
Return the wheel number.
Definition: DTChamberId.h:42
int sector() const
Definition: DTChamberId.h:52
const DTConfigManager * _conf_manager
Definition: DTTrig.h:243
bool _inputexist
Definition: DTTrig.h:249
DTSectCollThSegm * chSectCollThSegm(DTSectColl *unit, int step)
Return the theta track segment in req. chamber/step [SC step].
Definition: DTTrig.cc:376
DTSuperLayerId superlayerId() const
Return the corresponding SuperLayerId.
Definition: DTLayerId.h:48
const std::vector< const DTChamber * > & chambers() const
Return a vector of all Chamber.
Definition: DTGeometry.cc:84
step
Definition: StallMonitor.cc:83
edm::ESGetToken< DTGeometry, MuonGeometryRecord > dtGeomToken_
Definition: DTTrig.h:245
DTTSTheta * getTSTheta(int istat) const
Return TSTheta.
Definition: DTSectColl.h:72
void dumpGeom() const
Dump the geometry.
Definition: DTTrig.cc:390
void dumpLuts(short int lut_btic, const DTConfigManager *conf) const
Dump the LUT files.
Definition: DTTrig.cc:397
DTSectColl const * SCUnit(DTSectCollId scid) const
Return a SC unit - Muon numbering - const version.
Definition: DTTrig.cc:262
void createTUs(const edm::EventSetup &iSetup)
Create the trigger units and store them in the cache.
Definition: DTTrig.cc:63
void clearCache()
Clear cache vector.
Definition: DTCache.h:47
const_iterator begin() const
Get first cache element.
Definition: DTCache.h:38
const DTChamber * chamber(const DTChamberId &id) const
Return a DTChamber given its id.
Definition: DTGeometry.cc:90
unsigned long long _configid
Definition: DTTrig.h:251
SCcontainer::const_iterator SC_const_iterator
Definition: DTTrig.h:66
DTTSTheta * TSThTrigs() const
Return the chamber Trigger Server (Theta)
Definition: DTSCTrigUnit.h:93
#define LogDebug(id)
void clear()
Clear the trigger units cache.
Definition: DTTrig.cc:239
std::vector< DTSectCollPhSegm >::const_iterator beginPh() const
Return iterator to the beghin of Phi cache.
Definition: DTSectColl.h:148