Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-11-15 23:40:39

0001 // L1TGlobalUtil
0002 //
0003 // author: Brian Winer Ohio State
0004 //
0005 
0006 #include "L1Trigger/L1TGlobal/interface/L1TGlobalUtil.h"
0007 
0008 #include <fstream>
0009 #include <iostream>
0010 #include <memory>
0011 
0012 #include "FWCore/Framework/interface/Event.h"
0013 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0014 #include "FWCore/Utilities/interface/EDGetToken.h"
0015 #include "FWCore/Utilities/interface/InputTag.h"
0016 #include "FWCore/Utilities/interface/Transition.h"
0017 
0018 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0019 #include "FWCore/MessageLogger/interface/MessageDrop.h"
0020 
0021 // constructor
0022 l1t::L1TGlobalUtil::L1TGlobalUtil() : m_l1GtMenu(nullptr) {
0023   // initialize cached IDs
0024   m_l1GtMenuCacheID = 0ULL;
0025   m_l1GtPfAlgoCacheID = 0ULL;
0026   m_filledPrescales = false;
0027 
0028   edm::FileInPath f1("L1Trigger/L1TGlobal/data/Luminosity/startup/prescale_L1TGlobal.csv");
0029   m_preScaleFileName = f1.fullPath();
0030   m_numberPhysTriggers = 512;  //need to get this out of the EventSetup
0031   m_PreScaleColumn = 0;
0032   m_readPrescalesFromFile = false;
0033 
0034   m_prescaleFactorsAlgoTrig = nullptr;
0035   m_triggerMaskAlgoTrig = nullptr;
0036 }
0037 
0038 l1t::L1TGlobalUtil::L1TGlobalUtil(edm::ParameterSet const& pset,
0039                                   edm::ConsumesCollector&& iC,
0040                                   UseEventSetupIn useEventSetupIn)
0041     : L1TGlobalUtil(pset, iC, useEventSetupIn) {}
0042 
0043 l1t::L1TGlobalUtil::L1TGlobalUtil(edm::ParameterSet const& pset,
0044                                   edm::ConsumesCollector& iC,
0045                                   UseEventSetupIn useEventSetupIn)
0046     : L1TGlobalUtil() {
0047   m_l1tGlobalUtilHelper = std::make_unique<L1TGlobalUtilHelper>(pset, iC);
0048   m_readPrescalesFromFile = m_l1tGlobalUtilHelper->readPrescalesFromFile();
0049   eventSetupConsumes(iC, useEventSetupIn);
0050 }
0051 
0052 // destructor
0053 l1t::L1TGlobalUtil::~L1TGlobalUtil() {
0054   // empty
0055 }
0056 
0057 /// check that the L1TGlobalUtil has been properly initialised
0058 bool l1t::L1TGlobalUtil::valid() const { return m_l1GtMenuCacheID != 0ULL and m_l1GtMenu != nullptr; }
0059 
0060 void l1t::L1TGlobalUtil::OverridePrescalesAndMasks(std::string filename, unsigned int psColumn) {
0061   edm::FileInPath f1("L1Trigger/L1TGlobal/data/Luminosity/startup/" + filename);
0062   m_preScaleFileName = f1.fullPath();
0063   m_PreScaleColumn = psColumn;
0064 }
0065 
0066 void l1t::L1TGlobalUtil::retrieveL1(const edm::Event& iEvent, const edm::EventSetup& evSetup) {
0067   // typically, the L1T menu and prescale table (may change only between Runs)
0068   bool isRun = false;
0069   retrieveL1Setup(evSetup, isRun);
0070   // typically the prescale set index used and the event by event accept/reject info (changes between Events)
0071   retrieveL1Event(iEvent, evSetup);
0072 }
0073 
0074 void l1t::L1TGlobalUtil::retrieveL1(const edm::Event& iEvent,
0075                                     const edm::EventSetup& evSetup,
0076                                     edm::EDGetToken gtAlgToken) {
0077   // typically, the L1T menu and prescale table (may change only between Runs)
0078   bool isRun = false;
0079   retrieveL1Setup(evSetup, isRun);
0080   // typically the prescale set index used and the event by event accept/reject info (changes between Events)
0081   retrieveL1Event(iEvent, evSetup, gtAlgToken);
0082 }
0083 
0084 void l1t::L1TGlobalUtil::retrieveL1Setup(const edm::EventSetup& evSetup) {
0085   bool isRun = true;
0086   retrieveL1Setup(evSetup, isRun);
0087 }
0088 
0089 void l1t::L1TGlobalUtil::retrieveL1Setup(const edm::EventSetup& evSetup, bool isRun) {
0090   // get / update the trigger menu from the EventSetup
0091   // local cache & check on cacheIdentifier
0092   auto menuRcd = evSetup.get<L1TUtmTriggerMenuRcd>();
0093   unsigned long long l1GtMenuCacheID = menuRcd.cacheIdentifier();
0094 
0095   if (m_l1GtMenuCacheID != l1GtMenuCacheID) {
0096     if (isRun) {
0097       m_l1GtMenu = &menuRcd.get(m_L1TUtmTriggerMenuRunToken);
0098     } else {
0099       m_l1GtMenu = &menuRcd.get(m_L1TUtmTriggerMenuEventToken);
0100     }
0101 
0102     //std::cout << "Attempting to fill the map " << std::endl;
0103     m_algorithmMap = &(m_l1GtMenu->getAlgorithmMap());
0104 
0105     //reset vectors since we have new menu
0106     resetDecisionVectors();
0107 
0108     m_l1GtMenuCacheID = l1GtMenuCacheID;
0109   }
0110 
0111   if (!m_readPrescalesFromFile) {
0112     auto vetosRcd = evSetup.get<L1TGlobalPrescalesVetosFractRcd>();
0113     unsigned long long l1GtPfAlgoCacheID = vetosRcd.cacheIdentifier();
0114 
0115     if (m_l1GtPfAlgoCacheID != l1GtPfAlgoCacheID) {
0116       //std::cout << "Reading Prescales and Masks from dB" << std::endl;
0117 
0118       // clear and dimension
0119       resetPrescaleVectors();
0120       resetMaskVectors();
0121       m_PreScaleColumn = 0;
0122       m_numberOfPreScaleColumns = 0;
0123       m_numberPhysTriggers = 0;
0124 
0125       const L1TGlobalPrescalesVetosFract* es = nullptr;
0126       if (isRun) {
0127         es = &vetosRcd.get(m_L1TGlobalPrescalesVetosFractRunToken);
0128       } else {
0129         es = &vetosRcd.get(m_L1TGlobalPrescalesVetosFractEventToken);
0130       }
0131       m_l1GtPrescalesVetoes = PrescalesVetosFractHelper::readFromEventSetup(es);
0132 
0133       m_prescaleFactorsAlgoTrig = &(m_l1GtPrescalesVetoes->prescaleTable());
0134       m_numberOfPreScaleColumns = m_prescaleFactorsAlgoTrig->size();
0135       m_numberPhysTriggers =
0136           (*m_prescaleFactorsAlgoTrig)[0].size();  // assumes all prescale columns are the same length
0137 
0138       m_triggerMaskAlgoTrig = &(m_l1GtPrescalesVetoes->triggerAlgoBxMask());
0139 
0140       m_l1GtPfAlgoCacheID = l1GtPfAlgoCacheID;
0141     }
0142   } else {
0143     //Load the prescales from external file
0144     // clear and dimension
0145     if (!m_filledPrescales) {
0146       resetPrescaleVectors();
0147       resetMaskVectors();
0148 
0149       loadPrescalesAndMasks();
0150 
0151       // Set Prescale factors to initial
0152       m_prescaleFactorsAlgoTrig = &m_initialPrescaleFactorsAlgoTrig;
0153       m_triggerMaskAlgoTrig = &m_initialTriggerMaskAlgoTrig;
0154       m_filledPrescales = true;
0155     }
0156   }
0157 
0158   //Protect against poor prescale column choice (I don't think there is a way this happen as currently structured)
0159   if (m_PreScaleColumn > m_prescaleFactorsAlgoTrig->size()) {
0160     LogTrace("l1t|Global") << "\nNo Prescale Set: " << m_PreScaleColumn
0161                            << "\nMax Prescale Set value : " << m_prescaleFactorsAlgoTrig->size()
0162                            << "\nSetting prescale column to 0" << std::endl;
0163     m_PreScaleColumn = 0;
0164   }
0165   //std::cout << "Using prescale column: " << m_PreScaleColumn << std::endl;
0166   const std::vector<double>& prescaleSet = (*m_prescaleFactorsAlgoTrig)[m_PreScaleColumn];
0167 
0168   for (std::map<std::string, L1TUtmAlgorithm>::const_iterator itAlgo = m_algorithmMap->begin();
0169        itAlgo != m_algorithmMap->end();
0170        itAlgo++) {
0171     // Get the algorithm name
0172     std::string_view algName = itAlgo->first;
0173     int algBit = (itAlgo->second).getIndex();  //algoBitNumber();
0174 
0175     (m_prescales[algBit]).first = algName;
0176     if (size_t(algBit) < prescaleSet.size()) {
0177       (m_prescales[algBit]).second = prescaleSet[algBit];
0178     }
0179     LogDebug("l1t|Global") << "Number of bunch crossings stored: " << (*m_triggerMaskAlgoTrig).size() << endl;
0180 
0181     const std::map<int, std::vector<int> >* triggerAlgoMaskAlgoTrig = m_triggerMaskAlgoTrig;
0182     std::map<int, std::vector<int> >::const_iterator it = triggerAlgoMaskAlgoTrig->begin();
0183 
0184     std::vector<int> maskedBxs;
0185     (m_masks[algBit]).first = algName;
0186     (m_masks[algBit]).second = maskedBxs;
0187     while (it != triggerAlgoMaskAlgoTrig->end()) {
0188       std::vector<int> masks = it->second;
0189       //std::cout<< "BX: " << it->first<<" VecSize: "<< masks.size();
0190       //std::cout << "\tMasked algos: ";
0191       for (unsigned int imask = 0; imask < masks.size(); imask++) {
0192         if (masks.at(imask) == algBit)
0193           maskedBxs.push_back(it->first);
0194         // std::cout << "\t" << masks.at(imask);
0195       }
0196       //std::cout << "\n";
0197       it++;
0198     }
0199 
0200     if (!maskedBxs.empty()) {
0201       LogDebug("l1t|Global") << "i Algo: " << algBit << "\t" << algName << " masked\n";
0202       for (unsigned int ibx = 0; ibx < maskedBxs.size(); ibx++) {
0203         // std::cout << "\t" << maskedBxs.at(ibx);
0204         (m_masks[algBit]).second = maskedBxs;
0205       }
0206     }
0207   }
0208 }
0209 
0210 void l1t::L1TGlobalUtil::retrieveL1Event(const edm::Event& iEvent, const edm::EventSetup& evSetup) {
0211   retrieveL1Event(iEvent, evSetup, m_l1tGlobalUtilHelper->l1tAlgBlkToken());
0212 }
0213 
0214 void l1t::L1TGlobalUtil::retrieveL1Event(const edm::Event& iEvent,
0215                                          const edm::EventSetup& evSetup,
0216                                          edm::EDGetToken gtAlgToken) {
0217   // Get the Global Trigger Output Algorithm block
0218   iEvent.getByToken(gtAlgToken, m_uGtAlgBlk);
0219   m_finalOR = false;
0220 
0221   //Make sure we have a valid AlgBlk
0222   if (m_uGtAlgBlk.isValid()) {
0223     // get the GlabalAlgBlk (Stupid find better way) of BX=0
0224     std::vector<GlobalAlgBlk>::const_iterator algBlk = m_uGtAlgBlk->begin(0);
0225     if (algBlk != m_uGtAlgBlk->end(0)) {
0226       if (!m_readPrescalesFromFile) {
0227         m_PreScaleColumn = static_cast<unsigned int>(algBlk->getPreScColumn());
0228 
0229         // Fix for MC prescale column being set to index+1 in early versions of uGT emulator
0230         if (iEvent.run() == 1) {
0231           if (m_prescaleFactorsAlgoTrig->size() == 1 && m_PreScaleColumn == 1)
0232             m_PreScaleColumn = 0;
0233         }
0234 
0235         // add protection against out-of-bound index for prescale column
0236         if (m_PreScaleColumn >= m_prescaleFactorsAlgoTrig->size()) {
0237           LogDebug("l1t|Global") << "Prescale column extracted from GlobalAlgBlk too large: " << m_PreScaleColumn
0238                                  << "\tMaximum value allowed: " << m_prescaleFactorsAlgoTrig->size() - 1
0239                                  << "\tResetting prescale column to 0" << std::endl;
0240           m_PreScaleColumn = 0;
0241         }
0242       }
0243       const std::vector<double>& prescaleSet = (*m_prescaleFactorsAlgoTrig)[m_PreScaleColumn];
0244 
0245       // Grab the final OR from the AlgBlk,
0246       m_finalOR = algBlk->getFinalOR();
0247 
0248       // Make a map of the trigger name and whether it passed various stages (initial,prescale,final)
0249       // Note: might be able to improve performance by not full remaking map with names each time
0250       for (std::map<std::string, L1TUtmAlgorithm>::const_iterator itAlgo = m_algorithmMap->begin();
0251            itAlgo != m_algorithmMap->end();
0252            itAlgo++) {
0253         // Get the algorithm name
0254         std::string_view algName = itAlgo->first;
0255         int algBit = (itAlgo->second).getIndex();  //algoBitNumber();
0256 
0257         bool decisionInitial = algBlk->getAlgoDecisionInitial(algBit);
0258         (m_decisionsInitial[algBit]).first = algName;
0259         (m_decisionsInitial[algBit]).second = decisionInitial;
0260 
0261         bool decisionInterm = algBlk->getAlgoDecisionInterm(algBit);
0262         (m_decisionsInterm[algBit]).first = algName;
0263         (m_decisionsInterm[algBit]).second = decisionInterm;
0264 
0265         bool decisionFinal = algBlk->getAlgoDecisionFinal(algBit);
0266         (m_decisionsFinal[algBit]).first = algName;
0267         (m_decisionsFinal[algBit]).second = decisionFinal;
0268 
0269         (m_prescales[algBit]).first = algName;
0270         if (size_t(algBit) < prescaleSet.size()) {
0271           (m_prescales[algBit]).second = prescaleSet[algBit];
0272         }
0273       }
0274     } else {
0275       //cout << "Error empty AlgBlk recovered.\n";
0276     }
0277   } else {
0278     //cout<< "Error no valid uGT Algorithm Data with Token provided " << endl;
0279   }
0280 }
0281 
0282 void l1t::L1TGlobalUtil::loadPrescalesAndMasks() {
0283   std::ifstream inputPrescaleFile;
0284   //std::cout << "Reading prescales from file: " << m_preScaleFileName << std::endl;
0285   inputPrescaleFile.open(m_preScaleFileName);
0286 
0287   std::vector<std::vector<int> > vec;
0288   std::vector<std::vector<double> > prescale_vec;
0289 
0290   if (inputPrescaleFile) {
0291     std::string prefix1("#");
0292     std::string prefix2("-1");
0293 
0294     std::string line;
0295 
0296     bool first = true;
0297 
0298     while (getline(inputPrescaleFile, line)) {
0299       if (!line.compare(0, prefix1.size(), prefix1))
0300         continue;
0301       //if( !line.compare(0, prefix2.size(), prefix2) ) continue;
0302 
0303       istringstream split(line);
0304       int value;
0305       int col = 0;
0306       char sep;
0307 
0308       while (split >> value) {
0309         if (first) {
0310           // Each new value read on line 1 should create a new inner vector
0311           vec.push_back(std::vector<int>());
0312         }
0313 
0314         vec[col].push_back(value);
0315         ++col;
0316 
0317         // read past the separator
0318         split >> sep;
0319       }
0320 
0321       // Finished reading line 1 and creating as many inner
0322       // vectors as required
0323       first = false;
0324     }
0325 
0326     int NumPrescaleSets = 0;
0327     for (int iCol = 0; iCol < int(vec.size()); iCol++) {
0328       if (!vec[iCol].empty()) {
0329         int firstRow = vec[iCol][0];
0330 
0331         if (firstRow >= 0)
0332           NumPrescaleSets++;
0333         //else if( firstRow==-2 ) maskColumn = iCol;
0334         //else if( firstRow==-3 ) maskVetoColumn = iCol;
0335       }
0336     }
0337 
0338     //std::cout << "NumPrescaleSets= " << NumPrescaleSets << std::endl;
0339     if (NumPrescaleSets > 0) {
0340       // Fill default prescale set
0341       prescale_vec.reserve(NumPrescaleSets);
0342       for (int iSet = 0; iSet < NumPrescaleSets; iSet++) {
0343         prescale_vec.emplace_back();
0344         prescale_vec.back().reserve(m_numberPhysTriggers);
0345         for (unsigned int iBit = 0; iBit < m_numberPhysTriggers; ++iBit) {
0346           int inputDefaultPrescale = 1;
0347           prescale_vec[iSet].push_back(inputDefaultPrescale);
0348         }
0349       }
0350 
0351       // Fill non-trivial prescale set
0352       for (int iBit = 1; iBit < int(vec[0].size()); iBit++) {
0353         unsigned int algoBit = vec[0][iBit];
0354         // algoBit must be less than the number of triggers
0355         if (algoBit < m_numberPhysTriggers) {
0356           for (int iSet = 0; iSet < int(vec.size()); iSet++) {
0357             int useSet = -1;
0358             if (!vec[iSet].empty()) {
0359               useSet = vec[iSet][0];
0360             }
0361             useSet -= 1;
0362 
0363             if (useSet < 0)
0364               continue;
0365 
0366             int prescale = vec[iSet][iBit];
0367             prescale_vec[useSet][algoBit] = prescale;
0368           }
0369         } else {
0370           LogTrace("l1t|Global") << "\nPrescale file has algo bit: " << algoBit
0371                                  << "\nThis is larger than the number of triggers: " << m_numberPhysTriggers
0372                                  << "\nSomething is wrong. Ignoring." << std::endl;
0373         }
0374       }
0375     }
0376 
0377   } else {
0378     LogTrace("l1t|Global") << "\nCould not find file: " << m_preScaleFileName
0379                            << "\nFilling the prescale vectors with prescale 1"
0380                            << "\nSetting prescale set to 0" << std::endl;
0381 
0382     m_PreScaleColumn = 0;
0383 
0384     prescale_vec.reserve(1);
0385     for (int col = 0; col < 1; col++) {
0386       prescale_vec.emplace_back();
0387       prescale_vec.back().reserve(m_numberPhysTriggers);
0388       for (unsigned int iBit = 0; iBit < m_numberPhysTriggers; ++iBit) {
0389         int inputDefaultPrescale = 0;
0390         prescale_vec[col].push_back(inputDefaultPrescale);
0391       }
0392     }
0393   }
0394 
0395   inputPrescaleFile.close();
0396 
0397   m_initialPrescaleFactorsAlgoTrig = std::move(prescale_vec);
0398   // setting of bx masks from an input file not enabled; do not see a use case at the moment
0399   //std::map<int, std::vector<int> > m_initialTriggerMaskAlgoTrig;
0400 }
0401 
0402 void l1t::L1TGlobalUtil::eventSetupConsumes(edm::ConsumesCollector& iC, UseEventSetupIn useEventSetupIn) {
0403   if (useEventSetupIn == UseEventSetupIn::Run || useEventSetupIn == UseEventSetupIn::RunAndEvent) {
0404     m_L1TUtmTriggerMenuRunToken = iC.esConsumes<L1TUtmTriggerMenu, L1TUtmTriggerMenuRcd, edm::Transition::BeginRun>();
0405     if (!m_readPrescalesFromFile) {
0406       m_L1TGlobalPrescalesVetosFractRunToken =
0407           iC.esConsumes<L1TGlobalPrescalesVetosFract, L1TGlobalPrescalesVetosFractRcd, edm::Transition::BeginRun>();
0408     }
0409   }
0410   if (useEventSetupIn == UseEventSetupIn::Event || useEventSetupIn == UseEventSetupIn::RunAndEvent) {
0411     m_L1TUtmTriggerMenuEventToken = iC.esConsumes<L1TUtmTriggerMenu, L1TUtmTriggerMenuRcd>();
0412     if (!m_readPrescalesFromFile) {
0413       m_L1TGlobalPrescalesVetosFractEventToken =
0414           iC.esConsumes<L1TGlobalPrescalesVetosFract, L1TGlobalPrescalesVetosFractRcd>();
0415     }
0416   }
0417 }
0418 
0419 void l1t::L1TGlobalUtil::resetDecisionVectors() {
0420   // Reset all the vector contents with null information
0421   m_decisionsInitial.clear();
0422   m_decisionsInitial.resize(m_maxNumberPhysTriggers);
0423   m_decisionsInterm.clear();
0424   m_decisionsInterm.resize(m_maxNumberPhysTriggers);
0425   m_decisionsFinal.clear();
0426   m_decisionsFinal.resize(m_maxNumberPhysTriggers);
0427 
0428   for (unsigned int algBit = 0; algBit < m_maxNumberPhysTriggers; algBit++) {
0429     (m_decisionsInitial[algBit]).first = "NULL";
0430     (m_decisionsInitial[algBit]).second = false;
0431 
0432     (m_decisionsInterm[algBit]).first = "NULL";
0433     (m_decisionsInterm[algBit]).second = false;
0434 
0435     (m_decisionsFinal[algBit]).first = "NULL";
0436     (m_decisionsFinal[algBit]).second = false;
0437   }
0438 }
0439 
0440 void l1t::L1TGlobalUtil::resetPrescaleVectors() {
0441   // Reset all the vector contents with null information
0442   m_prescales.clear();
0443   m_prescales.resize(m_maxNumberPhysTriggers);
0444 
0445   for (unsigned int algBit = 0; algBit < m_maxNumberPhysTriggers; algBit++) {
0446     (m_prescales[algBit]).first = "NULL";
0447     (m_prescales[algBit]).second = 1;
0448   }
0449 }
0450 
0451 void l1t::L1TGlobalUtil::resetMaskVectors() {
0452   // Reset all the vector contents with null information
0453   m_masks.clear();
0454   m_masks.resize(m_maxNumberPhysTriggers);
0455 
0456   for (unsigned int algBit = 0; algBit < m_maxNumberPhysTriggers; algBit++) {
0457     (m_masks[algBit]).first = "NULL";
0458     // ccla (m_masks[algBit]).second = true;
0459   }
0460 }
0461 
0462 const bool l1t::L1TGlobalUtil::getAlgBitFromName(const std::string& algName, int& bit) const {
0463   std::map<std::string, L1TUtmAlgorithm>::const_iterator itAlgo = m_algorithmMap->find(algName);
0464   if (itAlgo != m_algorithmMap->end()) {
0465     bit = (itAlgo->second).getIndex();  //algoBitNumber();
0466     return true;
0467   }
0468 
0469   return false;  //did not find anything by that name
0470 }
0471 
0472 const bool l1t::L1TGlobalUtil::getAlgNameFromBit(int& bit, std::string_view& algName) const {
0473   // since we are just looking up the name, doesn't matter which vector we get it from
0474   if ((m_decisionsInitial[bit]).first != "NULL") {
0475     algName = (m_decisionsInitial[bit]).first;
0476     return true;
0477   }
0478   return false;  //No name associated with this bit
0479 }
0480 
0481 const bool l1t::L1TGlobalUtil::getInitialDecisionByBit(int& bit, bool& decision) const {
0482   /*
0483     for(std::vector<GlobalAlgBlk>::const_iterator algBlk = m_uGtAlgBlk->begin(0); algBlk != m_uGtAlgBlk->end(0); ++algBlk) {
0484       decision = algBlk->getAlgoDecisionFinal(bit);
0485     }
0486   */
0487   // Need some check that this is a valid bit
0488   if ((m_decisionsInitial[bit]).first != "NULL") {
0489     decision = (m_decisionsInitial[bit]).second;
0490     return true;
0491   }
0492 
0493   return false;  //couldn't get the information requested.
0494 }
0495 const bool l1t::L1TGlobalUtil::getIntermDecisionByBit(int& bit, bool& decision) const {
0496   // Need some check that this is a valid bit
0497   if ((m_decisionsInterm[bit]).first != "NULL") {
0498     decision = (m_decisionsInterm[bit]).second;
0499     return true;
0500   }
0501 
0502   return false;  //couldn't get the information requested.
0503 }
0504 const bool l1t::L1TGlobalUtil::getFinalDecisionByBit(int& bit, bool& decision) const {
0505   // Need some check that this is a valid bit
0506   if ((m_decisionsFinal[bit]).first != "NULL") {
0507     decision = (m_decisionsFinal[bit]).second;
0508     return true;
0509   }
0510 
0511   return false;  //couldn't get the information requested.
0512 }
0513 const bool l1t::L1TGlobalUtil::getPrescaleByBit(int& bit, double& prescale) const {
0514   // Need some check that this is a valid bit
0515   if ((m_prescales[bit]).first != "NULL") {
0516     prescale = (m_prescales[bit]).second;
0517     return true;
0518   }
0519 
0520   return false;  //couldn't get the information requested.
0521 }
0522 const bool l1t::L1TGlobalUtil::getMaskByBit(int& bit, std::vector<int>& mask) const {
0523   // Need some check that this is a valid bit
0524   if ((m_masks[bit]).first != "NULL") {
0525     mask = (m_masks[bit]).second;
0526     return true;
0527   }
0528 
0529   return false;  //couldn't get the information requested.
0530 }
0531 
0532 const bool l1t::L1TGlobalUtil::getInitialDecisionByName(const std::string& algName, bool& decision) const {
0533   int bit = -1;
0534   if (getAlgBitFromName(algName, bit)) {
0535     decision = (m_decisionsInitial[bit]).second;
0536     return true;
0537   }
0538 
0539   return false;  //trigger name was not the menu.
0540 }
0541 
0542 const bool l1t::L1TGlobalUtil::getIntermDecisionByName(const std::string& algName, bool& decision) const {
0543   int bit = -1;
0544   if (getAlgBitFromName(algName, bit)) {
0545     decision = (m_decisionsInterm[bit]).second;
0546     return true;
0547   }
0548 
0549   return false;  //trigger name was not the menu.
0550 }
0551 
0552 const bool l1t::L1TGlobalUtil::getFinalDecisionByName(const std::string& algName, bool& decision) const {
0553   int bit = -1;
0554   if (getAlgBitFromName(algName, bit)) {
0555     decision = (m_decisionsFinal[bit]).second;
0556     return true;
0557   }
0558 
0559   return false;  //trigger name was not the menu.
0560 }
0561 const bool l1t::L1TGlobalUtil::getPrescaleByName(const std::string& algName, double& prescale) const {
0562   int bit = -1;
0563   if (getAlgBitFromName(algName, bit)) {
0564     prescale = (m_prescales[bit]).second;
0565     return true;
0566   }
0567 
0568   return false;  //trigger name was not the menu.
0569 }
0570 const bool l1t::L1TGlobalUtil::getMaskByName(const std::string& algName, std::vector<int>& mask) const {
0571   int bit = -1;
0572   if (getAlgBitFromName(algName, bit)) {
0573     mask = (m_masks[bit]).second;
0574     return true;
0575   }
0576 
0577   return false;  //trigger name was not the menu.
0578 }