Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:20:03

0001 /**
0002  * \class L1GtTrigReport
0003  *
0004  *
0005  * Description: L1 Trigger report.
0006  *
0007  * Implementation:
0008  *    <TODO: enter implementation details>
0009  *
0010  * \author: Vasile Mihai Ghete - HEPHY Vienna
0011  *
0012  *
0013  */
0014 
0015 // this class header
0016 #include "L1Trigger/GlobalTriggerAnalyzer/interface/L1GtTrigReport.h"
0017 
0018 // system include files
0019 #include <memory>
0020 
0021 #include <iostream>
0022 #include <iomanip>
0023 
0024 #include <map>
0025 #include <set>
0026 #include <cmath>
0027 #include <string>
0028 
0029 // user include files
0030 #include "DataFormats/L1GlobalTrigger/interface/L1GlobalTriggerReadoutSetupFwd.h"
0031 
0032 #include "FWCore/Framework/interface/Frameworkfwd.h"
0033 
0034 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0035 
0036 #include "FWCore/Framework/interface/Event.h"
0037 #include "FWCore/Framework/interface/MakerMacros.h"
0038 
0039 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0040 #include "FWCore/Utilities/interface/InputTag.h"
0041 
0042 #include "FWCore/Framework/interface/EventSetup.h"
0043 
0044 #include "CondFormats/L1TObjects/interface/L1GtStableParameters.h"
0045 #include "CondFormats/DataRecord/interface/L1GtStableParametersRcd.h"
0046 
0047 #include "CondFormats/L1TObjects/interface/L1GtPrescaleFactors.h"
0048 #include "CondFormats/DataRecord/interface/L1GtPrescaleFactorsAlgoTrigRcd.h"
0049 #include "CondFormats/DataRecord/interface/L1GtPrescaleFactorsTechTrigRcd.h"
0050 
0051 #include "CondFormats/L1TObjects/interface/L1GtTriggerMask.h"
0052 #include "CondFormats/DataRecord/interface/L1GtTriggerMaskAlgoTrigRcd.h"
0053 #include "CondFormats/DataRecord/interface/L1GtTriggerMaskTechTrigRcd.h"
0054 
0055 #include "CondFormats/DataRecord/interface/L1GtTriggerMaskVetoAlgoTrigRcd.h"
0056 #include "CondFormats/DataRecord/interface/L1GtTriggerMaskVetoTechTrigRcd.h"
0057 
0058 #include "CondFormats/L1TObjects/interface/L1GtTriggerMenu.h"
0059 #include "CondFormats/DataRecord/interface/L1GtTriggerMenuRcd.h"
0060 
0061 #include "L1Trigger/GlobalTriggerAnalyzer/interface/L1GtTrigReportEntry.h"
0062 
0063 // constructor(s)
0064 L1GtTrigReport::L1GtTrigReport(const edm::ParameterSet& pSet)
0065     :
0066 
0067       // initialize cached IDs
0068 
0069       //
0070       m_l1GtStableParCacheID(0ULL),
0071 
0072       m_numberPhysTriggers(0),
0073       m_numberTechnicalTriggers(0),
0074       m_numberDaqPartitions(0),
0075       m_numberDaqPartitionsMax(0),
0076 
0077       //
0078       m_l1GtPfAlgoCacheID(0ULL),
0079       m_l1GtPfTechCacheID(0ULL),
0080 
0081       m_l1GtTmAlgoCacheID(0ULL),
0082       m_l1GtTmTechCacheID(0ULL),
0083 
0084       m_l1GtTmVetoAlgoCacheID(0ULL),
0085       m_l1GtTmVetoTechCacheID(0ULL),
0086 
0087       //
0088       m_l1GtMenuCacheID(0ULL),
0089 
0090       // boolean flag to select the input record
0091       // if true, it will use L1GlobalTriggerRecord
0092       m_useL1GlobalTriggerRecord(pSet.getParameter<bool>("UseL1GlobalTriggerRecord")),
0093 
0094       /// input tag for GT record (L1 GT DAQ record or L1 GT "lite" record):
0095       m_l1GtRecordInputTag(pSet.getParameter<edm::InputTag>("L1GtRecordInputTag")),
0096       m_l1GtRecordInputToken1(m_useL1GlobalTriggerRecord ? consumes<L1GlobalTriggerRecord>(m_l1GtRecordInputTag)
0097                                                          : edm::EDGetTokenT<L1GlobalTriggerRecord>()),
0098       m_l1GtRecordInputToken2(not m_useL1GlobalTriggerRecord
0099                                   ? consumes<L1GlobalTriggerReadoutRecord>(m_l1GtRecordInputTag)
0100                                   : edm::EDGetTokenT<L1GlobalTriggerReadoutRecord>()),
0101 
0102       m_stableParToken(esConsumes()),
0103       m_pfAlgoToken(esConsumes()),
0104       m_pfTechToken(esConsumes()),
0105       m_tmAlgoToken(esConsumes()),
0106       m_tmTechToken(esConsumes()),
0107       m_tmVetoAlgoToken(esConsumes()),
0108       m_tmVetoTechToken(esConsumes()),
0109       m_menuToken(esConsumes()),
0110 
0111       // print verbosity
0112       m_printVerbosity(pSet.getUntrackedParameter<int>("PrintVerbosity", 2)),
0113 
0114       // print output
0115       m_printOutput(pSet.getUntrackedParameter<int>("PrintOutput", 3)),
0116 
0117       // initialize global counters
0118 
0119       // number of events processed
0120       m_totalEvents(0),
0121 
0122       //
0123       m_entryList(),
0124       m_entryListTechTrig(),
0125 
0126       // set the index of physics DAQ partition TODO input parameter?
0127       m_physicsDaqPartition(0)
0128 
0129 {
0130   LogDebug("L1GtTrigReport") << "\n  Use L1GlobalTriggerRecord:   " << m_useL1GlobalTriggerRecord
0131                              << "\n   (if false: L1GtTrigReport uses L1GlobalTriggerReadoutRecord.)"
0132                              << "\n  Input tag for L1 GT record:  " << m_l1GtRecordInputTag.label() << " \n"
0133                              << "\n  Print verbosity level:           " << m_printVerbosity << " \n"
0134                              << "\n  Print output:                    " << m_printOutput << " \n"
0135                              << std::endl;
0136 }
0137 
0138 // destructor
0139 L1GtTrigReport::~L1GtTrigReport() {
0140   for (ItEntry itEntry = m_entryList.begin(); itEntry != m_entryList.end(); itEntry++) {
0141     if (*itEntry != nullptr) {
0142       delete *itEntry;
0143       *itEntry = nullptr;
0144     }
0145   }
0146 
0147   m_entryList.clear();
0148 
0149   for (ItEntry itEntry = m_entryListTechTrig.begin(); itEntry != m_entryListTechTrig.end(); itEntry++) {
0150     if (*itEntry != nullptr) {
0151       delete *itEntry;
0152       *itEntry = nullptr;
0153     }
0154   }
0155 
0156   m_entryListTechTrig.clear();
0157 }
0158 
0159 // member functions
0160 
0161 // method called once each job just before starting event loop
0162 void L1GtTrigReport::beginJob() {
0163   // empty
0164 }
0165 
0166 // analyze each event
0167 void L1GtTrigReport::analyze(const edm::Event& iEvent, const edm::EventSetup& evSetup) {
0168   // increase the number of processed events
0169   m_totalEvents++;
0170 
0171   // get / update the stable parameters from the EventSetup
0172   // local cache & check on cacheIdentifier
0173 
0174   unsigned long long l1GtStableParCacheID = evSetup.get<L1GtStableParametersRcd>().cacheIdentifier();
0175 
0176   if (m_l1GtStableParCacheID != l1GtStableParCacheID) {
0177     m_l1GtStablePar = &evSetup.getData(m_stableParToken);
0178 
0179     // number of physics triggers
0180     m_numberPhysTriggers = m_l1GtStablePar->gtNumberPhysTriggers();
0181 
0182     // number of technical triggers
0183     m_numberTechnicalTriggers = m_l1GtStablePar->gtNumberTechnicalTriggers();
0184 
0185     // number of DAQ partitions
0186     m_numberDaqPartitions = 8;  // FIXME add it to stable parameters
0187 
0188     if (m_numberDaqPartitionsMax < m_numberDaqPartitions) {
0189       int numberDaqPartitionsOld = m_numberDaqPartitionsMax;
0190       m_numberDaqPartitionsMax = m_numberDaqPartitions;
0191 
0192       m_globalNrErrors.reserve(m_numberDaqPartitionsMax);
0193       m_globalNrAccepts.reserve(m_numberDaqPartitionsMax);
0194 
0195       for (unsigned int iDaq = numberDaqPartitionsOld; iDaq < m_numberDaqPartitionsMax; ++iDaq) {
0196         m_globalNrErrors.push_back(0);
0197         m_globalNrAccepts.push_back(0);
0198       }
0199     }
0200 
0201     //
0202     m_l1GtStableParCacheID = l1GtStableParCacheID;
0203   }
0204 
0205   // get / update the prescale factors from the EventSetup
0206   // local cache & check on cacheIdentifier
0207 
0208   unsigned long long l1GtPfAlgoCacheID = evSetup.get<L1GtPrescaleFactorsAlgoTrigRcd>().cacheIdentifier();
0209 
0210   if (m_l1GtPfAlgoCacheID != l1GtPfAlgoCacheID) {
0211     m_l1GtPfAlgo = &evSetup.getData(m_pfAlgoToken);
0212 
0213     m_prescaleFactorsAlgoTrig = &(m_l1GtPfAlgo->gtPrescaleFactors());
0214 
0215     m_l1GtPfAlgoCacheID = l1GtPfAlgoCacheID;
0216   }
0217 
0218   unsigned long long l1GtPfTechCacheID = evSetup.get<L1GtPrescaleFactorsTechTrigRcd>().cacheIdentifier();
0219 
0220   if (m_l1GtPfTechCacheID != l1GtPfTechCacheID) {
0221     m_l1GtPfTech = &evSetup.getData(m_pfTechToken);
0222 
0223     m_prescaleFactorsTechTrig = &(m_l1GtPfTech->gtPrescaleFactors());
0224 
0225     m_l1GtPfTechCacheID = l1GtPfTechCacheID;
0226   }
0227 
0228   // get / update the trigger mask from the EventSetup
0229   // local cache & check on cacheIdentifier
0230 
0231   unsigned long long l1GtTmAlgoCacheID = evSetup.get<L1GtTriggerMaskAlgoTrigRcd>().cacheIdentifier();
0232 
0233   if (m_l1GtTmAlgoCacheID != l1GtTmAlgoCacheID) {
0234     m_l1GtTmAlgo = &evSetup.getData(m_tmAlgoToken);
0235 
0236     m_triggerMaskAlgoTrig = m_l1GtTmAlgo->gtTriggerMask();
0237 
0238     m_l1GtTmAlgoCacheID = l1GtTmAlgoCacheID;
0239   }
0240 
0241   unsigned long long l1GtTmTechCacheID = evSetup.get<L1GtTriggerMaskTechTrigRcd>().cacheIdentifier();
0242 
0243   if (m_l1GtTmTechCacheID != l1GtTmTechCacheID) {
0244     m_l1GtTmTech = &evSetup.getData(m_tmTechToken);
0245 
0246     m_triggerMaskTechTrig = m_l1GtTmTech->gtTriggerMask();
0247 
0248     m_l1GtTmTechCacheID = l1GtTmTechCacheID;
0249   }
0250 
0251   unsigned long long l1GtTmVetoAlgoCacheID = evSetup.get<L1GtTriggerMaskVetoAlgoTrigRcd>().cacheIdentifier();
0252 
0253   if (m_l1GtTmVetoAlgoCacheID != l1GtTmVetoAlgoCacheID) {
0254     m_l1GtTmVetoAlgo = &evSetup.getData(m_tmVetoAlgoToken);
0255 
0256     m_triggerMaskVetoAlgoTrig = m_l1GtTmVetoAlgo->gtTriggerMask();
0257 
0258     m_l1GtTmVetoAlgoCacheID = l1GtTmVetoAlgoCacheID;
0259   }
0260 
0261   unsigned long long l1GtTmVetoTechCacheID = evSetup.get<L1GtTriggerMaskVetoTechTrigRcd>().cacheIdentifier();
0262 
0263   if (m_l1GtTmVetoTechCacheID != l1GtTmVetoTechCacheID) {
0264     m_l1GtTmVetoTech = &evSetup.getData(m_tmVetoTechToken);
0265 
0266     m_triggerMaskVetoTechTrig = m_l1GtTmVetoTech->gtTriggerMask();
0267 
0268     m_l1GtTmVetoTechCacheID = l1GtTmVetoTechCacheID;
0269   }
0270 
0271   // get / update the trigger menu from the EventSetup
0272   // local cache & check on cacheIdentifier
0273 
0274   unsigned long long l1GtMenuCacheID = evSetup.get<L1GtTriggerMenuRcd>().cacheIdentifier();
0275 
0276   if (m_l1GtMenuCacheID != l1GtMenuCacheID) {
0277     m_l1GtMenu = &evSetup.getData(m_menuToken);
0278 
0279     m_l1GtMenuCacheID = l1GtMenuCacheID;
0280 
0281     LogDebug("L1GtTrigReport") << "\n  Changing L1 menu to : \n"
0282                                << m_l1GtMenu->gtTriggerMenuName() << "\n\n"
0283                                << std::endl;
0284   }
0285 
0286   const AlgorithmMap& algorithmMap = m_l1GtMenu->gtAlgorithmMap();
0287   const AlgorithmMap& technicalTriggerMap = m_l1GtMenu->gtTechnicalTriggerMap();
0288 
0289   const std::string& menuName = m_l1GtMenu->gtTriggerMenuName();
0290 
0291   // ... end EventSetup
0292 
0293   // get L1GlobalTriggerReadoutRecord or L1GlobalTriggerRecord
0294   // in L1GlobalTriggerRecord, only the physics partition is available
0295   edm::Handle<L1GlobalTriggerReadoutRecord> gtReadoutRecord;
0296   edm::Handle<L1GlobalTriggerRecord> gtRecord;
0297 
0298   if (m_useL1GlobalTriggerRecord) {
0299     iEvent.getByToken(m_l1GtRecordInputToken1, gtRecord);
0300   } else {
0301     iEvent.getByToken(m_l1GtRecordInputToken2, gtReadoutRecord);
0302   }
0303 
0304   bool validRecord = false;
0305 
0306   unsigned int pfIndexAlgo = 0;  // get them later from the record
0307   unsigned int pfIndexTech = 0;
0308 
0309   DecisionWord gtDecisionWordBeforeMask;
0310   DecisionWord gtDecisionWordAfterMask;
0311 
0312   TechnicalTriggerWord technicalTriggerWordBeforeMask;
0313   TechnicalTriggerWord technicalTriggerWordAfterMask;
0314 
0315   if (m_useL1GlobalTriggerRecord) {
0316     if (gtRecord.isValid()) {
0317       // get Global Trigger decision and the decision word
0318       bool gtDecision = gtRecord->decision();
0319 
0320       gtDecisionWordBeforeMask = gtRecord->decisionWordBeforeMask();
0321       gtDecisionWordAfterMask = gtRecord->decisionWord();
0322 
0323       technicalTriggerWordBeforeMask = gtRecord->technicalTriggerWordBeforeMask();
0324       technicalTriggerWordAfterMask = gtRecord->technicalTriggerWord();
0325 
0326       if (gtDecision) {
0327         m_globalNrAccepts[m_physicsDaqPartition]++;
0328       }
0329 
0330       pfIndexAlgo = gtRecord->gtPrescaleFactorIndexAlgo();
0331       pfIndexTech = gtRecord->gtPrescaleFactorIndexTech();
0332 
0333       validRecord = true;
0334 
0335     } else {
0336       m_globalNrErrors[m_physicsDaqPartition]++;
0337 
0338       edm::LogWarning("L1GtTrigReport") << "\n  L1GlobalTriggerRecord with input tag " << m_l1GtRecordInputTag.label()
0339                                         << " not found."
0340                                         << "\n  Event classified as Error\n\n"
0341                                         << std::endl;
0342     }
0343 
0344   } else {
0345     if (gtReadoutRecord.isValid()) {
0346       // check if the readout record has size greater than zero, then proceeds
0347       const std::vector<L1GtFdlWord>& fdlVec = gtReadoutRecord->gtFdlVector();
0348       size_t fdlVecSize = fdlVec.size();
0349 
0350       if (fdlVecSize > 0) {
0351         LogDebug("L1GtTrigReport") << "\n  L1GlobalTriggerReadoutRecord with input tag " << m_l1GtRecordInputTag.label()
0352                                    << " has gtFdlVector of size " << fdlVecSize << std::endl;
0353 
0354         // get Global Trigger finalOR and the decision word
0355         uint16_t gtFinalOR = gtReadoutRecord->finalOR();
0356 
0357         gtDecisionWordBeforeMask = gtReadoutRecord->decisionWord();
0358         technicalTriggerWordBeforeMask = gtReadoutRecord->technicalTriggerWord();
0359 
0360         for (unsigned int iDaqPartition = 0; iDaqPartition < m_numberDaqPartitions; ++iDaqPartition) {
0361           bool gtDecision = static_cast<bool>(gtFinalOR & (1 << iDaqPartition));
0362           if (gtDecision) {
0363             m_globalNrAccepts[iDaqPartition]++;
0364           }
0365         }
0366 
0367         pfIndexAlgo = static_cast<unsigned int>((gtReadoutRecord->gtFdlWord()).gtPrescaleFactorIndexAlgo());
0368         pfIndexTech = static_cast<unsigned int>((gtReadoutRecord->gtFdlWord()).gtPrescaleFactorIndexTech());
0369 
0370         validRecord = true;
0371 
0372       } else {
0373         for (unsigned int iDaqPartition = 0; iDaqPartition < m_numberDaqPartitions; ++iDaqPartition) {
0374           m_globalNrErrors[iDaqPartition]++;
0375         }
0376 
0377         edm::LogWarning("L1GtTrigReport")
0378             << "\n  L1GlobalTriggerReadoutRecord with input tag " << m_l1GtRecordInputTag.label()
0379             << " has gtFdlVector of size " << fdlVecSize << "\n  Invalid L1GlobalTriggerReadoutRecord!"
0380             << "\n  Event classified as Error\n\n"
0381             << std::endl;
0382 
0383         validRecord = false;
0384       }
0385 
0386     } else {
0387       for (unsigned int iDaqPartition = 0; iDaqPartition < m_numberDaqPartitions; ++iDaqPartition) {
0388         m_globalNrErrors[iDaqPartition]++;
0389       }
0390 
0391       edm::LogWarning("L1GtTrigReport") << "\n  L1GlobalTriggerReadoutRecord with input tag "
0392                                         << m_l1GtRecordInputTag.label() << " not found."
0393                                         << "\n  Event classified as Error\n\n"
0394                                         << std::endl;
0395     }
0396   }
0397 
0398   // get the prescale factor set used in the actual luminosity segment
0399   const std::vector<int>& prescaleFactorsAlgoTrig = (*m_prescaleFactorsAlgoTrig).at(pfIndexAlgo);
0400 
0401   const std::vector<int>& prescaleFactorsTechTrig = (*m_prescaleFactorsTechTrig).at(pfIndexTech);
0402 
0403   if (validRecord) {
0404     // loop over algorithms and increase the corresponding counters
0405     for (CItAlgo itAlgo = algorithmMap.begin(); itAlgo != algorithmMap.end(); itAlgo++) {
0406       std::string algName = itAlgo->first;
0407       int algBitNumber = (itAlgo->second).algoBitNumber();
0408 
0409       // the result before applying the trigger masks is available
0410       // in both L1GlobalTriggerReadoutRecord or L1GlobalTriggerRecord
0411       bool algResultBeforeMask = gtDecisionWordBeforeMask[algBitNumber];
0412 
0413       int prescaleFactor = prescaleFactorsAlgoTrig.at(algBitNumber);
0414 
0415       for (unsigned int iDaqPartition = 0; iDaqPartition < m_numberDaqPartitions; ++iDaqPartition) {
0416         unsigned int triggerMask = (m_triggerMaskAlgoTrig.at(algBitNumber)) & (1 << iDaqPartition);
0417 
0418         bool algResultAfterMask = false;
0419 
0420         if (m_useL1GlobalTriggerRecord) {
0421           if (iDaqPartition == m_physicsDaqPartition) {
0422             // result available already for physics DAQ partition
0423             // in lite record
0424             algResultAfterMask = gtDecisionWordAfterMask[algBitNumber];
0425           } else {
0426             // apply the masks for other partitions
0427             algResultAfterMask = algResultBeforeMask;
0428 
0429             if (triggerMask) {
0430               algResultAfterMask = false;
0431             }
0432           }
0433         } else {
0434           // apply the masks for L1GlobalTriggerReadoutRecord
0435           algResultAfterMask = algResultBeforeMask;
0436 
0437           if (triggerMask) {
0438             algResultAfterMask = false;
0439           }
0440         }
0441 
0442         L1GtTrigReportEntry* entryRep =
0443             new L1GtTrigReportEntry(menuName, algName, prescaleFactor, triggerMask, iDaqPartition);
0444 
0445         int iCount = 0;
0446 
0447         for (CItEntry itEntry = m_entryList.begin(); itEntry != m_entryList.end(); itEntry++) {
0448           if ((*entryRep) == *(*itEntry)) {
0449             iCount++;
0450             // increase the corresponding counter in the list entry
0451             (*itEntry)->addValidEntry(algResultAfterMask, algResultBeforeMask);
0452           }
0453         }
0454 
0455         if (iCount == 0) {
0456           // if entry not in the list, increase the corresponding counter
0457           // and push the entry in the list
0458           entryRep->addValidEntry(algResultAfterMask, algResultBeforeMask);
0459           m_entryList.push_back(entryRep);
0460         } else {
0461           delete entryRep;
0462         }
0463       }
0464     }
0465 
0466     // loop over technical triggers and increase the corresponding counters
0467     for (CItAlgo itAlgo = technicalTriggerMap.begin(); itAlgo != technicalTriggerMap.end(); itAlgo++) {
0468       //for (unsigned int iTechTrig = 0; iTechTrig < m_numberTechnicalTriggers; ++iTechTrig) {
0469 
0470       std::string ttName = itAlgo->first;
0471       int ttBitNumber = (itAlgo->second).algoBitNumber();
0472       // std::string ttName = std::to_string(iTechTrig);
0473       // int ttBitNumber = iTechTrig;
0474 
0475       // the result before applying the trigger masks is available
0476       // in both L1GlobalTriggerReadoutRecord or L1GlobalTriggerRecord
0477       bool ttResultBeforeMask = technicalTriggerWordBeforeMask[ttBitNumber];
0478 
0479       int prescaleFactor = prescaleFactorsTechTrig.at(ttBitNumber);
0480 
0481       for (unsigned int iDaqPartition = 0; iDaqPartition < m_numberDaqPartitions; ++iDaqPartition) {
0482         unsigned int triggerMask = (m_triggerMaskTechTrig.at(ttBitNumber)) & (1 << iDaqPartition);
0483 
0484         bool ttResultAfterMask = false;
0485 
0486         if (m_useL1GlobalTriggerRecord) {
0487           if (iDaqPartition == m_physicsDaqPartition) {
0488             // result available already for physics DAQ partition
0489             // in lite record
0490             ttResultAfterMask = technicalTriggerWordAfterMask[ttBitNumber];
0491           } else {
0492             // apply the masks for other partitions
0493             ttResultAfterMask = ttResultBeforeMask;
0494 
0495             if (triggerMask) {
0496               ttResultAfterMask = false;
0497             }
0498           }
0499         } else {
0500           // apply the masks for L1GlobalTriggerReadoutRecord
0501           ttResultAfterMask = ttResultBeforeMask;
0502 
0503           if (triggerMask) {
0504             ttResultAfterMask = false;
0505           }
0506         }
0507 
0508         L1GtTrigReportEntry* entryRep =
0509             new L1GtTrigReportEntry(menuName, ttName, prescaleFactor, triggerMask, iDaqPartition);
0510 
0511         int iCount = 0;
0512 
0513         for (CItEntry itEntry = m_entryListTechTrig.begin(); itEntry != m_entryListTechTrig.end(); itEntry++) {
0514           if ((*entryRep) == *(*itEntry)) {
0515             iCount++;
0516             // increase the corresponding counter in the list entry
0517             (*itEntry)->addValidEntry(ttResultAfterMask, ttResultBeforeMask);
0518           }
0519         }
0520 
0521         if (iCount == 0) {
0522           // if entry not in the list, increase the corresponding counter
0523           // and push the entry in the list
0524           entryRep->addValidEntry(ttResultAfterMask, ttResultBeforeMask);
0525           m_entryListTechTrig.push_back(entryRep);
0526         } else {
0527           delete entryRep;
0528         }
0529       }
0530     }
0531 
0532   } else {
0533     // loop over algorithms and increase the error counters
0534     for (CItAlgo itAlgo = algorithmMap.begin(); itAlgo != algorithmMap.end(); itAlgo++) {
0535       std::string algName = itAlgo->first;
0536       int algBitNumber = (itAlgo->second).algoBitNumber();
0537 
0538       int prescaleFactor = prescaleFactorsAlgoTrig.at(algBitNumber);
0539 
0540       for (unsigned int iDaqPartition = 0; iDaqPartition < m_numberDaqPartitions; ++iDaqPartition) {
0541         unsigned int triggerMask = (m_triggerMaskAlgoTrig.at(algBitNumber)) & (1 << iDaqPartition);
0542 
0543         L1GtTrigReportEntry* entryRep =
0544             new L1GtTrigReportEntry(menuName, algName, prescaleFactor, triggerMask, iDaqPartition);
0545 
0546         int iCount = 0;
0547 
0548         for (CItEntry itEntry = m_entryList.begin(); itEntry != m_entryList.end(); itEntry++) {
0549           if ((*entryRep) == *(*itEntry)) {
0550             iCount++;
0551             // increase the corresponding counter in the list entry
0552             (*itEntry)->addErrorEntry();
0553           }
0554         }
0555 
0556         if (iCount == 0) {
0557           // if entry not in the list, increase the corresponding counter
0558           // and push the entry in the list
0559           entryRep->addErrorEntry();
0560           m_entryList.push_back(entryRep);
0561         } else {
0562           delete entryRep;
0563         }
0564       }
0565     }
0566 
0567     // loop over technical triggers and increase the error counters
0568     // FIXME move to names when technical triggers available in menu
0569     //for (CItAlgo itAlgo = technicalTriggerMap.begin(); itAlgo != technicalTriggerMap.end(); itAlgo++) {
0570     for (unsigned int iTechTrig = 0; iTechTrig < m_numberTechnicalTriggers; ++iTechTrig) {
0571       //std::string ttName = itAlgo->first;
0572       //int ttBitNumber = ( itAlgo->second ).algoBitNumber();
0573       std::string ttName = std::to_string(iTechTrig);
0574       int ttBitNumber = iTechTrig;
0575 
0576       int prescaleFactor = prescaleFactorsTechTrig.at(ttBitNumber);
0577 
0578       for (unsigned int iDaqPartition = 0; iDaqPartition < m_numberDaqPartitions; ++iDaqPartition) {
0579         unsigned int triggerMask = (m_triggerMaskTechTrig.at(ttBitNumber)) & (1 << iDaqPartition);
0580 
0581         L1GtTrigReportEntry* entryRep =
0582             new L1GtTrigReportEntry(menuName, ttName, prescaleFactor, triggerMask, iDaqPartition);
0583 
0584         int iCount = 0;
0585 
0586         for (CItEntry itEntry = m_entryListTechTrig.begin(); itEntry != m_entryListTechTrig.end(); itEntry++) {
0587           if ((*entryRep) == *(*itEntry)) {
0588             iCount++;
0589             // increase the corresponding counter in the list entry
0590             (*itEntry)->addErrorEntry();
0591           }
0592         }
0593 
0594         if (iCount == 0) {
0595           // if entry not in the list, increase the corresponding counter
0596           // and push the entry in the list
0597           entryRep->addErrorEntry();
0598           m_entryListTechTrig.push_back(entryRep);
0599         } else {
0600           delete entryRep;
0601         }
0602       }
0603     }
0604   }
0605 }
0606 
0607 // method called once each job just after ending the event loop
0608 void L1GtTrigReport::endJob() {
0609   // define an output stream to print into
0610   // it can then be directed to whatever log level is desired
0611   std::ostringstream myCout;
0612 
0613   myCout << std::dec << std::endl;
0614   myCout << "L1T-Report "
0615          << "----------       Event Summary       ----------\n";
0616   myCout << "L1T-Report "
0617          << "Total number of events processed: " << m_totalEvents << "\n";
0618   myCout << "L1T-Report\n";
0619 
0620   myCout << "\n"
0621          << "   DAQ partition "
0622          << "           Total "
0623          << " Passed[finalOR] "
0624          << "        Rejected "
0625          << "          Errors "
0626          << "\n"
0627          << std::endl;
0628 
0629   for (unsigned int iDaqPartition = 0; iDaqPartition < m_numberDaqPartitions; ++iDaqPartition) {
0630     int rejectedEvents = m_totalEvents - m_globalNrErrors[iDaqPartition] - m_globalNrAccepts[iDaqPartition];
0631 
0632     if (m_useL1GlobalTriggerRecord && (iDaqPartition != m_physicsDaqPartition)) {
0633       continue;
0634     } else {
0635       myCout << std::right << std::setw(16) << iDaqPartition << " " << std::right << std::setw(16) << m_totalEvents
0636              << " " << std::right << std::setw(16) << m_globalNrAccepts[iDaqPartition] << " " << std::right
0637              << std::setw(16) << rejectedEvents << " " << std::right << std::setw(16) << m_globalNrErrors[iDaqPartition]
0638              << std::endl;
0639     }
0640   }
0641 
0642   // get the list of menus for the sample analyzed
0643   //
0644   std::set<std::string> menuList;
0645   typedef std::set<std::string>::const_iterator CItL1Menu;
0646 
0647   for (CItEntry itEntry = m_entryList.begin(); itEntry != m_entryList.end(); itEntry++) {
0648     menuList.insert((*itEntry)->gtTriggerMenuName());
0649   }
0650 
0651   myCout << "\nThe following L1 menus were used for this sample: " << std::endl;
0652   for (CItL1Menu itMenu = menuList.begin(); itMenu != menuList.end(); itMenu++) {
0653     myCout << "  " << (*itMenu) << std::endl;
0654   }
0655   myCout << "\n" << std::endl;
0656 
0657   switch (m_printVerbosity) {
0658     case 0: {
0659       myCout << "\nL1T-Report "
0660              << "---------- L1 Trigger Global Summary - DAQ Partition " << m_physicsDaqPartition << "----------\n\n";
0661 
0662       myCout << "\n\n Number of events written after applying L1 prescale factors"
0663              << " and trigger masks\n"
0664              << " if not explicitly mentioned.\n\n";
0665 
0666       for (CItL1Menu itMenu = menuList.begin(); itMenu != menuList.end(); itMenu++) {
0667         myCout << "\nReport for L1 menu:  " << (*itMenu) << "\n" << std::endl;
0668 
0669         myCout << std::right << std::setw(45) << "Algorithm Key"
0670                << " " << std::right << std::setw(10) << "Passed"
0671                << " " << std::right << std::setw(10) << "Rejected"
0672                << " " << std::right << std::setw(10) << "Error"
0673                << "\n";
0674 
0675         for (CItEntry itEntry = m_entryList.begin(); itEntry != m_entryList.end(); itEntry++) {
0676           if (((*itEntry)->gtDaqPartition() == m_physicsDaqPartition) && ((*itEntry)->gtTriggerMenuName() == *itMenu)) {
0677             myCout << std::right << std::setw(45) << (*itEntry)->gtAlgoName() << " " << std::right << std::setw(10)
0678                    << (*itEntry)->gtNrEventsAccept() << " " << std::right << std::setw(10)
0679                    << (*itEntry)->gtNrEventsReject() << " " << std::right << std::setw(10)
0680                    << (*itEntry)->gtNrEventsError() << "\n";
0681           }
0682         }
0683 
0684         myCout << "\n\n"
0685                << std::right << std::setw(45) << "Technical Trigger Key"
0686                << " " << std::right << std::setw(10) << "Passed"
0687                << " " << std::right << std::setw(10) << "Rejected"
0688                << " " << std::right << std::setw(10) << "Error"
0689                << "\n";
0690 
0691         for (CItEntry itEntry = m_entryListTechTrig.begin(); itEntry != m_entryListTechTrig.end(); itEntry++) {
0692           if (((*itEntry)->gtDaqPartition() == m_physicsDaqPartition) && ((*itEntry)->gtTriggerMenuName() == *itMenu)) {
0693             myCout << std::right << std::setw(45) << (*itEntry)->gtAlgoName() << " " << std::right << std::setw(10)
0694                    << (*itEntry)->gtNrEventsAccept() << " " << std::right << std::setw(10)
0695                    << (*itEntry)->gtNrEventsReject() << " " << std::right << std::setw(10)
0696                    << (*itEntry)->gtNrEventsError() << "\n";
0697           }
0698         }
0699       }
0700 
0701     }
0702 
0703     break;
0704     case 1: {
0705       myCout << "\nL1T-Report "
0706              << "---------- L1 Trigger Global Summary - DAQ Partition " << m_physicsDaqPartition << "----------\n\n";
0707 
0708       myCout << "\n\n Number of events written after applying L1 prescale factors"
0709              << " and trigger masks\n"
0710              << " if not explicitly mentioned.\n\n";
0711 
0712       for (CItL1Menu itMenu = menuList.begin(); itMenu != menuList.end(); itMenu++) {
0713         myCout << "\nReport for L1 menu:  " << (*itMenu) << "\n" << std::endl;
0714         myCout << std::right << std::setw(45) << "Algorithm Key"
0715                << " " << std::right << std::setw(10) << "Prescale"
0716                << " " << std::right << std::setw(5) << "Mask"
0717                << " " << std::right << std::setw(10) << "Passed"
0718                << " " << std::right << std::setw(10) << "Rejected"
0719                << " " << std::right << std::setw(10) << "Error" << std::setw(2) << " "
0720                << "\n";
0721 
0722         for (CItEntry itEntry = m_entryList.begin(); itEntry != m_entryList.end(); itEntry++) {
0723           if (((*itEntry)->gtDaqPartition() == m_physicsDaqPartition) && ((*itEntry)->gtTriggerMenuName() == *itMenu)) {
0724             myCout << std::right << std::setw(45) << (*itEntry)->gtAlgoName() << " " << std::right << std::setw(10)
0725                    << (*itEntry)->gtPrescaleFactor() << "    " << std::right << std::setw(2)  //<< std::setfill('0')
0726                    << std::hex << (*itEntry)->gtTriggerMask()                                 //<< std::setfill(' ')
0727                    << std::dec << " " << std::right << std::setw(10) << (*itEntry)->gtNrEventsAccept() << " "
0728                    << std::right << std::setw(10) << (*itEntry)->gtNrEventsReject() << " " << std::right
0729                    << std::setw(10) << (*itEntry)->gtNrEventsError() << std::setw(2) << " "
0730                    << "\n";
0731           }
0732         }
0733 
0734         myCout << "\n\n"
0735                << std::right << std::setw(45) << "Technical Trigger Key"
0736                << " " << std::right << std::setw(10) << "Prescale"
0737                << " " << std::right << std::setw(5) << "Mask"
0738                << " " << std::right << std::setw(10) << "Passed"
0739                << " " << std::right << std::setw(10) << "Rejected"
0740                << " " << std::right << std::setw(10) << "Error" << std::setw(2) << " "
0741                << "\n";
0742 
0743         for (CItEntry itEntry = m_entryListTechTrig.begin(); itEntry != m_entryListTechTrig.end(); itEntry++) {
0744           if (((*itEntry)->gtDaqPartition() == m_physicsDaqPartition) && ((*itEntry)->gtTriggerMenuName() == *itMenu)) {
0745             myCout << std::right << std::setw(45) << (*itEntry)->gtAlgoName() << " " << std::right << std::setw(10)
0746                    << (*itEntry)->gtPrescaleFactor() << "    " << std::right << std::setw(2)  //<< std::setfill('0')
0747                    << std::hex << (*itEntry)->gtTriggerMask()                                 //<< std::setfill(' ')
0748                    << std::dec << " " << std::right << std::setw(10) << (*itEntry)->gtNrEventsAccept() << " "
0749                    << std::right << std::setw(10) << (*itEntry)->gtNrEventsReject() << " " << std::right
0750                    << std::setw(10) << (*itEntry)->gtNrEventsError() << std::setw(2) << " "
0751                    << "\n";
0752           }
0753         }
0754       }
0755 
0756     }
0757 
0758     break;
0759     case 2: {
0760       for (CItL1Menu itMenu = menuList.begin(); itMenu != menuList.end(); itMenu++) {
0761         myCout << "\nReport for L1 menu:  " << (*itMenu) << "\n" << std::endl;
0762 
0763         myCout << std::right << std::setw(45) << "Algorithm Key"
0764                << " " << std::right << std::setw(10) << "Passed"
0765                << " " << std::right << std::setw(10) << "Rejected"
0766                << " " << std::right << std::setw(10) << "Error"
0767                << "\n";
0768 
0769         for (CItEntry itEntry = m_entryList.begin(); itEntry != m_entryList.end(); itEntry++) {
0770           if (((*itEntry)->gtDaqPartition() == m_physicsDaqPartition) && ((*itEntry)->gtTriggerMenuName() == *itMenu)) {
0771             int nrEventsAccept = (*itEntry)->gtNrEventsAccept();
0772             int nrEventsReject = (*itEntry)->gtNrEventsReject();
0773             int nrEventsError = (*itEntry)->gtNrEventsError();
0774 
0775             myCout << std::right << std::setw(45) << ((*itEntry)->gtAlgoName()) << " " << std::right << std::setw(10)
0776                    << nrEventsAccept << " " << std::right << std::setw(10) << nrEventsReject << " " << std::right
0777                    << std::setw(10) << nrEventsError << "\n";
0778           }
0779         }
0780 
0781         // efficiency and its statistical error
0782 
0783         myCout << "\n\n"
0784                << std::right << std::setw(45) << "Algorithm Key"
0785                << "    " << std::right << std::setw(10) << "Efficiency "
0786                << " " << std::right << std::setw(10) << "Stat error (%)"
0787                << "\n";
0788 
0789         for (CItEntry itEntry = m_entryList.begin(); itEntry != m_entryList.end(); itEntry++) {
0790           if (((*itEntry)->gtDaqPartition() == 0) && ((*itEntry)->gtTriggerMenuName() == *itMenu)) {
0791             int nrEventsAccept = (*itEntry)->gtNrEventsAccept();
0792             int nrEventsReject = (*itEntry)->gtNrEventsReject();
0793             int nrEventsError = (*itEntry)->gtNrEventsError();
0794 
0795             int totalEvents = nrEventsAccept + nrEventsReject + nrEventsError;
0796 
0797             // efficiency and their statistical error
0798             float eff = 0.;
0799             float statErrEff = 0.;
0800 
0801             if (totalEvents != 0) {
0802               eff = static_cast<float>(nrEventsAccept) / static_cast<float>(totalEvents);
0803               statErrEff = sqrt(eff * (1.0 - eff) / static_cast<float>(totalEvents));
0804             }
0805 
0806             myCout << std::right << std::setw(45) << ((*itEntry)->gtAlgoName()) << " " << std::right << std::setw(10)
0807                    << std::fixed << std::setprecision(2) << 100. * eff << " +- " << std::right << std::setw(10)
0808                    << std::setprecision(2) << 100. * statErrEff << "\n";
0809           }
0810         }
0811 
0812         myCout << "\n\n"
0813                << std::right << std::setw(45) << "Technical Trigger Key"
0814                << " " << std::right << std::setw(10) << "Passed"
0815                << " " << std::right << std::setw(10) << "Rejected"
0816                << " " << std::right << std::setw(10) << "Error"
0817                << "\n";
0818 
0819         for (CItEntry itEntry = m_entryListTechTrig.begin(); itEntry != m_entryListTechTrig.end(); itEntry++) {
0820           if (((*itEntry)->gtDaqPartition() == m_physicsDaqPartition) && ((*itEntry)->gtTriggerMenuName() == *itMenu)) {
0821             int nrEventsAccept = (*itEntry)->gtNrEventsAccept();
0822             int nrEventsReject = (*itEntry)->gtNrEventsReject();
0823             int nrEventsError = (*itEntry)->gtNrEventsError();
0824 
0825             myCout << std::right << std::setw(45) << ((*itEntry)->gtAlgoName()) << " " << std::right << std::setw(10)
0826                    << nrEventsAccept << " " << std::right << std::setw(10) << nrEventsReject << " " << std::right
0827                    << std::setw(10) << nrEventsError << "\n";
0828           }
0829         }
0830 
0831         // efficiency and its statistical error
0832 
0833         myCout << "\n\n"
0834                << std::right << std::setw(45) << "Technical Trigger Key"
0835                << "    " << std::right << std::setw(10) << "Efficiency "
0836                << " " << std::right << std::setw(10) << "Stat error (%)"
0837                << "\n";
0838 
0839         for (CItEntry itEntry = m_entryListTechTrig.begin(); itEntry != m_entryListTechTrig.end(); itEntry++) {
0840           if (((*itEntry)->gtDaqPartition() == 0) && ((*itEntry)->gtTriggerMenuName() == *itMenu)) {
0841             int nrEventsAccept = (*itEntry)->gtNrEventsAccept();
0842             int nrEventsReject = (*itEntry)->gtNrEventsReject();
0843             int nrEventsError = (*itEntry)->gtNrEventsError();
0844 
0845             int totalEvents = nrEventsAccept + nrEventsReject + nrEventsError;
0846 
0847             // efficiency and their statistical error
0848             float eff = 0.;
0849             float statErrEff = 0.;
0850 
0851             if (totalEvents != 0) {
0852               eff = static_cast<float>(nrEventsAccept) / static_cast<float>(totalEvents);
0853               statErrEff = sqrt(eff * (1.0 - eff) / static_cast<float>(totalEvents));
0854             }
0855 
0856             myCout << std::right << std::setw(45) << ((*itEntry)->gtAlgoName()) << " " << std::right << std::setw(10)
0857                    << std::fixed << std::setprecision(2) << 100. * eff << " +- " << std::right << std::setw(10)
0858                    << std::setprecision(2) << 100. * statErrEff << "\n";
0859           }
0860         }
0861       }
0862 
0863     } break;
0864 
0865     case 10: {
0866       myCout << "\nL1T-Report "
0867              << "---------- L1 Trigger Global Summary - DAQ Partition " << m_physicsDaqPartition << "----------\n\n";
0868 
0869       for (CItL1Menu itMenu = menuList.begin(); itMenu != menuList.end(); itMenu++) {
0870         myCout << "\nReport for L1 menu:  " << (*itMenu) << "\n" << std::endl;
0871         myCout << std::right << std::setw(45) << "Algorithm Key"
0872                << " " << std::right << std::setw(10) << "Prescale"
0873                << "  " << std::right << std::setw(5) << "Mask"
0874                << " " << std::right << std::setw(25) << "Before Mask"
0875                << " " << std::right << std::setw(30) << "After Mask"
0876                << " " << std::right << std::setw(22) << "Error"
0877                << "\n"
0878                << std::right << std::setw(64) << " " << std::setw(15) << "Passed" << std::right << std::setw(1) << " "
0879                << std::setw(15) << "Rejected" << std::right << std::setw(1) << " " << std::setw(15) << "Passed"
0880                << std::right << std::setw(1) << " " << std::setw(15) << "Rejected"
0881                << "\n";
0882 
0883         for (CItEntry itEntry = m_entryList.begin(); itEntry != m_entryList.end(); itEntry++) {
0884           if (((*itEntry)->gtDaqPartition() == m_physicsDaqPartition) && ((*itEntry)->gtTriggerMenuName() == *itMenu)) {
0885             myCout << std::right << std::setw(45) << (*itEntry)->gtAlgoName() << " " << std::right << std::setw(10)
0886                    << (*itEntry)->gtPrescaleFactor() << " " << std::right << std::setw(5) << " " << std::hex
0887                    << (*itEntry)->gtTriggerMask() << std::dec << " " << std::right << std::setw(15)
0888                    << (*itEntry)->gtNrEventsAcceptBeforeMask() << " " << std::right << std::setw(15)
0889                    << (*itEntry)->gtNrEventsRejectBeforeMask() << " " << std::right << std::setw(15)
0890                    << (*itEntry)->gtNrEventsAccept() << " " << std::right << std::setw(15)
0891                    << (*itEntry)->gtNrEventsReject() << " " << std::right << std::setw(15)
0892                    << (*itEntry)->gtNrEventsError() << "\n";
0893           }
0894         }
0895 
0896         myCout << "\n\n"
0897                << std::right << std::setw(45) << "Technical Trigger Key"
0898                << " " << std::right << std::setw(10) << "Prescale"
0899                << "  " << std::right << std::setw(5) << "Mask"
0900                << " " << std::right << std::setw(25) << "Before Mask"
0901                << " " << std::right << std::setw(30) << "After Mask"
0902                << " " << std::right << std::setw(22) << "Error"
0903                << "\n"
0904                << std::right << std::setw(64) << " " << std::setw(15) << "Passed" << std::right << std::setw(1) << " "
0905                << std::setw(15) << "Rejected" << std::right << std::setw(1) << " " << std::setw(15) << "Passed"
0906                << std::right << std::setw(1) << " " << std::setw(15) << "Rejected"
0907                << "\n";
0908 
0909         for (CItEntry itEntry = m_entryListTechTrig.begin(); itEntry != m_entryListTechTrig.end(); itEntry++) {
0910           if (((*itEntry)->gtDaqPartition() == m_physicsDaqPartition) && ((*itEntry)->gtTriggerMenuName() == *itMenu)) {
0911             myCout << std::right << std::setw(45) << (*itEntry)->gtAlgoName() << " " << std::right << std::setw(10)
0912                    << (*itEntry)->gtPrescaleFactor() << " " << std::right << std::setw(5) << " " << std::hex
0913                    << (*itEntry)->gtTriggerMask() << std::dec << " " << std::right << std::setw(15)
0914                    << (*itEntry)->gtNrEventsAcceptBeforeMask() << " " << std::right << std::setw(15)
0915                    << (*itEntry)->gtNrEventsRejectBeforeMask() << " " << std::right << std::setw(15)
0916                    << (*itEntry)->gtNrEventsAccept() << " " << std::right << std::setw(15)
0917                    << (*itEntry)->gtNrEventsReject() << " " << std::right << std::setw(15)
0918                    << (*itEntry)->gtNrEventsError() << "\n";
0919           }
0920         }
0921       }
0922     }
0923 
0924     break;
0925     case 100: {
0926       for (unsigned int iDaqPartition = 0; iDaqPartition < m_numberDaqPartitions; ++iDaqPartition) {
0927         myCout << "\nL1T-Report "
0928                << "---------- L1 Trigger Global Summary - DAQ Partition " << iDaqPartition << " "
0929                << "----------\n\n";
0930 
0931         myCout << std::right << std::setw(45) << "Algorithm Key"
0932                << " " << std::right << std::setw(10) << "Passed"
0933                << " " << std::right << std::setw(10) << "Rejected"
0934                << " " << std::right << std::setw(10) << "Error" << std::setw(2) << " "
0935                << "\n";
0936 
0937         for (CItEntry itEntry = m_entryList.begin(); itEntry != m_entryList.end(); itEntry++) {
0938           if (((*itEntry)->gtDaqPartition() == 0)) {
0939             myCout << std::right << std::setw(45) << (*itEntry)->gtAlgoName() << " " << std::right << std::setw(10)
0940                    << (*itEntry)->gtNrEventsAccept() << " " << std::right << std::setw(10)
0941                    << (*itEntry)->gtNrEventsReject() << " " << std::right << std::setw(10)
0942                    << (*itEntry)->gtNrEventsError() << std::setw(2) << " "
0943                    << "\n";
0944           }
0945         }
0946 
0947         myCout << "\n\n"
0948                << std::right << std::setw(45) << "Technical Trigger Key"
0949                << " " << std::right << std::setw(10) << "Passed"
0950                << " " << std::right << std::setw(10) << "Rejected"
0951                << " " << std::right << std::setw(10) << "Error" << std::setw(2) << " "
0952                << "\n";
0953 
0954         for (CItEntry itEntry = m_entryListTechTrig.begin(); itEntry != m_entryListTechTrig.end(); itEntry++) {
0955           if ((*itEntry)->gtDaqPartition() == 0) {
0956             myCout << std::right << std::setw(45) << (*itEntry)->gtAlgoName() << " " << std::right << std::setw(10)
0957                    << (*itEntry)->gtNrEventsAccept() << " " << std::right << std::setw(10)
0958                    << (*itEntry)->gtNrEventsReject() << " " << std::right << std::setw(10)
0959                    << (*itEntry)->gtNrEventsError() << std::setw(2) << " " << std::right << std::setw(20)
0960                    << (*itEntry)->gtTriggerMenuName() << "\n";
0961           }
0962         }
0963       }
0964     }
0965 
0966     break;
0967     case 101: {
0968       for (unsigned int iDaqPartition = 0; iDaqPartition < m_numberDaqPartitions; ++iDaqPartition) {
0969         myCout << "\nL1T-Report "
0970                << "---------- L1 Trigger Global Summary - DAQ Partition " << iDaqPartition << " "
0971                << "----------\n\n";
0972 
0973         myCout << std::right << std::setw(45) << "Algorithm Key"
0974                << " " << std::right << std::setw(10) << "Prescale"
0975                << " " << std::right << std::setw(5) << "Mask"
0976                << " " << std::right << std::setw(10) << "Passed"
0977                << " " << std::right << std::setw(10) << "Rejected"
0978                << " " << std::right << std::setw(10) << "Error" << std::setw(2) << " "
0979                << "\n";
0980 
0981         for (CItEntry itEntry = m_entryList.begin(); itEntry != m_entryList.end(); itEntry++) {
0982           if ((*itEntry)->gtDaqPartition() == 0) {
0983             myCout << std::right << std::setw(45) << (*itEntry)->gtAlgoName() << " " << std::right << std::setw(10)
0984                    << (*itEntry)->gtPrescaleFactor() << "   " << std::right << std::setw(2)  //<< std::setfill('0')
0985                    << std::hex << (*itEntry)->gtTriggerMask()                                //<< std::setfill(' ')
0986                    << std::dec << " " << std::right << std::setw(10) << (*itEntry)->gtNrEventsAccept() << " "
0987                    << std::right << std::setw(10) << (*itEntry)->gtNrEventsReject() << " " << std::right
0988                    << std::setw(10) << (*itEntry)->gtNrEventsError() << std::setw(2) << " "
0989                    << "\n";
0990           }
0991         }
0992 
0993         myCout << "\n\n"
0994                << std::right << std::setw(45) << "Technical Trigger Key"
0995                << " " << std::right << std::setw(10) << "Prescale"
0996                << " " << std::right << std::setw(5) << "Mask"
0997                << " " << std::right << std::setw(10) << "Passed"
0998                << " " << std::right << std::setw(10) << "Rejected"
0999                << " " << std::right << std::setw(10) << "Error" << std::setw(2) << " "
1000                << "\n";
1001 
1002         for (CItEntry itEntry = m_entryListTechTrig.begin(); itEntry != m_entryListTechTrig.end(); itEntry++) {
1003           if ((*itEntry)->gtDaqPartition() == 0) {
1004             myCout << std::right << std::setw(45) << (*itEntry)->gtAlgoName() << " " << std::right << std::setw(10)
1005                    << (*itEntry)->gtPrescaleFactor() << "   " << std::right << std::setw(2)  //<< std::setfill('0')
1006                    << std::hex << (*itEntry)->gtTriggerMask()                                //<< std::setfill(' ')
1007                    << std::dec << " " << std::right << std::setw(10) << (*itEntry)->gtNrEventsAccept() << " "
1008                    << std::right << std::setw(10) << (*itEntry)->gtNrEventsReject() << " " << std::right
1009                    << std::setw(10) << (*itEntry)->gtNrEventsError() << std::setw(2) << " "
1010                    << "\n";
1011           }
1012         }
1013       }
1014     }
1015 
1016     break;
1017     default: {
1018       myCout << "\n\nL1GtTrigReport: Error - no print verbosity level = " << m_printVerbosity
1019              << " defined! \nCheck available values in the cfi file."
1020              << "\n";
1021     }
1022 
1023     break;
1024   }
1025 
1026   // TODO for other verbosity levels
1027   // print the trigger menu, the prescale factors and the trigger mask, etc
1028 
1029   myCout << std::endl;
1030   myCout << "L1T-Report end!" << std::endl;
1031   myCout << std::endl;
1032 
1033   switch (m_printOutput) {
1034     case 0: {
1035       std::cout << myCout.str() << std::endl;
1036 
1037     }
1038 
1039     break;
1040     case 1: {
1041       LogTrace("L1GtTrigReport") << myCout.str() << std::endl;
1042 
1043     } break;
1044 
1045     case 2: {
1046       edm::LogVerbatim("L1GtTrigReport") << myCout.str() << std::endl;
1047 
1048     }
1049 
1050     break;
1051     case 3: {
1052       edm::LogInfo("L1GtTrigReport") << myCout.str();
1053 
1054     }
1055 
1056     break;
1057     default: {
1058       std::cout << "\n\n  L1GtTrigReport: Error - no print output = " << m_printOutput
1059                 << " defined! \n  Check available values in the cfi file."
1060                 << "\n"
1061                 << std::endl;
1062 
1063     } break;
1064   }
1065 }