Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-10-08 05:11:42

0001 /**
0002  * \class L1GlobalTriggerFDL
0003  *
0004  *
0005  * Description: Final Decision Logic board.
0006  *
0007  * Implementation:
0008  *    <TODO: enter implementation details>
0009  *
0010  * \author: M. Fierro            - HEPHY Vienna - ORCA version
0011  * \author: Vasile Mihai Ghete   - HEPHY Vienna - CMSSW version
0012  *
0013  *
0014  */
0015 
0016 // this class header
0017 #include "L1Trigger/GlobalTrigger/interface/L1GlobalTriggerFDL.h"
0018 
0019 // system include files
0020 #include <iostream>
0021 
0022 // user include files
0023 #include "DataFormats/L1GlobalTrigger/interface/L1GlobalTriggerReadoutSetup.h"
0024 #include "DataFormats/L1GlobalTrigger/interface/L1GlobalTriggerReadoutSetupFwd.h"
0025 #include "DataFormats/L1GlobalTrigger/interface/L1GtFdlWord.h"
0026 
0027 #include "DataFormats/L1GlobalTrigger/interface/L1GlobalTriggerEvmReadoutRecord.h"
0028 #include "DataFormats/L1GlobalTrigger/interface/L1GlobalTriggerReadoutRecord.h"
0029 
0030 #include "L1Trigger/GlobalTrigger/interface/L1GlobalTriggerGTL.h"
0031 #include "L1Trigger/GlobalTrigger/interface/L1GlobalTriggerPSB.h"
0032 
0033 #include "FWCore/MessageLogger/interface/MessageDrop.h"
0034 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0035 
0036 #include "FWCore/Framework/interface/Event.h"
0037 
0038 #include "CondFormats/L1TObjects/interface/L1GtBoard.h"
0039 #include "CondFormats/L1TObjects/interface/L1GtFwd.h"
0040 
0041 // forward declarations
0042 
0043 // constructor
0044 L1GlobalTriggerFDL::L1GlobalTriggerFDL()
0045     :  // logical switches
0046       m_firstEv(true),
0047       m_firstEvLumiSegment(true),
0048       m_firstEvRun(true),
0049       m_isDebugEnabled(edm::isDebugEnabled()) {
0050   // create empty FDL word
0051   m_gtFdlWord = new L1GtFdlWord();
0052 
0053   // can not reserve memory here for prescale counters - no access to EventSetup
0054 }
0055 
0056 // destructor
0057 L1GlobalTriggerFDL::~L1GlobalTriggerFDL() {
0058   reset();
0059   delete m_gtFdlWord;
0060 }
0061 
0062 // Operations
0063 
0064 // run FDL
0065 void L1GlobalTriggerFDL::run(edm::Event &iEvent,
0066                              const std::vector<int> &prescaleFactorsAlgoTrig,
0067                              const std::vector<int> &prescaleFactorsTechTrig,
0068                              const std::vector<unsigned int> &triggerMaskAlgoTrig,
0069                              const std::vector<unsigned int> &triggerMaskTechTrig,
0070                              const std::vector<unsigned int> &triggerMaskVetoAlgoTrig,
0071                              const std::vector<unsigned int> &triggerMaskVetoTechTrig,
0072                              const std::vector<L1GtBoard> &boardMaps,
0073                              const int totalBxInEvent,
0074                              const int iBxInEvent,
0075                              const unsigned int numberPhysTriggers,
0076                              const unsigned int numberTechnicalTriggers,
0077                              const unsigned int numberDaqPartitions,
0078                              const L1GlobalTriggerGTL *ptrGTL,
0079                              const L1GlobalTriggerPSB *ptrPSB,
0080                              const int pfAlgoSetIndex,
0081                              const int pfTechSetIndex,
0082                              const bool algorithmTriggersUnprescaled,
0083                              const bool algorithmTriggersUnmasked,
0084                              const bool technicalTriggersUnprescaled,
0085                              const bool technicalTriggersUnmasked,
0086                              const bool technicalTriggersVetoUnmasked) {
0087   // FIXME get rid of bitset in GTL in order to use only EventSetup
0088   const unsigned int numberPhysTriggersSet = L1GlobalTriggerReadoutSetup::NumberPhysTriggers;
0089 
0090   // get gtlDecisionWord from GTL
0091   std::bitset<numberPhysTriggersSet> gtlDecisionWord = ptrGTL->getAlgorithmOR();
0092 
0093   // convert decision word from std::bitset to std::vector<bool>
0094   DecisionWord algoDecisionWord(numberPhysTriggers);
0095 
0096   for (unsigned int iBit = 0; iBit < numberPhysTriggers; ++iBit) {
0097     bool bitValue = gtlDecisionWord.test(iBit);
0098     algoDecisionWord[iBit] = bitValue;
0099   }
0100 
0101   // prescale counters are reset at the beginning of the luminosity segment
0102 
0103   if (m_firstEv) {
0104     // prescale counters: numberPhysTriggers counters per bunch cross
0105     m_prescaleCounterAlgoTrig.reserve(numberPhysTriggers * totalBxInEvent);
0106 
0107     for (int iBxInEvent = 0; iBxInEvent <= totalBxInEvent; ++iBxInEvent) {
0108       m_prescaleCounterAlgoTrig.push_back(prescaleFactorsAlgoTrig);
0109     }
0110 
0111     // prescale counters: numberTechnicalTriggers counters per bunch cross
0112     m_prescaleCounterTechTrig.reserve(numberTechnicalTriggers * totalBxInEvent);
0113 
0114     for (int iBxInEvent = 0; iBxInEvent <= totalBxInEvent; ++iBxInEvent) {
0115       m_prescaleCounterTechTrig.push_back(prescaleFactorsTechTrig);
0116     }
0117 
0118     m_firstEv = false;
0119   }
0120 
0121   // TODO FIXME find the beginning of the luminosity segment
0122   if (m_firstEvLumiSegment) {
0123     m_prescaleCounterAlgoTrig.clear();
0124     for (int iBxInEvent = 0; iBxInEvent <= totalBxInEvent; ++iBxInEvent) {
0125       m_prescaleCounterAlgoTrig.push_back(prescaleFactorsAlgoTrig);
0126     }
0127 
0128     m_prescaleCounterTechTrig.clear();
0129     for (int iBxInEvent = 0; iBxInEvent <= totalBxInEvent; ++iBxInEvent) {
0130       m_prescaleCounterTechTrig.push_back(prescaleFactorsTechTrig);
0131     }
0132 
0133     m_firstEvLumiSegment = false;
0134   }
0135 
0136   // prescale the algorithm, if necessary
0137 
0138   // iBxInEvent is ... -2 -1 0 1 2 ... while counters are 0 1 2 3 4 ...
0139   int inBxInEvent = totalBxInEvent / 2 + iBxInEvent;
0140 
0141   for (unsigned int iBit = 0; iBit < numberPhysTriggers; ++iBit) {
0142     if ((!algorithmTriggersUnprescaled) && (prescaleFactorsAlgoTrig.at(iBit) != 1)) {
0143       bool bitValue = algoDecisionWord.at(iBit);
0144       if (bitValue) {
0145         (m_prescaleCounterAlgoTrig.at(inBxInEvent).at(iBit))--;
0146         if (m_prescaleCounterAlgoTrig.at(inBxInEvent).at(iBit) == 0) {
0147           // bit already true in algoDecisionWord, just reset counter
0148           m_prescaleCounterAlgoTrig.at(inBxInEvent).at(iBit) = prescaleFactorsAlgoTrig.at(iBit);
0149 
0150           // LogTrace("L1GlobalTrigger")
0151           //<< "\nPrescaled algorithm: " << iBit << ". Reset counter to "
0152           //<< prescaleFactorsAlgoTrig.at(iBit) << "\n"
0153           //<< std::endl;
0154 
0155         } else {
0156           // change bit to false
0157           algoDecisionWord[iBit] = false;
0158           ;
0159 
0160           // LogTrace("L1GlobalTrigger")
0161           //<< "\nPrescaled algorithm: " << iBit << ". Result set to false"
0162           //<< std::endl;
0163         }
0164       }
0165     }
0166   }
0167 
0168   // algo decision word written in the FDL readout before the trigger mask
0169   // in order to allow multiple DAQ partitions
0170 
0171   //
0172   // technical triggers
0173   //
0174 
0175   std::vector<bool> techDecisionWord = *(ptrPSB->getGtTechnicalTriggers());
0176 
0177   // prescale the technical trigger, if necessary
0178 
0179   for (unsigned int iBit = 0; iBit < numberTechnicalTriggers; ++iBit) {
0180     if ((!technicalTriggersUnprescaled) && (prescaleFactorsTechTrig.at(iBit) != 1)) {
0181       bool bitValue = techDecisionWord.at(iBit);
0182       if (bitValue) {
0183         (m_prescaleCounterTechTrig.at(inBxInEvent).at(iBit))--;
0184         if (m_prescaleCounterTechTrig.at(inBxInEvent).at(iBit) == 0) {
0185           // bit already true in techDecisionWord, just reset counter
0186           m_prescaleCounterTechTrig.at(inBxInEvent).at(iBit) = prescaleFactorsTechTrig.at(iBit);
0187 
0188           // LogTrace("L1GlobalTrigger")
0189           //<< "\nPrescaled algorithm: " << iBit << ". Reset counter to "
0190           //<< prescaleFactorsTechTrig.at(iBit) << "\n"
0191           //<< std::endl;
0192 
0193         } else {
0194           // change bit to false
0195           techDecisionWord[iBit] = false;
0196 
0197           // LogTrace("L1GlobalTrigger")
0198           //<< "\nPrescaled technical trigger: " << iBit << ". Result set to
0199           // false"
0200           //<< std::endl;
0201         }
0202       }
0203     }
0204   }
0205 
0206   // technical trigger decision word written in the FDL readout before the
0207   // trigger mask in order to allow multiple DAQ partitions
0208 
0209   //
0210   // compute the final decision word per DAQ partition
0211   //
0212 
0213   uint16_t finalOrValue = 0;
0214 
0215   for (unsigned int iDaq = 0; iDaq < numberDaqPartitions; ++iDaq) {
0216     bool daqPartitionFinalOR = false;
0217 
0218     // starts with technical trigger veto mask to minimize computation
0219     // no algorithm trigger veto mask is implemented up to now in hardware,
0220     // therefore do not implement it here
0221     bool vetoTechTrig = false;
0222 
0223     // vetoTechTrig can change only when using trigger veto masks
0224     if (!technicalTriggersVetoUnmasked) {
0225       for (unsigned int iBit = 0; iBit < numberTechnicalTriggers; ++iBit) {
0226         int triggerMaskVetoTechTrigBit = triggerMaskVetoTechTrig[iBit] & (1 << iDaq);
0227         // LogTrace("L1GlobalTrigger")
0228         //<< "\nTechnical trigger bit: " << iBit
0229         //<< " mask = " << triggerMaskVetoTechTrigBit
0230         //<< " DAQ partition " << iDaq
0231         //<< std::endl;
0232 
0233         if (triggerMaskVetoTechTrigBit && techDecisionWord[iBit]) {
0234           vetoTechTrig = true;
0235 
0236           // LogTrace("L1GlobalTrigger")
0237           //<< "\nVeto mask technical trigger: " << iBit
0238           // << ". FinalOR for DAQ partition " << iDaq << " set to false"
0239           //<< std::endl;
0240 
0241           break;
0242         }
0243       }
0244     }
0245 
0246     // apply algorithm and technical trigger masks only if no veto from
0247     // technical trigger
0248     if (!vetoTechTrig) {
0249       // algorithm trigger mask
0250       bool algoFinalOr = false;
0251 
0252       for (unsigned int iBit = 0; iBit < numberPhysTriggers; ++iBit) {
0253         bool iBitDecision = false;
0254 
0255         int triggerMaskAlgoTrigBit = -1;
0256 
0257         if (algorithmTriggersUnmasked) {
0258           triggerMaskAlgoTrigBit = 0;
0259         } else {
0260           triggerMaskAlgoTrigBit = triggerMaskAlgoTrig[iBit] & (1 << iDaq);
0261         }
0262         // LogTrace("L1GlobalTrigger")
0263         //<< "\nAlgorithm trigger bit: " << iBit
0264         //<< " mask = " << triggerMaskAlgoTrigBit
0265         //<< " DAQ partition " << iDaq
0266         //<< std::endl;
0267 
0268         if (triggerMaskAlgoTrigBit) {
0269           iBitDecision = false;
0270 
0271           // LogTrace("L1GlobalTrigger")
0272           //<< "\nMasked algorithm trigger: " << iBit << ". Result set to false"
0273           //<< std::endl;
0274         } else {
0275           iBitDecision = algoDecisionWord[iBit];
0276         }
0277 
0278         algoFinalOr = algoFinalOr || iBitDecision;
0279       }
0280 
0281       // set the technical trigger mask: block the corresponding algorithm if
0282       // bit value is 1
0283 
0284       bool techFinalOr = false;
0285 
0286       for (unsigned int iBit = 0; iBit < numberTechnicalTriggers; ++iBit) {
0287         bool iBitDecision = false;
0288 
0289         int triggerMaskTechTrigBit = -1;
0290 
0291         if (technicalTriggersUnmasked) {
0292           triggerMaskTechTrigBit = 0;
0293         } else {
0294           triggerMaskTechTrigBit = triggerMaskTechTrig[iBit] & (1 << iDaq);
0295         }
0296         // LogTrace("L1GlobalTrigger")
0297         //<< "\nTechnical trigger bit: " << iBit
0298         //<< " mask = " << triggerMaskTechTrigBit
0299         //<< std::endl;
0300 
0301         if (triggerMaskTechTrigBit) {
0302           iBitDecision = false;
0303 
0304           // LogTrace("L1GlobalTrigger")
0305           //<< "\nMasked technical trigger: " << iBit << ". Result set to false"
0306           //<< std::endl;
0307         } else {
0308           iBitDecision = techDecisionWord[iBit];
0309         }
0310 
0311         techFinalOr = techFinalOr || iBitDecision;
0312       }
0313 
0314       daqPartitionFinalOR = algoFinalOr || techFinalOr;
0315 
0316     } else {
0317       daqPartitionFinalOR = false;  // vetoTechTrig
0318     }
0319 
0320     // push it in finalOrValue
0321     uint16_t daqPartitionFinalORValue = static_cast<uint16_t>(daqPartitionFinalOR);
0322 
0323     finalOrValue = finalOrValue | (daqPartitionFinalORValue << iDaq);
0324   }
0325 
0326   // fill everything we know in the L1GtFdlWord
0327 
0328   typedef std::vector<L1GtBoard>::const_iterator CItBoardMaps;
0329   for (CItBoardMaps itBoard = boardMaps.begin(); itBoard != boardMaps.end(); ++itBoard) {
0330     if ((itBoard->gtBoardType() == FDL)) {
0331       m_gtFdlWord->setBoardId(itBoard->gtBoardId());
0332 
0333       // BxInEvent
0334       m_gtFdlWord->setBxInEvent(iBxInEvent);
0335 
0336       // bunch crossing
0337 
0338       // fill in emulator the same bunch crossing (12 bits - hardwired number of
0339       // bits...) and the same local bunch crossing for all boards
0340       int bxCross = iEvent.bunchCrossing();
0341       uint16_t bxCrossHw = 0;
0342       if ((bxCross & 0xFFF) == bxCross) {
0343         bxCrossHw = static_cast<uint16_t>(bxCross);
0344       } else {
0345         bxCrossHw = 0;  // Bx number too large, set to 0!
0346         if (m_verbosity) {
0347           LogDebug("L1GlobalTrigger") << "\nBunch cross number [hex] = " << std::hex << bxCross
0348                                       << "\n  larger than 12 bits. Set to 0! \n"
0349                                       << std::dec << std::endl;
0350         }
0351       }
0352 
0353       m_gtFdlWord->setBxNr(bxCrossHw);
0354 
0355       // set event number since last L1 reset generated in FDL
0356       m_gtFdlWord->setEventNr(static_cast<uint32_t>(iEvent.id().event()));
0357 
0358       // technical trigger decision word
0359       m_gtFdlWord->setGtTechnicalTriggerWord(techDecisionWord);
0360 
0361       // algorithm trigger decision word
0362       m_gtFdlWord->setGtDecisionWord(algoDecisionWord);
0363 
0364       // index of prescale factor set - technical triggers and algo
0365       m_gtFdlWord->setGtPrescaleFactorIndexTech(static_cast<uint16_t>(pfTechSetIndex));
0366       m_gtFdlWord->setGtPrescaleFactorIndexAlgo(static_cast<uint16_t>(pfAlgoSetIndex));
0367 
0368       // NoAlgo bit FIXME
0369 
0370       // finalOR
0371       m_gtFdlWord->setFinalOR(finalOrValue);
0372 
0373       // orbit number
0374       m_gtFdlWord->setOrbitNr(static_cast<uint32_t>(iEvent.orbitNumber()));
0375 
0376       // luminosity segment number
0377       m_gtFdlWord->setLumiSegmentNr(static_cast<uint16_t>(iEvent.luminosityBlock()));
0378 
0379       // local bunch crossing - set identical with absolute BxNr
0380       m_gtFdlWord->setLocalBxNr(bxCrossHw);
0381     }
0382   }
0383 }
0384 
0385 // fill the FDL block in the L1 GT DAQ record for iBxInEvent
0386 void L1GlobalTriggerFDL::fillDaqFdlBlock(const int iBxInEvent,
0387                                          const uint16_t &activeBoardsGtDaq,
0388                                          const int recordLength0,
0389                                          const int recordLength1,
0390                                          const unsigned int altNrBxBoardDaq,
0391                                          const std::vector<L1GtBoard> &boardMaps,
0392                                          L1GlobalTriggerReadoutRecord *gtDaqReadoutRecord) {
0393   typedef std::vector<L1GtBoard>::const_iterator CItBoardMaps;
0394   for (CItBoardMaps itBoard = boardMaps.begin(); itBoard != boardMaps.end(); ++itBoard) {
0395     int iPosition = itBoard->gtPositionDaqRecord();
0396     if (iPosition > 0) {
0397       int iActiveBit = itBoard->gtBitDaqActiveBoards();
0398       bool activeBoard = false;
0399       bool writeBoard = false;
0400 
0401       int recLength = -1;
0402 
0403       if (iActiveBit >= 0) {
0404         activeBoard = activeBoardsGtDaq & (1 << iActiveBit);
0405 
0406         int altNrBxBoard = (altNrBxBoardDaq & (1 << iActiveBit)) >> iActiveBit;
0407 
0408         if (altNrBxBoard == 1) {
0409           recLength = recordLength1;
0410         } else {
0411           recLength = recordLength0;
0412         }
0413 
0414         int lowBxInEvent = (recLength + 1) / 2 - recLength;
0415         int uppBxInEvent = (recLength + 1) / 2 - 1;
0416 
0417         if ((iBxInEvent >= lowBxInEvent) && (iBxInEvent <= uppBxInEvent)) {
0418           writeBoard = true;
0419         }
0420       }
0421 
0422       if (activeBoard && writeBoard && (itBoard->gtBoardType() == FDL)) {
0423         gtDaqReadoutRecord->setGtFdlWord(*m_gtFdlWord);
0424       }
0425     }
0426   }
0427 }
0428 
0429 // fill the FDL block in the L1 GT EVM record for iBxInEvent
0430 void L1GlobalTriggerFDL::fillEvmFdlBlock(const int iBxInEvent,
0431                                          const uint16_t &activeBoardsGtEvm,
0432                                          const int recordLength0,
0433                                          const int recordLength1,
0434                                          const unsigned int altNrBxBoardEvm,
0435                                          const std::vector<L1GtBoard> &boardMaps,
0436                                          L1GlobalTriggerEvmReadoutRecord *gtEvmReadoutRecord) {
0437   typedef std::vector<L1GtBoard>::const_iterator CItBoardMaps;
0438   for (CItBoardMaps itBoard = boardMaps.begin(); itBoard != boardMaps.end(); ++itBoard) {
0439     int iPosition = itBoard->gtPositionEvmRecord();
0440     if (iPosition > 0) {
0441       int iActiveBit = itBoard->gtBitEvmActiveBoards();
0442       bool activeBoard = false;
0443 
0444       if (iActiveBit >= 0) {
0445         activeBoard = activeBoardsGtEvm & (1 << iActiveBit);
0446       }
0447 
0448       if (activeBoard && (itBoard->gtBoardType() == FDL)) {
0449         gtEvmReadoutRecord->setGtFdlWord(*m_gtFdlWord);
0450       }
0451     }
0452   }
0453 }
0454 
0455 // clear FDL
0456 void L1GlobalTriggerFDL::reset() {
0457   m_gtFdlWord->reset();
0458 
0459   // do NOT reset the prescale counters
0460 }