CMS 3D CMS Logo

L1TGlobalUtil.cc
Go to the documentation of this file.
1 // L1TGlobalUtil
2 //
3 // author: Brian Winer Ohio State
4 //
5 
7 
8 #include <fstream>
9 #include <iostream>
10 #include <memory>
11 
17 
20 
21 // constructor
22 l1t::L1TGlobalUtil::L1TGlobalUtil() : m_l1GtMenu(nullptr) {
23  // initialize cached IDs
24  m_l1GtMenuCacheID = 0ULL;
25  m_l1GtPfAlgoCacheID = 0ULL;
26  m_filledPrescales = false;
27 
28  edm::FileInPath f1("L1Trigger/L1TGlobal/data/Luminosity/startup/prescale_L1TGlobal.csv");
29  m_preScaleFileName = f1.fullPath();
30  m_numberPhysTriggers = 512; //need to get this out of the EventSetup
31  m_PreScaleColumn = 0;
33 
34  m_prescaleFactorsAlgoTrig = nullptr;
35  m_triggerMaskAlgoTrig = nullptr;
36 }
37 
40  UseEventSetupIn useEventSetupIn)
41  : L1TGlobalUtil(pset, iC, useEventSetupIn) {}
42 
45  UseEventSetupIn useEventSetupIn)
46  : L1TGlobalUtil() {
47  m_l1tGlobalUtilHelper = std::make_unique<L1TGlobalUtilHelper>(pset, iC);
48  m_readPrescalesFromFile = m_l1tGlobalUtilHelper->readPrescalesFromFile();
49  eventSetupConsumes(iC, useEventSetupIn);
50 }
51 
52 // destructor
54  // empty
55 }
56 
58 bool l1t::L1TGlobalUtil::valid() const { return m_l1GtMenuCacheID != 0ULL and m_l1GtMenu != nullptr; }
59 
61  edm::FileInPath f1("L1Trigger/L1TGlobal/data/Luminosity/startup/" + filename);
62  m_preScaleFileName = f1.fullPath();
63  m_PreScaleColumn = psColumn;
64 }
65 
67  // typically, the L1T menu and prescale table (may change only between Runs)
68  bool isRun = false;
69  retrieveL1Setup(evSetup, isRun);
70  // typically the prescale set index used and the event by event accept/reject info (changes between Events)
71  retrieveL1Event(iEvent, evSetup);
72 }
73 
75  const edm::EventSetup& evSetup,
76  edm::EDGetToken gtAlgToken) {
77  // typically, the L1T menu and prescale table (may change only between Runs)
78  bool isRun = false;
79  retrieveL1Setup(evSetup, isRun);
80  // typically the prescale set index used and the event by event accept/reject info (changes between Events)
81  retrieveL1Event(iEvent, evSetup, gtAlgToken);
82 }
83 
85  bool isRun = true;
86  retrieveL1Setup(evSetup, isRun);
87 }
88 
89 void l1t::L1TGlobalUtil::retrieveL1Setup(const edm::EventSetup& evSetup, bool isRun) {
90  // get / update the trigger menu from the EventSetup
91  // local cache & check on cacheIdentifier
92  auto menuRcd = evSetup.get<L1TUtmTriggerMenuRcd>();
93  unsigned long long l1GtMenuCacheID = menuRcd.cacheIdentifier();
94 
95  if (m_l1GtMenuCacheID != l1GtMenuCacheID) {
96  if (isRun) {
97  m_l1GtMenu = &menuRcd.get(m_L1TUtmTriggerMenuRunToken);
98  } else {
99  m_l1GtMenu = &menuRcd.get(m_L1TUtmTriggerMenuEventToken);
100  }
101 
102  //std::cout << "Attempting to fill the map " << std::endl;
103  m_algorithmMap = &(m_l1GtMenu->getAlgorithmMap());
104 
105  //reset vectors since we have new menu
106  resetDecisionVectors();
107 
108  m_l1GtMenuCacheID = l1GtMenuCacheID;
109  }
110 
111  if (!m_readPrescalesFromFile) {
112  auto vetosRcd = evSetup.get<L1TGlobalPrescalesVetosFractRcd>();
113  unsigned long long l1GtPfAlgoCacheID = vetosRcd.cacheIdentifier();
114 
115  if (m_l1GtPfAlgoCacheID != l1GtPfAlgoCacheID) {
116  //std::cout << "Reading Prescales and Masks from dB" << std::endl;
117 
118  // clear and dimension
119  resetPrescaleVectors();
120  resetMaskVectors();
121  m_PreScaleColumn = 0;
122  m_numberOfPreScaleColumns = 0;
123  m_numberPhysTriggers = 0;
124 
125  const L1TGlobalPrescalesVetosFract* es = nullptr;
126  if (isRun) {
127  es = &vetosRcd.get(m_L1TGlobalPrescalesVetosFractRunToken);
128  } else {
129  es = &vetosRcd.get(m_L1TGlobalPrescalesVetosFractEventToken);
130  }
131  m_l1GtPrescalesVetoes = PrescalesVetosFractHelper::readFromEventSetup(es);
132 
133  m_prescaleFactorsAlgoTrig = &(m_l1GtPrescalesVetoes->prescaleTable());
134  m_numberOfPreScaleColumns = m_prescaleFactorsAlgoTrig->size();
135  m_numberPhysTriggers =
136  (*m_prescaleFactorsAlgoTrig)[0].size(); // assumes all prescale columns are the same length
137 
138  m_triggerMaskAlgoTrig = &(m_l1GtPrescalesVetoes->triggerAlgoBxMask());
139 
140  m_l1GtPfAlgoCacheID = l1GtPfAlgoCacheID;
141  }
142  } else {
143  //Load the prescales from external file
144  // clear and dimension
145  if (!m_filledPrescales) {
146  resetPrescaleVectors();
147  resetMaskVectors();
148 
149  loadPrescalesAndMasks();
150 
151  // Set Prescale factors to initial
152  m_prescaleFactorsAlgoTrig = &m_initialPrescaleFactorsAlgoTrig;
153  m_triggerMaskAlgoTrig = &m_initialTriggerMaskAlgoTrig;
154  m_filledPrescales = true;
155  }
156  }
157 
158  //Protect against poor prescale column choice (I don't think there is a way this happen as currently structured)
159  if (m_PreScaleColumn > m_prescaleFactorsAlgoTrig->size()) {
160  LogTrace("l1t|Global") << "\nNo Prescale Set: " << m_PreScaleColumn
161  << "\nMax Prescale Set value : " << m_prescaleFactorsAlgoTrig->size()
162  << "\nSetting prescale column to 0" << std::endl;
163  m_PreScaleColumn = 0;
164  }
165  //std::cout << "Using prescale column: " << m_PreScaleColumn << std::endl;
166  const std::vector<double>& prescaleSet = (*m_prescaleFactorsAlgoTrig)[m_PreScaleColumn];
167 
168  for (std::map<std::string, L1TUtmAlgorithm>::const_iterator itAlgo = m_algorithmMap->begin();
169  itAlgo != m_algorithmMap->end();
170  itAlgo++) {
171  // Get the algorithm name
172  std::string algName = itAlgo->first;
173  int algBit = (itAlgo->second).getIndex(); //algoBitNumber();
174 
175  (m_prescales[algBit]).first = algName;
176  if (size_t(algBit) < prescaleSet.size()) {
177  (m_prescales[algBit]).second = prescaleSet[algBit];
178  }
179  LogDebug("l1t|Global") << "Number of bunch crossings stored: " << (*m_triggerMaskAlgoTrig).size() << endl;
180 
181  const std::map<int, std::vector<int> >* triggerAlgoMaskAlgoTrig = m_triggerMaskAlgoTrig;
182  std::map<int, std::vector<int> >::const_iterator it = triggerAlgoMaskAlgoTrig->begin();
183 
184  std::vector<int> maskedBxs;
185  (m_masks[algBit]).first = algName;
186  (m_masks[algBit]).second = maskedBxs;
187  while (it != triggerAlgoMaskAlgoTrig->end()) {
188  std::vector<int> masks = it->second;
189  //std::cout<< "BX: " << it->first<<" VecSize: "<< masks.size();
190  //std::cout << "\tMasked algos: ";
191  for (unsigned int imask = 0; imask < masks.size(); imask++) {
192  if (masks.at(imask) == algBit)
193  maskedBxs.push_back(it->first);
194  // std::cout << "\t" << masks.at(imask);
195  }
196  //std::cout << "\n";
197  it++;
198  }
199 
200  if (!maskedBxs.empty()) {
201  LogDebug("l1t|Global") << "i Algo: " << algBit << "\t" << algName << " masked\n";
202  for (unsigned int ibx = 0; ibx < maskedBxs.size(); ibx++) {
203  // std::cout << "\t" << maskedBxs.at(ibx);
204  (m_masks[algBit]).second = maskedBxs;
205  }
206  }
207  }
208 }
209 
211  retrieveL1Event(iEvent, evSetup, m_l1tGlobalUtilHelper->l1tAlgBlkToken());
212 }
213 
215  const edm::EventSetup& evSetup,
216  edm::EDGetToken gtAlgToken) {
217  // Get the Global Trigger Output Algorithm block
218  iEvent.getByToken(gtAlgToken, m_uGtAlgBlk);
219  m_finalOR = false;
220 
221  //Make sure we have a valid AlgBlk
222  if (m_uGtAlgBlk.isValid()) {
223  // get the GlabalAlgBlk (Stupid find better way) of BX=0
224  std::vector<GlobalAlgBlk>::const_iterator algBlk = m_uGtAlgBlk->begin(0);
225  if (algBlk != m_uGtAlgBlk->end(0)) {
226  if (!m_readPrescalesFromFile) {
227  m_PreScaleColumn = static_cast<unsigned int>(algBlk->getPreScColumn());
228 
229  // Fix for MC prescale column being set to index+1 in early versions of uGT emulator
230  if (iEvent.run() == 1) {
231  if (m_prescaleFactorsAlgoTrig->size() == 1 && m_PreScaleColumn == 1)
232  m_PreScaleColumn = 0;
233  }
234 
235  // add protection against out-of-bound index for prescale column
236  if (m_PreScaleColumn >= m_prescaleFactorsAlgoTrig->size()) {
237  LogDebug("l1t|Global") << "Prescale column extracted from GlobalAlgBlk too large: " << m_PreScaleColumn
238  << "\tMaximum value allowed: " << m_prescaleFactorsAlgoTrig->size() - 1
239  << "\tResetting prescale column to 0" << std::endl;
240  m_PreScaleColumn = 0;
241  }
242  }
243  const std::vector<double>& prescaleSet = (*m_prescaleFactorsAlgoTrig)[m_PreScaleColumn];
244 
245  // Grab the final OR from the AlgBlk,
246  m_finalOR = algBlk->getFinalOR();
247 
248  // Make a map of the trigger name and whether it passed various stages (initial,prescale,final)
249  // Note: might be able to improve performance by not full remaking map with names each time
250  for (std::map<std::string, L1TUtmAlgorithm>::const_iterator itAlgo = m_algorithmMap->begin();
251  itAlgo != m_algorithmMap->end();
252  itAlgo++) {
253  // Get the algorithm name
254  std::string algName = itAlgo->first;
255  int algBit = (itAlgo->second).getIndex(); //algoBitNumber();
256 
257  bool decisionInitial = algBlk->getAlgoDecisionInitial(algBit);
258  (m_decisionsInitial[algBit]).first = algName;
259  (m_decisionsInitial[algBit]).second = decisionInitial;
260 
261  bool decisionInterm = algBlk->getAlgoDecisionInterm(algBit);
262  (m_decisionsInterm[algBit]).first = algName;
263  (m_decisionsInterm[algBit]).second = decisionInterm;
264 
265  bool decisionFinal = algBlk->getAlgoDecisionFinal(algBit);
266  (m_decisionsFinal[algBit]).first = algName;
267  (m_decisionsFinal[algBit]).second = decisionFinal;
268 
269  (m_prescales[algBit]).first = algName;
270  if (size_t(algBit) < prescaleSet.size()) {
271  (m_prescales[algBit]).second = prescaleSet[algBit];
272  }
273  }
274  } else {
275  //cout << "Error empty AlgBlk recovered.\n";
276  }
277  } else {
278  //cout<< "Error no valid uGT Algorithm Data with Token provided " << endl;
279  }
280 }
281 
283  std::ifstream inputPrescaleFile;
284  //std::cout << "Reading prescales from file: " << m_preScaleFileName << std::endl;
285  inputPrescaleFile.open(m_preScaleFileName);
286 
287  std::vector<std::vector<int> > vec;
288  std::vector<std::vector<double> > prescale_vec;
289 
290  if (inputPrescaleFile) {
291  std::string prefix1("#");
292  std::string prefix2("-1");
293 
295 
296  bool first = true;
297 
298  while (getline(inputPrescaleFile, line)) {
299  if (!line.compare(0, prefix1.size(), prefix1))
300  continue;
301  //if( !line.compare(0, prefix2.size(), prefix2) ) continue;
302 
303  istringstream split(line);
304  int value;
305  int col = 0;
306  char sep;
307 
308  while (split >> value) {
309  if (first) {
310  // Each new value read on line 1 should create a new inner vector
311  vec.push_back(std::vector<int>());
312  }
313 
314  vec[col].push_back(value);
315  ++col;
316 
317  // read past the separator
318  split >> sep;
319  }
320 
321  // Finished reading line 1 and creating as many inner
322  // vectors as required
323  first = false;
324  }
325 
326  int NumPrescaleSets = 0;
327  for (int iCol = 0; iCol < int(vec.size()); iCol++) {
328  if (!vec[iCol].empty()) {
329  int firstRow = vec[iCol][0];
330 
331  if (firstRow >= 0)
332  NumPrescaleSets++;
333  //else if( firstRow==-2 ) maskColumn = iCol;
334  //else if( firstRow==-3 ) maskVetoColumn = iCol;
335  }
336  }
337 
338  //std::cout << "NumPrescaleSets= " << NumPrescaleSets << std::endl;
339  if (NumPrescaleSets > 0) {
340  // Fill default prescale set
341  for (int iSet = 0; iSet < NumPrescaleSets; iSet++) {
342  prescale_vec.push_back(std::vector<double>());
343  for (unsigned int iBit = 0; iBit < m_numberPhysTriggers; ++iBit) {
344  int inputDefaultPrescale = 1;
345  prescale_vec[iSet].push_back(inputDefaultPrescale);
346  }
347  }
348 
349  // Fill non-trivial prescale set
350  for (int iBit = 1; iBit < int(vec[0].size()); iBit++) {
351  unsigned int algoBit = vec[0][iBit];
352  // algoBit must be less than the number of triggers
353  if (algoBit < m_numberPhysTriggers) {
354  for (int iSet = 0; iSet < int(vec.size()); iSet++) {
355  int useSet = -1;
356  if (!vec[iSet].empty()) {
357  useSet = vec[iSet][0];
358  }
359  useSet -= 1;
360 
361  if (useSet < 0)
362  continue;
363 
364  int prescale = vec[iSet][iBit];
365  prescale_vec[useSet][algoBit] = prescale;
366  }
367  } else {
368  LogTrace("l1t|Global") << "\nPrescale file has algo bit: " << algoBit
369  << "\nThis is larger than the number of triggers: " << m_numberPhysTriggers
370  << "\nSomething is wrong. Ignoring." << std::endl;
371  }
372  }
373  }
374 
375  } else {
376  LogTrace("l1t|Global") << "\nCould not find file: " << m_preScaleFileName
377  << "\nFilling the prescale vectors with prescale 1"
378  << "\nSetting prescale set to 0" << std::endl;
379 
380  m_PreScaleColumn = 0;
381 
382  for (int col = 0; col < 1; col++) {
383  prescale_vec.push_back(std::vector<double>());
384  for (unsigned int iBit = 0; iBit < m_numberPhysTriggers; ++iBit) {
385  int inputDefaultPrescale = 0;
386  prescale_vec[col].push_back(inputDefaultPrescale);
387  }
388  }
389  }
390 
391  inputPrescaleFile.close();
392 
393  m_initialPrescaleFactorsAlgoTrig = prescale_vec;
394  // setting of bx masks from an input file not enabled; do not see a use case at the moment
395  std::map<int, std::vector<int> > m_initialTriggerMaskAlgoTrig;
396 }
397 
399  if (useEventSetupIn == UseEventSetupIn::Run || useEventSetupIn == UseEventSetupIn::RunAndEvent) {
400  m_L1TUtmTriggerMenuRunToken = iC.esConsumes<L1TUtmTriggerMenu, L1TUtmTriggerMenuRcd, edm::Transition::BeginRun>();
401  if (!m_readPrescalesFromFile) {
402  m_L1TGlobalPrescalesVetosFractRunToken =
404  }
405  }
406  if (useEventSetupIn == UseEventSetupIn::Event || useEventSetupIn == UseEventSetupIn::RunAndEvent) {
407  m_L1TUtmTriggerMenuEventToken = iC.esConsumes<L1TUtmTriggerMenu, L1TUtmTriggerMenuRcd>();
408  if (!m_readPrescalesFromFile) {
409  m_L1TGlobalPrescalesVetosFractEventToken =
411  }
412  }
413 }
414 
416  // Reset all the vector contents with null information
417  m_decisionsInitial.clear();
418  m_decisionsInitial.resize(m_maxNumberPhysTriggers);
419  m_decisionsInterm.clear();
420  m_decisionsInterm.resize(m_maxNumberPhysTriggers);
421  m_decisionsFinal.clear();
422  m_decisionsFinal.resize(m_maxNumberPhysTriggers);
423 
424  for (unsigned int algBit = 0; algBit < m_maxNumberPhysTriggers; algBit++) {
425  (m_decisionsInitial[algBit]).first = "NULL";
426  (m_decisionsInitial[algBit]).second = false;
427 
428  (m_decisionsInterm[algBit]).first = "NULL";
429  (m_decisionsInterm[algBit]).second = false;
430 
431  (m_decisionsFinal[algBit]).first = "NULL";
432  (m_decisionsFinal[algBit]).second = false;
433  }
434 }
435 
437  // Reset all the vector contents with null information
438  m_prescales.clear();
439  m_prescales.resize(m_maxNumberPhysTriggers);
440 
441  for (unsigned int algBit = 0; algBit < m_maxNumberPhysTriggers; algBit++) {
442  (m_prescales[algBit]).first = "NULL";
443  (m_prescales[algBit]).second = 1;
444  }
445 }
446 
448  // Reset all the vector contents with null information
449  m_masks.clear();
450  m_masks.resize(m_maxNumberPhysTriggers);
451 
452  for (unsigned int algBit = 0; algBit < m_maxNumberPhysTriggers; algBit++) {
453  (m_masks[algBit]).first = "NULL";
454  // ccla (m_masks[algBit]).second = true;
455  }
456 }
457 
458 const bool l1t::L1TGlobalUtil::getAlgBitFromName(const std::string& algName, int& bit) const {
459  std::map<std::string, L1TUtmAlgorithm>::const_iterator itAlgo = m_algorithmMap->find(algName);
460  if (itAlgo != m_algorithmMap->end()) {
461  bit = (itAlgo->second).getIndex(); //algoBitNumber();
462  return true;
463  }
464 
465  return false; //did not find anything by that name
466 }
467 
468 const bool l1t::L1TGlobalUtil::getAlgNameFromBit(int& bit, std::string& algName) const {
469  // since we are just looking up the name, doesn't matter which vector we get it from
470  if ((m_decisionsInitial[bit]).first != "NULL") {
471  algName = (m_decisionsInitial[bit]).first;
472  return true;
473  }
474  return false; //No name associated with this bit
475 }
476 
477 const bool l1t::L1TGlobalUtil::getInitialDecisionByBit(int& bit, bool& decision) const {
478  /*
479  for(std::vector<GlobalAlgBlk>::const_iterator algBlk = m_uGtAlgBlk->begin(0); algBlk != m_uGtAlgBlk->end(0); ++algBlk) {
480  decision = algBlk->getAlgoDecisionFinal(bit);
481  }
482  */
483  // Need some check that this is a valid bit
484  if ((m_decisionsInitial[bit]).first != "NULL") {
485  decision = (m_decisionsInitial[bit]).second;
486  return true;
487  }
488 
489  return false; //couldn't get the information requested.
490 }
491 const bool l1t::L1TGlobalUtil::getIntermDecisionByBit(int& bit, bool& decision) const {
492  // Need some check that this is a valid bit
493  if ((m_decisionsInterm[bit]).first != "NULL") {
494  decision = (m_decisionsInterm[bit]).second;
495  return true;
496  }
497 
498  return false; //couldn't get the information requested.
499 }
500 const bool l1t::L1TGlobalUtil::getFinalDecisionByBit(int& bit, bool& decision) const {
501  // Need some check that this is a valid bit
502  if ((m_decisionsFinal[bit]).first != "NULL") {
503  decision = (m_decisionsFinal[bit]).second;
504  return true;
505  }
506 
507  return false; //couldn't get the information requested.
508 }
509 const bool l1t::L1TGlobalUtil::getPrescaleByBit(int& bit, double& prescale) const {
510  // Need some check that this is a valid bit
511  if ((m_prescales[bit]).first != "NULL") {
512  prescale = (m_prescales[bit]).second;
513  return true;
514  }
515 
516  return false; //couldn't get the information requested.
517 }
518 const bool l1t::L1TGlobalUtil::getMaskByBit(int& bit, std::vector<int>& mask) const {
519  // Need some check that this is a valid bit
520  if ((m_masks[bit]).first != "NULL") {
521  mask = (m_masks[bit]).second;
522  return true;
523  }
524 
525  return false; //couldn't get the information requested.
526 }
527 
528 const bool l1t::L1TGlobalUtil::getInitialDecisionByName(const std::string& algName, bool& decision) const {
529  int bit = -1;
530  if (getAlgBitFromName(algName, bit)) {
531  decision = (m_decisionsInitial[bit]).second;
532  return true;
533  }
534 
535  return false; //trigger name was not the menu.
536 }
537 
538 const bool l1t::L1TGlobalUtil::getIntermDecisionByName(const std::string& algName, bool& decision) const {
539  int bit = -1;
540  if (getAlgBitFromName(algName, bit)) {
541  decision = (m_decisionsInterm[bit]).second;
542  return true;
543  }
544 
545  return false; //trigger name was not the menu.
546 }
547 
548 const bool l1t::L1TGlobalUtil::getFinalDecisionByName(const std::string& algName, bool& decision) const {
549  int bit = -1;
550  if (getAlgBitFromName(algName, bit)) {
551  decision = (m_decisionsFinal[bit]).second;
552  return true;
553  }
554 
555  return false; //trigger name was not the menu.
556 }
557 const bool l1t::L1TGlobalUtil::getPrescaleByName(const std::string& algName, double& prescale) const {
558  int bit = -1;
559  if (getAlgBitFromName(algName, bit)) {
560  prescale = (m_prescales[bit]).second;
561  return true;
562  }
563 
564  return false; //trigger name was not the menu.
565 }
566 const bool l1t::L1TGlobalUtil::getMaskByName(const std::string& algName, std::vector<int>& mask) const {
567  int bit = -1;
568  if (getAlgBitFromName(algName, bit)) {
569  mask = (m_masks[bit]).second;
570  return true;
571  }
572 
573  return false; //trigger name was not the menu.
574 }
size
Write out results.
std::string m_preScaleFileName
const bool getInitialDecisionByBit(int &bit, bool &decision) const
const bool getMaskByName(const std::string &algName, std::vector< int > &mask) const
void retrieveL1Event(const edm::Event &iEvent, const edm::EventSetup &evSetup)
const bool getAlgBitFromName(const std::string &AlgName, int &bit) const
unsigned long long m_l1GtPfAlgoCacheID
void OverridePrescalesAndMasks(std::string filename, unsigned int psColumn=1.)
static const PrescalesVetosFractHelper * readFromEventSetup(const L1TGlobalPrescalesVetosFract *es)
unsigned int m_PreScaleColumn
std::unique_ptr< L1TGlobalUtilHelper > m_l1tGlobalUtilHelper
void resetDecisionVectors()
clear decision vectors on a menu change
unsigned long long cacheIdentifier() const
#define LogTrace(id)
const bool getIntermDecisionByBit(int &bit, bool &decision) const
U second(std::pair< T, U > const &p)
int iEvent
Definition: GenABIO.cc:224
const std::map< int, std::vector< int > > * m_triggerMaskAlgoTrig
void eventSetupConsumes(edm::ConsumesCollector &iC, UseEventSetupIn useEventSetupIn)
const bool getFinalDecisionByBit(int &bit, bool &decision) const
const bool getInitialDecisionByName(const std::string &algName, bool &decision) const
T get() const
Definition: EventSetup.h:82
Definition: value.py:1
const std::vector< std::vector< double > > * m_prescaleFactorsAlgoTrig
const bool getPrescaleByName(const std::string &algName, double &prescale) const
unsigned long long m_l1GtMenuCacheID
bool valid() const
check that the L1TGlobalUtil has been properly initialised
const bool getIntermDecisionByName(const std::string &algName, bool &decision) const
const bool getPrescaleByBit(int &bit, double &prescale) const
const bool getAlgNameFromBit(int &bit, std::string &AlgName) const
col
Definition: cuy.py:1009
const bool getFinalDecisionByName(const std::string &algName, bool &decision) const
const bool getMaskByBit(int &bit, std::vector< int > &mask) const
UseEventSetupIn
Definition: L1TGlobalUtil.h:42
virtual ~L1TGlobalUtil()
destructor
void retrieveL1Setup(const edm::EventSetup &evSetup)
unsigned int m_numberPhysTriggers
void retrieveL1(const edm::Event &iEvent, const edm::EventSetup &evSetup)
initialize the class (mainly reserve)
#define LogDebug(id)