CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_4_5_patch3/src/L1Trigger/GlobalTriggerAnalyzer/src/L1GtTrigReport.cc

Go to the documentation of this file.
00001 
00017 // this class header
00018 #include "L1Trigger/GlobalTriggerAnalyzer/interface/L1GtTrigReport.h"
00019 
00020 // system include files
00021 #include <memory>
00022 
00023 #include <iostream>
00024 #include <iomanip>
00025 
00026 #include<map>
00027 #include<set>
00028 #include <cmath>
00029 #include <string>
00030 
00031 #include "boost/lexical_cast.hpp"
00032 
00033 // user include files
00034 #include "DataFormats/L1GlobalTrigger/interface/L1GlobalTriggerReadoutSetupFwd.h"
00035 #include "DataFormats/L1GlobalTrigger/interface/L1GlobalTriggerReadoutRecord.h"
00036 #include "DataFormats/L1GlobalTrigger/interface/L1GlobalTriggerRecord.h"
00037 
00038 #include "FWCore/Framework/interface/Frameworkfwd.h"
00039 #include "FWCore/Framework/interface/EDAnalyzer.h"
00040 
00041 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00042 
00043 #include "FWCore/Framework/interface/Event.h"
00044 #include "FWCore/Framework/interface/MakerMacros.h"
00045 
00046 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00047 #include "FWCore/Utilities/interface/InputTag.h"
00048 
00049 #include "FWCore/Framework/interface/EventSetup.h"
00050 #include "FWCore/Framework/interface/ESHandle.h"
00051 
00052 #include "CondFormats/L1TObjects/interface/L1GtStableParameters.h"
00053 #include "CondFormats/DataRecord/interface/L1GtStableParametersRcd.h"
00054 
00055 #include "CondFormats/L1TObjects/interface/L1GtPrescaleFactors.h"
00056 #include "CondFormats/DataRecord/interface/L1GtPrescaleFactorsAlgoTrigRcd.h"
00057 #include "CondFormats/DataRecord/interface/L1GtPrescaleFactorsTechTrigRcd.h"
00058 
00059 #include "CondFormats/L1TObjects/interface/L1GtTriggerMask.h"
00060 #include "CondFormats/DataRecord/interface/L1GtTriggerMaskAlgoTrigRcd.h"
00061 #include "CondFormats/DataRecord/interface/L1GtTriggerMaskTechTrigRcd.h"
00062 
00063 #include "CondFormats/DataRecord/interface/L1GtTriggerMaskVetoAlgoTrigRcd.h"
00064 #include "CondFormats/DataRecord/interface/L1GtTriggerMaskVetoTechTrigRcd.h"
00065 
00066 #include "CondFormats/L1TObjects/interface/L1GtTriggerMenu.h"
00067 #include "CondFormats/DataRecord/interface/L1GtTriggerMenuRcd.h"
00068 
00069 #include "L1Trigger/GlobalTriggerAnalyzer/interface/L1GtTrigReportEntry.h"
00070 
00071 // constructor(s)
00072 L1GtTrigReport::L1GtTrigReport(const edm::ParameterSet& pSet) {
00073 
00074     // boolean flag to select the input record
00075     // if true, it will use L1GlobalTriggerRecord
00076     m_useL1GlobalTriggerRecord = pSet.getParameter<bool>("UseL1GlobalTriggerRecord");
00077 
00079     m_l1GtRecordInputTag = pSet.getParameter<edm::InputTag>("L1GtRecordInputTag");
00080 
00081     // print verbosity
00082     m_printVerbosity = pSet.getUntrackedParameter<int>("PrintVerbosity", 2);
00083 
00084     // print output
00085     m_printOutput = pSet.getUntrackedParameter<int>("PrintOutput", 3);
00086 
00087     LogDebug("L1GtTrigReport") << "\n  Use L1GlobalTriggerRecord:   " << m_useL1GlobalTriggerRecord
00088             << "\n   (if false: L1GtTrigReport uses L1GlobalTriggerReadoutRecord.)"
00089             << "\n  Input tag for L1 GT record:  " << m_l1GtRecordInputTag.label() << " \n"
00090             << "\n  Print verbosity level:           " << m_printVerbosity << " \n"
00091             << "\n  Print output:                    " << m_printOutput << " \n" << std::endl;
00092 
00093     // initialize cached IDs
00094 
00095     //
00096     m_l1GtStableParCacheID = 0ULL;
00097 
00098     m_numberPhysTriggers = 0;
00099     m_numberTechnicalTriggers = 0;
00100     m_numberDaqPartitions = 0;
00101     m_numberDaqPartitionsMax = 0;
00102 
00103     //
00104     m_l1GtPfAlgoCacheID = 0ULL;
00105     m_l1GtPfTechCacheID = 0ULL;
00106 
00107     m_l1GtTmAlgoCacheID = 0ULL;
00108     m_l1GtTmTechCacheID = 0ULL;
00109 
00110     m_l1GtTmVetoAlgoCacheID = 0ULL;
00111     m_l1GtTmVetoTechCacheID = 0ULL;
00112 
00113     //
00114     m_l1GtMenuCacheID = 0ULL;
00115 
00116     // initialize global counters
00117 
00118     // number of events processed
00119     m_totalEvents = 0;
00120 
00121     //
00122     m_entryList.clear();
00123     m_entryListTechTrig.clear();
00124 
00125     // set the index of physics DAQ partition TODO input parameter?
00126     m_physicsDaqPartition = 0;
00127 
00128 }
00129 
00130 // destructor
00131 L1GtTrigReport::~L1GtTrigReport() {
00132 
00133     for (ItEntry itEntry = m_entryList.begin(); itEntry != m_entryList.end(); itEntry++) {
00134         if (*itEntry != 0) {
00135             delete *itEntry;
00136             *itEntry = 0;
00137         }
00138     }
00139 
00140     m_entryList.clear();
00141 
00142     for (ItEntry itEntry = m_entryListTechTrig.begin(); itEntry != m_entryListTechTrig.end(); itEntry++) {
00143         if (*itEntry != 0) {
00144             delete *itEntry;
00145             *itEntry = 0;
00146         }
00147     }
00148 
00149     m_entryListTechTrig.clear();
00150 
00151 }
00152 
00153 // member functions
00154 
00155 
00156 // method called once each job just before starting event loop
00157 void L1GtTrigReport::beginJob() {
00158 
00159     // empty
00160 
00161 }
00162 
00163 // analyze each event
00164 void L1GtTrigReport::analyze(const edm::Event& iEvent, const edm::EventSetup& evSetup) {
00165 
00166     // increase the number of processed events
00167     m_totalEvents++;
00168 
00169     // get / update the stable parameters from the EventSetup
00170     // local cache & check on cacheIdentifier
00171 
00172     unsigned long long l1GtStableParCacheID =
00173             evSetup.get<L1GtStableParametersRcd>().cacheIdentifier();
00174 
00175     if (m_l1GtStableParCacheID != l1GtStableParCacheID) {
00176 
00177         edm::ESHandle<L1GtStableParameters> l1GtStablePar;
00178         evSetup.get<L1GtStableParametersRcd>().get(l1GtStablePar);
00179         m_l1GtStablePar = l1GtStablePar.product();
00180 
00181         // number of physics triggers
00182         m_numberPhysTriggers = m_l1GtStablePar->gtNumberPhysTriggers();
00183 
00184         // number of technical triggers
00185         m_numberTechnicalTriggers = m_l1GtStablePar->gtNumberTechnicalTriggers();
00186 
00187         // number of DAQ partitions
00188         m_numberDaqPartitions = 8; // FIXME add it to stable parameters
00189 
00190         if (m_numberDaqPartitionsMax < m_numberDaqPartitions) {
00191 
00192             int numberDaqPartitionsOld = m_numberDaqPartitionsMax;
00193             m_numberDaqPartitionsMax = m_numberDaqPartitions;
00194 
00195             m_globalNrErrors.reserve(m_numberDaqPartitionsMax);
00196             m_globalNrAccepts.reserve(m_numberDaqPartitionsMax);
00197 
00198             for (unsigned int iDaq = numberDaqPartitionsOld; iDaq < m_numberDaqPartitionsMax; ++iDaq) {
00199 
00200                 m_globalNrErrors.push_back(0);
00201                 m_globalNrAccepts.push_back(0);
00202 
00203             }
00204 
00205         }
00206 
00207         //
00208         m_l1GtStableParCacheID = l1GtStableParCacheID;
00209 
00210     }
00211 
00212     // get / update the prescale factors from the EventSetup
00213     // local cache & check on cacheIdentifier
00214 
00215     unsigned long long l1GtPfAlgoCacheID =
00216             evSetup.get<L1GtPrescaleFactorsAlgoTrigRcd>().cacheIdentifier();
00217 
00218     if (m_l1GtPfAlgoCacheID != l1GtPfAlgoCacheID) {
00219 
00220         edm::ESHandle<L1GtPrescaleFactors> l1GtPfAlgo;
00221         evSetup.get<L1GtPrescaleFactorsAlgoTrigRcd>().get(l1GtPfAlgo);
00222         m_l1GtPfAlgo = l1GtPfAlgo.product();
00223 
00224         m_prescaleFactorsAlgoTrig = & ( m_l1GtPfAlgo->gtPrescaleFactors() );
00225 
00226         m_l1GtPfAlgoCacheID = l1GtPfAlgoCacheID;
00227 
00228     }
00229 
00230     unsigned long long l1GtPfTechCacheID =
00231             evSetup.get<L1GtPrescaleFactorsTechTrigRcd>().cacheIdentifier();
00232 
00233     if (m_l1GtPfTechCacheID != l1GtPfTechCacheID) {
00234 
00235         edm::ESHandle<L1GtPrescaleFactors> l1GtPfTech;
00236         evSetup.get<L1GtPrescaleFactorsTechTrigRcd>().get(l1GtPfTech);
00237         m_l1GtPfTech = l1GtPfTech.product();
00238 
00239         m_prescaleFactorsTechTrig = & ( m_l1GtPfTech->gtPrescaleFactors() );
00240 
00241         m_l1GtPfTechCacheID = l1GtPfTechCacheID;
00242 
00243     }
00244 
00245     // get / update the trigger mask from the EventSetup
00246     // local cache & check on cacheIdentifier
00247 
00248     unsigned long long l1GtTmAlgoCacheID =
00249             evSetup.get<L1GtTriggerMaskAlgoTrigRcd>().cacheIdentifier();
00250 
00251     if (m_l1GtTmAlgoCacheID != l1GtTmAlgoCacheID) {
00252 
00253         edm::ESHandle<L1GtTriggerMask> l1GtTmAlgo;
00254         evSetup.get<L1GtTriggerMaskAlgoTrigRcd>().get(l1GtTmAlgo);
00255         m_l1GtTmAlgo = l1GtTmAlgo.product();
00256 
00257         m_triggerMaskAlgoTrig = m_l1GtTmAlgo->gtTriggerMask();
00258 
00259         m_l1GtTmAlgoCacheID = l1GtTmAlgoCacheID;
00260 
00261     }
00262 
00263     unsigned long long l1GtTmTechCacheID =
00264             evSetup.get<L1GtTriggerMaskTechTrigRcd>().cacheIdentifier();
00265 
00266     if (m_l1GtTmTechCacheID != l1GtTmTechCacheID) {
00267 
00268         edm::ESHandle<L1GtTriggerMask> l1GtTmTech;
00269         evSetup.get<L1GtTriggerMaskTechTrigRcd>().get(l1GtTmTech);
00270         m_l1GtTmTech = l1GtTmTech.product();
00271 
00272         m_triggerMaskTechTrig = m_l1GtTmTech->gtTriggerMask();
00273 
00274         m_l1GtTmTechCacheID = l1GtTmTechCacheID;
00275 
00276     }
00277 
00278     unsigned long long l1GtTmVetoAlgoCacheID =
00279             evSetup.get<L1GtTriggerMaskVetoAlgoTrigRcd>().cacheIdentifier();
00280 
00281     if (m_l1GtTmVetoAlgoCacheID != l1GtTmVetoAlgoCacheID) {
00282 
00283         edm::ESHandle<L1GtTriggerMask> l1GtTmVetoAlgo;
00284         evSetup.get<L1GtTriggerMaskVetoAlgoTrigRcd>().get(l1GtTmVetoAlgo);
00285         m_l1GtTmVetoAlgo = l1GtTmVetoAlgo.product();
00286 
00287         m_triggerMaskVetoAlgoTrig = m_l1GtTmVetoAlgo->gtTriggerMask();
00288 
00289         m_l1GtTmVetoAlgoCacheID = l1GtTmVetoAlgoCacheID;
00290 
00291     }
00292 
00293     unsigned long long l1GtTmVetoTechCacheID =
00294             evSetup.get<L1GtTriggerMaskVetoTechTrigRcd>().cacheIdentifier();
00295 
00296     if (m_l1GtTmVetoTechCacheID != l1GtTmVetoTechCacheID) {
00297 
00298         edm::ESHandle<L1GtTriggerMask> l1GtTmVetoTech;
00299         evSetup.get<L1GtTriggerMaskVetoTechTrigRcd>().get(l1GtTmVetoTech);
00300         m_l1GtTmVetoTech = l1GtTmVetoTech.product();
00301 
00302         m_triggerMaskVetoTechTrig = m_l1GtTmVetoTech->gtTriggerMask();
00303 
00304         m_l1GtTmVetoTechCacheID = l1GtTmVetoTechCacheID;
00305 
00306     }
00307 
00308     // get / update the trigger menu from the EventSetup
00309     // local cache & check on cacheIdentifier
00310 
00311     unsigned long long l1GtMenuCacheID =
00312             evSetup.get<L1GtTriggerMenuRcd>().cacheIdentifier();
00313 
00314     if (m_l1GtMenuCacheID != l1GtMenuCacheID) {
00315 
00316         edm::ESHandle<L1GtTriggerMenu> l1GtMenu;
00317         evSetup.get<L1GtTriggerMenuRcd>().get(l1GtMenu);
00318         m_l1GtMenu = l1GtMenu.product();
00319 
00320         m_l1GtMenuCacheID = l1GtMenuCacheID;
00321 
00322         LogDebug("L1GtTrigReport") << "\n  Changing L1 menu to : \n"
00323                 << m_l1GtMenu->gtTriggerMenuName() << "\n\n" << std::endl;
00324 
00325     }
00326 
00327 
00328     const AlgorithmMap& algorithmMap = m_l1GtMenu->gtAlgorithmMap();
00329     const AlgorithmMap& technicalTriggerMap = m_l1GtMenu->gtTechnicalTriggerMap();
00330 
00331     const std::string& menuName = m_l1GtMenu->gtTriggerMenuName();
00332 
00333     // ... end EventSetup
00334 
00335     // get L1GlobalTriggerReadoutRecord or L1GlobalTriggerRecord
00336     // in L1GlobalTriggerRecord, only the physics partition is available
00337     edm::Handle<L1GlobalTriggerReadoutRecord> gtReadoutRecord;
00338     edm::Handle<L1GlobalTriggerRecord> gtRecord;
00339 
00340     if (m_useL1GlobalTriggerRecord) {
00341         iEvent.getByLabel(m_l1GtRecordInputTag, gtRecord);
00342     } else {
00343         iEvent.getByLabel(m_l1GtRecordInputTag, gtReadoutRecord);
00344     }
00345 
00346     bool validRecord = false;
00347 
00348     unsigned int pfIndexAlgo = 0; // get them later from the record
00349     unsigned int pfIndexTech = 0;
00350 
00351     DecisionWord gtDecisionWordBeforeMask;
00352     DecisionWord gtDecisionWordAfterMask;
00353 
00354     TechnicalTriggerWord technicalTriggerWordBeforeMask;
00355     TechnicalTriggerWord technicalTriggerWordAfterMask;
00356 
00357     if (m_useL1GlobalTriggerRecord) {
00358 
00359         if (gtRecord.isValid()) {
00360 
00361             // get Global Trigger decision and the decision word
00362             bool gtDecision = gtRecord->decision();
00363 
00364             gtDecisionWordBeforeMask = gtRecord->decisionWordBeforeMask();
00365             gtDecisionWordAfterMask = gtRecord->decisionWord();
00366 
00367             technicalTriggerWordBeforeMask = gtRecord->technicalTriggerWordBeforeMask();
00368             technicalTriggerWordAfterMask = gtRecord->technicalTriggerWord();
00369 
00370             if (gtDecision) {
00371                 m_globalNrAccepts[m_physicsDaqPartition]++;
00372             }
00373 
00374             pfIndexAlgo = gtRecord->gtPrescaleFactorIndexAlgo();
00375             pfIndexTech = gtRecord->gtPrescaleFactorIndexTech();
00376 
00377             validRecord = true;
00378 
00379         } else {
00380 
00381             m_globalNrErrors[m_physicsDaqPartition]++;
00382 
00383             edm::LogWarning("L1GtTrigReport") << "\n  L1GlobalTriggerRecord with input tag "
00384                     << m_l1GtRecordInputTag.label() << " not found."
00385                     << "\n  Event classified as Error\n\n"
00386                     << std::endl;
00387 
00388         }
00389 
00390     } else {
00391         if (gtReadoutRecord.isValid()) {
00392 
00393             // check if the readout record has size greater than zero, then proceeds
00394             const std::vector<L1GtFdlWord>& fdlVec = gtReadoutRecord->gtFdlVector();
00395             size_t fdlVecSize = fdlVec.size();
00396 
00397             if (fdlVecSize > 0) {
00398 
00399                 LogDebug("L1GtTrigReport") << "\n  L1GlobalTriggerReadoutRecord with input tag "
00400                         << m_l1GtRecordInputTag.label() << " has gtFdlVector of size " << fdlVecSize
00401                         << std::endl;
00402 
00403                 // get Global Trigger finalOR and the decision word
00404                 boost::uint16_t gtFinalOR = gtReadoutRecord->finalOR();
00405 
00406                 gtDecisionWordBeforeMask = gtReadoutRecord->decisionWord();
00407                 technicalTriggerWordBeforeMask = gtReadoutRecord->technicalTriggerWord();
00408 
00409                 for (unsigned int iDaqPartition = 0; iDaqPartition < m_numberDaqPartitions; ++iDaqPartition) {
00410 
00411                     bool gtDecision = static_cast<bool>(gtFinalOR & ( 1 << iDaqPartition ));
00412                     if (gtDecision) {
00413                         m_globalNrAccepts[iDaqPartition]++;
00414                     }
00415 
00416                 }
00417 
00418                 pfIndexAlgo
00419                         = static_cast<unsigned int>( ( gtReadoutRecord->gtFdlWord() ).gtPrescaleFactorIndexAlgo());
00420                 pfIndexTech
00421                         = static_cast<unsigned int>( ( gtReadoutRecord->gtFdlWord() ).gtPrescaleFactorIndexTech());
00422 
00423                 validRecord = true;
00424 
00425             } else {
00426 
00427                 for (unsigned int iDaqPartition = 0; iDaqPartition < m_numberDaqPartitions; ++iDaqPartition) {
00428                     m_globalNrErrors[iDaqPartition]++;
00429                 }
00430 
00431                 edm::LogWarning("L1GtTrigReport") << "\n  L1GlobalTriggerReadoutRecord with input tag "
00432                         << m_l1GtRecordInputTag.label() << " has gtFdlVector of size " << fdlVecSize
00433                         << "\n  Invalid L1GlobalTriggerReadoutRecord!"
00434                         << "\n  Event classified as Error\n\n"
00435                         << std::endl;
00436 
00437                 validRecord = false;
00438 
00439             }
00440 
00441         } else {
00442 
00443             for (unsigned int iDaqPartition = 0; iDaqPartition < m_numberDaqPartitions; ++iDaqPartition) {
00444                 m_globalNrErrors[iDaqPartition]++;
00445             }
00446 
00447             edm::LogWarning("L1GtTrigReport") << "\n  L1GlobalTriggerReadoutRecord with input tag "
00448                     << m_l1GtRecordInputTag.label() << " not found."
00449                     << "\n  Event classified as Error\n\n"
00450                     << std::endl;
00451 
00452         }
00453 
00454     }
00455 
00456     // get the prescale factor set used in the actual luminosity segment
00457     const std::vector<int>& prescaleFactorsAlgoTrig =
00458             ( *m_prescaleFactorsAlgoTrig ).at(pfIndexAlgo);
00459 
00460     const std::vector<int>& prescaleFactorsTechTrig =
00461             ( *m_prescaleFactorsTechTrig ).at(pfIndexTech);
00462 
00463 
00464     if (validRecord) {
00465 
00466         // loop over algorithms and increase the corresponding counters
00467         for (CItAlgo itAlgo = algorithmMap.begin(); itAlgo != algorithmMap.end(); itAlgo++) {
00468 
00469             std::string algName = itAlgo->first;
00470             int algBitNumber = ( itAlgo->second ).algoBitNumber();
00471 
00472             // the result before applying the trigger masks is available
00473             // in both L1GlobalTriggerReadoutRecord or L1GlobalTriggerRecord
00474             bool algResultBeforeMask = gtDecisionWordBeforeMask[algBitNumber];
00475 
00476             int prescaleFactor = prescaleFactorsAlgoTrig.at(algBitNumber);
00477 
00478             for (unsigned int iDaqPartition = 0; iDaqPartition < m_numberDaqPartitions; ++iDaqPartition) {
00479 
00480                 unsigned int triggerMask = ( m_triggerMaskAlgoTrig.at(algBitNumber) ) & ( 1
00481                         << iDaqPartition );
00482 
00483                 bool algResultAfterMask = false;
00484 
00485                 if (m_useL1GlobalTriggerRecord) {
00486                     if (iDaqPartition == m_physicsDaqPartition) {
00487                         // result available already for physics DAQ partition
00488                         // in lite record
00489                         algResultAfterMask = gtDecisionWordAfterMask[algBitNumber];
00490                     } else {
00491                         // apply the masks for other partitions
00492                         algResultAfterMask = algResultBeforeMask;
00493 
00494                         if (triggerMask) {
00495                             algResultAfterMask = false;
00496                         }
00497                     }
00498                 } else {
00499                     // apply the masks for L1GlobalTriggerReadoutRecord
00500                     algResultAfterMask = algResultBeforeMask;
00501 
00502                     if (triggerMask) {
00503                         algResultAfterMask = false;
00504                     }
00505                 }
00506 
00507                 L1GtTrigReportEntry* entryRep = new L1GtTrigReportEntry(
00508                         menuName, algName, prescaleFactor, triggerMask, iDaqPartition);
00509 
00510                 int iCount = 0;
00511 
00512                 for (CItEntry itEntry = m_entryList.begin(); itEntry != m_entryList.end(); itEntry++) {
00513                     if ( ( *entryRep ) == * ( *itEntry )) {
00514                         iCount++;
00515                         // increase the corresponding counter in the list entry
00516                         ( *itEntry )->addValidEntry(algResultAfterMask, algResultBeforeMask);
00517                     }
00518                 }
00519 
00520                 if (iCount == 0) {
00521                     // if entry not in the list, increase the corresponding counter
00522                     // and push the entry in the list
00523                     entryRep->addValidEntry(algResultAfterMask, algResultBeforeMask);
00524                     m_entryList.push_back(entryRep);
00525                 } else {
00526                     delete entryRep;
00527                 }
00528             }
00529         }
00530 
00531         // loop over technical triggers and increase the corresponding counters
00532         for (CItAlgo itAlgo = technicalTriggerMap.begin(); itAlgo != technicalTriggerMap.end(); itAlgo++) {
00533         //for (unsigned int iTechTrig = 0; iTechTrig < m_numberTechnicalTriggers; ++iTechTrig) {
00534 
00535             std::string ttName = itAlgo->first;
00536             int ttBitNumber = ( itAlgo->second ).algoBitNumber();
00537             // std::string ttName = boost::lexical_cast<std::string>(iTechTrig);
00538             // int ttBitNumber = iTechTrig;
00539 
00540             // the result before applying the trigger masks is available
00541             // in both L1GlobalTriggerReadoutRecord or L1GlobalTriggerRecord
00542             bool ttResultBeforeMask = technicalTriggerWordBeforeMask[ttBitNumber];
00543 
00544             int prescaleFactor = prescaleFactorsTechTrig.at(ttBitNumber);
00545 
00546             for (unsigned int iDaqPartition = 0; iDaqPartition < m_numberDaqPartitions; ++iDaqPartition) {
00547 
00548                 unsigned int triggerMask = ( m_triggerMaskTechTrig.at(ttBitNumber) ) & ( 1
00549                         << iDaqPartition );
00550 
00551                 bool ttResultAfterMask = false;
00552 
00553                 if (m_useL1GlobalTriggerRecord) {
00554                     if (iDaqPartition == m_physicsDaqPartition) {
00555                         // result available already for physics DAQ partition
00556                         // in lite record
00557                         ttResultAfterMask = technicalTriggerWordAfterMask[ttBitNumber];
00558                     } else {
00559                         // apply the masks for other partitions
00560                         ttResultAfterMask = ttResultBeforeMask;
00561 
00562                         if (triggerMask) {
00563                             ttResultAfterMask = false;
00564                         }
00565                     }
00566                 } else {
00567                     // apply the masks for L1GlobalTriggerReadoutRecord
00568                     ttResultAfterMask = ttResultBeforeMask;
00569 
00570                     if (triggerMask) {
00571                         ttResultAfterMask = false;
00572                     }
00573                 }
00574 
00575                 L1GtTrigReportEntry* entryRep = new L1GtTrigReportEntry(
00576                         menuName, ttName, prescaleFactor, triggerMask, iDaqPartition);
00577 
00578                 int iCount = 0;
00579 
00580                 for (CItEntry itEntry = m_entryListTechTrig.begin(); itEntry
00581                         != m_entryListTechTrig.end(); itEntry++) {
00582                     if ( ( *entryRep ) == * ( *itEntry )) {
00583                         iCount++;
00584                         // increase the corresponding counter in the list entry
00585                         ( *itEntry )->addValidEntry(ttResultAfterMask, ttResultBeforeMask);
00586                     }
00587                 }
00588 
00589                 if (iCount == 0) {
00590                     // if entry not in the list, increase the corresponding counter
00591                     // and push the entry in the list
00592                     entryRep->addValidEntry(ttResultAfterMask, ttResultBeforeMask);
00593                     m_entryListTechTrig.push_back(entryRep);
00594                 } else {
00595                     delete entryRep;
00596                 }
00597             }
00598         }
00599 
00600     } else {
00601 
00602         // loop over algorithms and increase the error counters
00603         for (CItAlgo itAlgo = algorithmMap.begin(); itAlgo != algorithmMap.end(); itAlgo++) {
00604 
00605             std::string algName = itAlgo->first;
00606             int algBitNumber = ( itAlgo->second ).algoBitNumber();
00607 
00608             int prescaleFactor = prescaleFactorsAlgoTrig.at(algBitNumber);
00609 
00610             for (unsigned int iDaqPartition = 0; iDaqPartition < m_numberDaqPartitions; ++iDaqPartition) {
00611 
00612                 unsigned int triggerMask = ( m_triggerMaskAlgoTrig.at(algBitNumber) ) & ( 1
00613                         << iDaqPartition );
00614 
00615                 L1GtTrigReportEntry* entryRep = new L1GtTrigReportEntry(
00616                         menuName, algName, prescaleFactor, triggerMask, iDaqPartition);
00617 
00618                 int iCount = 0;
00619 
00620                 for (CItEntry itEntry = m_entryList.begin(); itEntry != m_entryList.end(); itEntry++) {
00621 
00622                     if ( ( *entryRep ) == * ( *itEntry )) {
00623                         iCount++;
00624                         // increase the corresponding counter in the list entry
00625                         ( *itEntry )->addErrorEntry();
00626                     }
00627                 }
00628 
00629                 if (iCount == 0) {
00630                     // if entry not in the list, increase the corresponding counter
00631                     // and push the entry in the list
00632                     entryRep->addErrorEntry();
00633                     m_entryList.push_back(entryRep);
00634                 } else {
00635                     delete entryRep;
00636                 }
00637             }
00638 
00639         }
00640 
00641         // loop over technical triggers and increase the error counters
00642         // FIXME move to names when technical triggers available in menu
00643         //for (CItAlgo itAlgo = technicalTriggerMap.begin(); itAlgo != technicalTriggerMap.end(); itAlgo++) {
00644         for (unsigned int iTechTrig = 0; iTechTrig < m_numberTechnicalTriggers; ++iTechTrig) {
00645 
00646             //std::string ttName = itAlgo->first;
00647             //int ttBitNumber = ( itAlgo->second ).algoBitNumber();
00648             std::string ttName = boost::lexical_cast<std::string>(iTechTrig);
00649             int ttBitNumber = iTechTrig;
00650 
00651             int prescaleFactor = prescaleFactorsTechTrig.at(ttBitNumber);
00652 
00653             for (unsigned int iDaqPartition = 0; iDaqPartition < m_numberDaqPartitions; ++iDaqPartition) {
00654 
00655                 unsigned int triggerMask = ( m_triggerMaskTechTrig.at(ttBitNumber) ) & ( 1
00656                         << iDaqPartition );
00657 
00658                 L1GtTrigReportEntry* entryRep = new L1GtTrigReportEntry(
00659                         menuName, ttName, prescaleFactor, triggerMask, iDaqPartition);
00660 
00661                 int iCount = 0;
00662 
00663                 for (CItEntry itEntry = m_entryListTechTrig.begin(); itEntry
00664                         != m_entryListTechTrig.end(); itEntry++) {
00665 
00666                     if ( ( *entryRep ) == * ( *itEntry )) {
00667                         iCount++;
00668                         // increase the corresponding counter in the list entry
00669                         ( *itEntry )->addErrorEntry();
00670                     }
00671                 }
00672 
00673                 if (iCount == 0) {
00674                     // if entry not in the list, increase the corresponding counter
00675                     // and push the entry in the list
00676                     entryRep->addErrorEntry();
00677                     m_entryListTechTrig.push_back(entryRep);
00678                 } else {
00679                     delete entryRep;
00680                 }
00681             }
00682 
00683         }
00684 
00685     }
00686 
00687 }
00688 
00689 // method called once each job just after ending the event loop
00690 void L1GtTrigReport::endJob() {
00691 
00692     // define an output stream to print into
00693     // it can then be directed to whatever log level is desired
00694     std::ostringstream myCout;
00695 
00696     myCout << std::dec << std::endl;
00697     myCout << "L1T-Report " << "----------       Event Summary       ----------\n";
00698     myCout << "L1T-Report " << "Total number of events processed: " << m_totalEvents << "\n";
00699     myCout << "L1T-Report\n";
00700 
00701     myCout
00702         << "\n"
00703         << "   DAQ partition "
00704         << "           Total "
00705         << " Passed[finalOR] "
00706         << "        Rejected "
00707         << "          Errors "
00708         << "\n" << std::endl;
00709 
00710     for (unsigned int iDaqPartition = 0; iDaqPartition < m_numberDaqPartitions; ++iDaqPartition) {
00711 
00712         int rejectedEvents = m_totalEvents - m_globalNrErrors[iDaqPartition]
00713                 - m_globalNrAccepts[iDaqPartition];
00714 
00715         if (m_useL1GlobalTriggerRecord && ( iDaqPartition != m_physicsDaqPartition )) {
00716             continue;
00717         } else {
00718 
00719             myCout
00720                 << std::right << std::setw(16) << iDaqPartition << " "
00721                 << std::right << std::setw(16) << m_totalEvents << " "
00722                 << std::right << std::setw(16) << m_globalNrAccepts[iDaqPartition] << " "
00723                 << std::right << std::setw(16) << rejectedEvents << " "
00724                 << std::right << std::setw(16) << m_globalNrErrors[iDaqPartition] << std::endl;
00725 
00726         }
00727 
00728     }
00729 
00730     // get the list of menus for the sample analyzed
00731     //
00732     std::set<std::string> menuList;
00733     typedef std::set<std::string>::const_iterator CItL1Menu;
00734 
00735     for (CItEntry itEntry = m_entryList.begin(); itEntry != m_entryList.end(); itEntry++) {
00736         menuList.insert( ( *itEntry )->gtTriggerMenuName());
00737     }
00738 
00739     myCout
00740             << "\nThe following L1 menus were used for this sample: " << std::endl;
00741     for (CItL1Menu itMenu = menuList.begin(); itMenu != menuList.end(); itMenu++) {
00742         myCout << "  " << ( *itMenu ) << std::endl;
00743     }
00744     myCout << "\n" << std::endl;
00745 
00746     switch (m_printVerbosity) {
00747         case 0: {
00748 
00749             myCout
00750                 << "\nL1T-Report " << "---------- L1 Trigger Global Summary - DAQ Partition "
00751                 << m_physicsDaqPartition << "----------\n\n";
00752 
00753             myCout
00754                 << "\n\n Number of events written after applying L1 prescale factors"
00755                 << " and trigger masks\n" << " if not explicitly mentioned.\n\n";
00756 
00757             for (CItL1Menu itMenu = menuList.begin(); itMenu != menuList.end(); itMenu++) {
00758 
00759                 myCout << "Report for L1 menu:  " << (*itMenu) << "\n"
00760                         << std::endl;
00761 
00762                 myCout
00763                     << std::right << std::setw(35) << "Algorithm Key" << " "
00764                     << std::right << std::setw(10) << "Passed" << " "
00765                     << std::right << std::setw(10) << "Rejected" << " "
00766                     << std::right << std::setw(10) << "Error"
00767                     << "\n";
00768 
00769                 for (CItEntry itEntry = m_entryList.begin(); itEntry
00770                         != m_entryList.end(); itEntry++) {
00771 
00772                     if ((*itEntry)->gtDaqPartition() == m_physicsDaqPartition) {
00773 
00774                         myCout
00775                             << std::right << std::setw(35) << (*itEntry)->gtAlgoName() << " "
00776                             << std::right << std::setw(10) << (*itEntry)->gtNrEventsAccept() << " "
00777                             << std::right << std::setw(10) << (*itEntry)->gtNrEventsReject() << " "
00778                             << std::right << std::setw(10) << (*itEntry)->gtNrEventsError()
00779                             << "\n";
00780                     }
00781 
00782                 }
00783 
00784                 myCout
00785                     << "\n\n"
00786                     << std::right << std::setw(35) << "Technical Trigger Key" << " "
00787                     << std::right << std::setw(10) << "Passed" << " "
00788                     << std::right << std::setw(10) << "Rejected" << " "
00789                     << std::right << std::setw(10) << "Error"
00790                     << "\n";
00791 
00792                 for (CItEntry itEntry = m_entryListTechTrig.begin(); itEntry
00793                         != m_entryListTechTrig.end(); itEntry++) {
00794 
00795                     if ((*itEntry)->gtDaqPartition() == m_physicsDaqPartition) {
00796 
00797                         myCout
00798                             << std::right << std::setw(35) << (*itEntry)->gtAlgoName() << " "
00799                             << std::right << std::setw(10) << (*itEntry)->gtNrEventsAccept() << " "
00800                             << std::right << std::setw(10) << (*itEntry)->gtNrEventsReject() << " "
00801                             << std::right << std::setw(10) << (*itEntry)->gtNrEventsError()
00802                             << "\n";
00803                     }
00804 
00805                 }
00806             }
00807 
00808         }
00809 
00810             break;
00811         case 1: {
00812 
00813             myCout << "\nL1T-Report " << "---------- L1 Trigger Global Summary - DAQ Partition "
00814                     << m_physicsDaqPartition << "----------\n\n";
00815 
00816             myCout << "\n\n Number of events written after applying L1 prescale factors"
00817                     << " and trigger masks\n" << " if not explicitly mentioned.\n\n";
00818 
00819             for (CItL1Menu itMenu = menuList.begin(); itMenu != menuList.end(); itMenu++) {
00820 
00821                 myCout << "Report for L1 menu:  " << (*itMenu) << "\n"
00822                         << std::endl;
00823                 myCout
00824                     << std::right << std::setw(35) << "Algorithm Key" << " "
00825                     << std::right << std::setw(10) << "Prescale" << " "
00826                     << std::right << std::setw(5)  << "Mask" << " "
00827                     << std::right << std::setw(10) << "Passed" << " "
00828                     << std::right << std::setw(10) << "Rejected" << " "
00829                     << std::right << std::setw(10) << "Error" << std::setw(2) << " "
00830                     << "\n";
00831 
00832                 for (CItEntry itEntry = m_entryList.begin(); itEntry != m_entryList.end(); itEntry++) {
00833 
00834                     if ( ( *itEntry )->gtDaqPartition() == m_physicsDaqPartition) {
00835                         myCout
00836                             << std::right << std::setw(35) << ( *itEntry )->gtAlgoName() << " "
00837                             << std::right << std::setw(10) << ( *itEntry )->gtPrescaleFactor() << "    "
00838                             << std::right << std::setw(2) //<< std::setfill('0')
00839                             << std::hex << ( *itEntry )->gtTriggerMask() //<< std::setfill(' ')
00840                             << std::dec << " "
00841                             << std::right << std::setw(10) << ( *itEntry )->gtNrEventsAccept() << " "
00842                             << std::right << std::setw(10) << ( *itEntry )->gtNrEventsReject() << " "
00843                             << std::right << std::setw(10) << ( *itEntry )->gtNrEventsError() << std::setw(2) << " "
00844                             << "\n";
00845                     }
00846                 }
00847 
00848                 myCout
00849                     << "\n\n"
00850                     << std::right << std::setw(35) << "Technical Trigger Key" << " "
00851                     << std::right << std::setw(10) << "Prescale" << " "
00852                     << std::right << std::setw(5)  << "Mask" << " "
00853                     << std::right << std::setw(10) << "Passed" << " "
00854                     << std::right << std::setw(10) << "Rejected" << " "
00855                     << std::right << std::setw(10) << "Error" << std::setw(2) << " "
00856                     << "\n";
00857 
00858                 for (CItEntry itEntry = m_entryListTechTrig.begin(); itEntry != m_entryListTechTrig.end(); itEntry++) {
00859 
00860                     if ( ( *itEntry )->gtDaqPartition() == m_physicsDaqPartition) {
00861                         myCout
00862                             << std::right << std::setw(35) << ( *itEntry )->gtAlgoName() << " "
00863                             << std::right << std::setw(10) << ( *itEntry )->gtPrescaleFactor() << "    "
00864                             << std::right << std::setw(2) //<< std::setfill('0')
00865                             << std::hex << ( *itEntry )->gtTriggerMask() //<< std::setfill(' ')
00866                             << std::dec << " "
00867                             << std::right << std::setw(10) << ( *itEntry )->gtNrEventsAccept() << " "
00868                             << std::right << std::setw(10) << ( *itEntry )->gtNrEventsReject() << " "
00869                             << std::right << std::setw(10) << ( *itEntry )->gtNrEventsError() << std::setw(2) << " "
00870                             << "\n";
00871                     }
00872                 }
00873             }
00874 
00875         }
00876 
00877             break;
00878         case 2: {
00879 
00880 
00881             for (CItL1Menu itMenu = menuList.begin(); itMenu != menuList.end(); itMenu++) {
00882 
00883                 myCout << "Report for L1 menu:  " << ( *itMenu ) << "\n"
00884                         << std::endl;
00885 
00886                 myCout
00887                     << std::right << std::setw(35) << "Algorithm Key" << " "
00888                     << std::right << std::setw(10) << "Passed" << " "
00889                     << std::right << std::setw(10) << "Rejected" << " "
00890                     << std::right << std::setw(10) << "Error" << "\n";
00891 
00892                 for (CItEntry itEntry = m_entryList.begin(); itEntry != m_entryList.end(); itEntry++) {
00893 
00894                     if ( ( ( *itEntry )->gtDaqPartition() == m_physicsDaqPartition )
00895                             && ( ( *itEntry )->gtTriggerMenuName() == *itMenu )) {
00896 
00897                         int nrEventsAccept = ( *itEntry )->gtNrEventsAccept();
00898                         int nrEventsReject = ( *itEntry )->gtNrEventsReject();
00899                         int nrEventsError = ( *itEntry )->gtNrEventsError();
00900 
00901                         myCout
00902                             << std::right << std::setw(35) << (( *itEntry )->gtAlgoName()) << " "
00903                             << std::right << std::setw(10) << nrEventsAccept << " "
00904                             << std::right << std::setw(10) << nrEventsReject << " "
00905                             << std::right << std::setw(10) << nrEventsError << "\n";
00906 
00907                     }
00908                 }
00909 
00910                 // efficiency and its statistical error
00911 
00912                 myCout << "\n\n"
00913                     << std::right << std::setw(35) << "Algorithm Key" << "    "
00914                     << std::right << std::setw(10) << "Efficiency " << " "
00915                     << std::right << std::setw(10) << "Stat error (%)" << "\n";
00916 
00917                 for (CItEntry itEntry = m_entryList.begin(); itEntry != m_entryList.end(); itEntry++) {
00918 
00919                     if ( ( ( *itEntry )->gtDaqPartition() == 0 ) && ( ( *itEntry )->gtTriggerMenuName()
00920                             == *itMenu )) {
00921 
00922                         int nrEventsAccept = ( *itEntry )->gtNrEventsAccept();
00923                         int nrEventsReject = ( *itEntry )->gtNrEventsReject();
00924                         int nrEventsError = ( *itEntry )->gtNrEventsError();
00925 
00926                         int totalEvents = nrEventsAccept + nrEventsReject + nrEventsError;
00927 
00928                         // efficiency and their statistical error
00929                         float eff = 0.;
00930                         float statErrEff = 0.;
00931 
00932                         if (totalEvents != 0) {
00933                             eff = static_cast<float> (nrEventsAccept)
00934                                     / static_cast<float> (totalEvents);
00935                             statErrEff = sqrt(eff * ( 1.0 - eff )
00936                                     / static_cast<float> (totalEvents));
00937 
00938                         }
00939 
00940                         myCout
00941                             << std::right << std::setw(35) << (( *itEntry )->gtAlgoName()) << " "
00942                             << std::right << std::setw(10) << std::fixed << std::setprecision(2)
00943                             << 100.*eff << " +- "
00944                             << std::right << std::setw(10) << std::setprecision(2)
00945                             << 100.*statErrEff << "\n";
00946 
00947 
00948                     }
00949 
00950                 }
00951 
00952                 myCout
00953                     << "\n\n"
00954                     << std::right << std::setw(35) << "Technical Trigger Key" << " "
00955                     << std::right << std::setw(10) << "Passed" << " "
00956                     << std::right << std::setw(10) << "Rejected" << " "
00957                     << std::right << std::setw(10) << "Error" << "\n";
00958 
00959                 for (CItEntry itEntry = m_entryListTechTrig.begin(); itEntry
00960                         != m_entryListTechTrig.end(); itEntry++) {
00961 
00962                     if ( ( ( *itEntry )->gtDaqPartition() == m_physicsDaqPartition )
00963                             && ( ( *itEntry )->gtTriggerMenuName() == *itMenu )) {
00964 
00965                         int nrEventsAccept = ( *itEntry )->gtNrEventsAccept();
00966                         int nrEventsReject = ( *itEntry )->gtNrEventsReject();
00967                         int nrEventsError = ( *itEntry )->gtNrEventsError();
00968 
00969                         myCout
00970                             << std::right << std::setw(35) << (( *itEntry )->gtAlgoName()) << " "
00971                             << std::right << std::setw(10) << nrEventsAccept << " "
00972                             << std::right << std::setw(10) << nrEventsReject << " "
00973                             << std::right << std::setw(10) << nrEventsError << "\n";
00974 
00975                     }
00976                 }
00977 
00978                 // efficiency and its statistical error
00979 
00980                 myCout << "\n\n"
00981                     << std::right << std::setw(35) << "Technical Trigger Key" << "    "
00982                     << std::right << std::setw(10) << "Efficiency " << " "
00983                     << std::right << std::setw(10) << "Stat error (%)" << "\n";
00984 
00985                 for (CItEntry itEntry = m_entryListTechTrig.begin(); itEntry
00986                         != m_entryListTechTrig.end(); itEntry++) {
00987 
00988                     if ( ( ( *itEntry )->gtDaqPartition() == 0 ) && ( ( *itEntry )->gtTriggerMenuName()
00989                             == *itMenu )) {
00990 
00991                         int nrEventsAccept = ( *itEntry )->gtNrEventsAccept();
00992                         int nrEventsReject = ( *itEntry )->gtNrEventsReject();
00993                         int nrEventsError = ( *itEntry )->gtNrEventsError();
00994 
00995                         int totalEvents = nrEventsAccept + nrEventsReject + nrEventsError;
00996 
00997                         // efficiency and their statistical error
00998                         float eff = 0.;
00999                         float statErrEff = 0.;
01000 
01001                         if (totalEvents != 0) {
01002                             eff = static_cast<float> (nrEventsAccept)
01003                                     / static_cast<float> (totalEvents);
01004                             statErrEff = sqrt(eff * ( 1.0 - eff )
01005                                     / static_cast<float> (totalEvents));
01006 
01007                         }
01008 
01009                         myCout
01010                             << std::right << std::setw(35) << (( *itEntry )->gtAlgoName()) << " "
01011                             << std::right << std::setw(10) << std::fixed << std::setprecision(2)
01012                             << 100.*eff << " +- "
01013                             << std::right << std::setw(10) << std::setprecision(2)
01014                             << 100.*statErrEff << "\n";
01015 
01016 
01017                     }
01018 
01019                 }
01020 
01021             }
01022 
01023 
01024         }
01025             break;
01026 
01027         case 10: {
01028 
01029             myCout << "\nL1T-Report " << "---------- L1 Trigger Global Summary - DAQ Partition "
01030                     << m_physicsDaqPartition << "----------\n\n";
01031 
01032             for (CItL1Menu itMenu = menuList.begin(); itMenu != menuList.end(); itMenu++) {
01033 
01034                 myCout << "Report for L1 menu:  " << ( *itMenu ) << "\n"
01035                         << std::endl;
01036                 myCout
01037                     << std::right << std::setw(45) << "Algorithm Key" << " "
01038                     << std::right << std::setw(10) << "Prescale" << "  "
01039                     << std::right << std::setw(5)  << "Mask" << " "
01040                     << std::right << std::setw(25) << "Before Mask" << " "
01041                     << std::right << std::setw(30) << "After Mask" << " "
01042                     << std::right << std::setw(22) << "Error"
01043                     << "\n"
01044                     << std::right << std::setw(64) << " " << std::setw(15) << "Passed"
01045                     << std::right << std::setw(1)  << " " << std::setw(15) << "Rejected"
01046                     << std::right << std::setw(1)  << " " << std::setw(15) << "Passed"
01047                     << std::right << std::setw(1)  << " " << std::setw(15) << "Rejected"
01048                     << "\n";
01049 
01050                 for (CItEntry itEntry = m_entryList.begin(); itEntry != m_entryList.end(); itEntry++) {
01051 
01052                     if ( ( *itEntry )->gtDaqPartition() == m_physicsDaqPartition) {
01053                         myCout
01054                             << std::right << std::setw(45) << ( *itEntry )->gtAlgoName() << " "
01055                             << std::right << std::setw(10) << ( *itEntry )->gtPrescaleFactor() << " "
01056                             << std::right << std::setw(5)  << " " << std::hex << ( *itEntry )->gtTriggerMask() << std::dec << " "
01057                             << std::right << std::setw(15) << ( *itEntry )->gtNrEventsAcceptBeforeMask() << " "
01058                             << std::right << std::setw(15) << ( *itEntry )->gtNrEventsRejectBeforeMask() << " "
01059                             << std::right << std::setw(15) << ( *itEntry )->gtNrEventsAccept() << " "
01060                             << std::right << std::setw(15) << ( *itEntry )->gtNrEventsReject() << " "
01061                             << std::right << std::setw(15) << ( *itEntry )->gtNrEventsError()
01062                             << "\n";
01063                     }
01064                 }
01065 
01066                 myCout
01067                     << "\n\n"
01068                     << std::right << std::setw(45) << "Technical Trigger Key" << " "
01069                     << std::right << std::setw(10) << "Prescale" << "  "
01070                     << std::right << std::setw(5)  << "Mask" << " "
01071                     << std::right << std::setw(25) << "Before Mask" << " "
01072                     << std::right << std::setw(30) << "After Mask" << " "
01073                     << std::right << std::setw(22) << "Error"
01074                     << "\n"
01075                     << std::right << std::setw(64) << " " << std::setw(15) << "Passed"
01076                     << std::right << std::setw(1)  << " " << std::setw(15) << "Rejected"
01077                     << std::right << std::setw(1)  << " " << std::setw(15) << "Passed"
01078                     << std::right << std::setw(1)  << " " << std::setw(15) << "Rejected"
01079                     << "\n";
01080 
01081                 for (CItEntry itEntry = m_entryListTechTrig.begin(); itEntry != m_entryListTechTrig.end(); itEntry++) {
01082 
01083                     if ( ( *itEntry )->gtDaqPartition() == m_physicsDaqPartition) {
01084                         myCout
01085                             << std::right << std::setw(45) << ( *itEntry )->gtAlgoName() << " "
01086                             << std::right << std::setw(10) << ( *itEntry )->gtPrescaleFactor() << " "
01087                             << std::right << std::setw(5)  << " " << std::hex << ( *itEntry )->gtTriggerMask() << std::dec << " "
01088                             << std::right << std::setw(15) << ( *itEntry )->gtNrEventsAcceptBeforeMask() << " "
01089                             << std::right << std::setw(15) << ( *itEntry )->gtNrEventsRejectBeforeMask() << " "
01090                             << std::right << std::setw(15) << ( *itEntry )->gtNrEventsAccept() << " "
01091                             << std::right << std::setw(15) << ( *itEntry )->gtNrEventsReject() << " "
01092                             << std::right << std::setw(15) << ( *itEntry )->gtNrEventsError()
01093                             << "\n";
01094                     }
01095                 }
01096             }
01097         }
01098 
01099             break;
01100         case 100: {
01101 
01102             for (unsigned int iDaqPartition = 0; iDaqPartition < m_numberDaqPartitions; ++iDaqPartition) {
01103 
01104                 myCout << "\nL1T-Report "
01105                         << "---------- L1 Trigger Global Summary - DAQ Partition " << iDaqPartition
01106                         << " " << "----------\n\n";
01107 
01108                 myCout
01109                     << std::right << std::setw(35) << "Algorithm Key" << " "
01110                     << std::right << std::setw(10) << "Passed" << " "
01111                     << std::right << std::setw(10) << "Rejected" << " "
01112                     << std::right << std::setw(10) << "Error" << std::setw(2) << " "
01113                     << "\n";
01114 
01115                 for (CItEntry itEntry = m_entryList.begin(); itEntry != m_entryList.end(); itEntry++) {
01116 
01117                     if ( ( *itEntry )->gtDaqPartition() == 0) {
01118 
01119                         myCout
01120                             << std::right << std::setw(35) << ( *itEntry )->gtAlgoName() << " "
01121                             << std::right << std::setw(10) << ( *itEntry )->gtNrEventsAccept() << " "
01122                             << std::right << std::setw(10) << ( *itEntry )->gtNrEventsReject() << " "
01123                             << std::right << std::setw(10) << ( *itEntry )->gtNrEventsError() << std::setw(2) << " "
01124                             << "\n";
01125                     }
01126 
01127                 }
01128 
01129                 myCout
01130                     << "\n\n"
01131                     << std::right << std::setw(35) << "Technical Trigger Key" << " "
01132                     << std::right << std::setw(10) << "Passed" << " "
01133                     << std::right << std::setw(10) << "Rejected" << " "
01134                     << std::right << std::setw(10) << "Error" << std::setw(2) << " "
01135                     << "\n";
01136 
01137                 for (CItEntry itEntry = m_entryListTechTrig.begin(); itEntry != m_entryListTechTrig.end(); itEntry++) {
01138 
01139                     if ( ( *itEntry )->gtDaqPartition() == 0) {
01140 
01141                         myCout
01142                             << std::right << std::setw(35) << ( *itEntry )->gtAlgoName() << " "
01143                             << std::right << std::setw(10) << ( *itEntry )->gtNrEventsAccept() << " "
01144                             << std::right << std::setw(10) << ( *itEntry )->gtNrEventsReject() << " "
01145                             << std::right << std::setw(10) << ( *itEntry )->gtNrEventsError() << std::setw(2) << " "
01146                             << std::right << std::setw(20) << ( *itEntry )->gtTriggerMenuName()
01147                             << "\n";
01148                     }
01149 
01150                 }
01151 
01152             }
01153         }
01154 
01155             break;
01156         case 101: {
01157 
01158             for (unsigned int iDaqPartition = 0; iDaqPartition < m_numberDaqPartitions; ++iDaqPartition) {
01159 
01160                 myCout << "\nL1T-Report "
01161                         << "---------- L1 Trigger Global Summary - DAQ Partition " << iDaqPartition
01162                         << " " << "----------\n\n";
01163 
01164                 myCout
01165                     << std::right << std::setw(35) << "Algorithm Key" << " "
01166                     << std::right << std::setw(10) << "Prescale" << " "
01167                     << std::right << std::setw(5)  << "Mask" << " "
01168                     << std::right << std::setw(10) << "Passed" << " "
01169                     << std::right << std::setw(10) << "Rejected" << " "
01170                     << std::right << std::setw(10) << "Error" << std::setw(2) << " "
01171                     << "\n";
01172 
01173                 for (CItEntry itEntry = m_entryList.begin(); itEntry != m_entryList.end(); itEntry++) {
01174 
01175                     if ( ( *itEntry )->gtDaqPartition() == 0) {
01176                         myCout
01177                             << std::right << std::setw(35) << ( *itEntry )->gtAlgoName() << " "
01178                             << std::right << std::setw(10) << ( *itEntry )->gtPrescaleFactor() << "   "
01179                             << std::right << std::setw(2) //<< std::setfill('0')
01180                             << std::hex << ( *itEntry )->gtTriggerMask() //<< std::setfill(' ')
01181                             << std::dec << " "
01182                             << std::right << std::setw(10) << ( *itEntry )->gtNrEventsAccept() << " "
01183                             << std::right << std::setw(10) << ( *itEntry )->gtNrEventsReject() << " "
01184                             << std::right << std::setw(10) << ( *itEntry )->gtNrEventsError() << std::setw(2) << " "
01185                             << "\n";
01186                     }
01187                 }
01188 
01189                 myCout
01190                     << "\n\n"
01191                     << std::right << std::setw(35) << "Technical Trigger Key" << " "
01192                     << std::right << std::setw(10) << "Prescale" << " "
01193                     << std::right << std::setw(5)  << "Mask" << " "
01194                     << std::right << std::setw(10) << "Passed" << " "
01195                     << std::right << std::setw(10) << "Rejected" << " "
01196                     << std::right << std::setw(10) << "Error" << std::setw(2) << " "
01197                     << "\n";
01198 
01199                 for (CItEntry itEntry = m_entryListTechTrig.begin(); itEntry != m_entryListTechTrig.end(); itEntry++) {
01200 
01201                     if ( ( *itEntry )->gtDaqPartition() == 0) {
01202                         myCout
01203                             << std::right << std::setw(35) << ( *itEntry )->gtAlgoName() << " "
01204                             << std::right << std::setw(10) << ( *itEntry )->gtPrescaleFactor() << "   "
01205                             << std::right << std::setw(2) //<< std::setfill('0')
01206                             << std::hex << ( *itEntry )->gtTriggerMask() //<< std::setfill(' ')
01207                             << std::dec << " "
01208                             << std::right << std::setw(10) << ( *itEntry )->gtNrEventsAccept() << " "
01209                             << std::right << std::setw(10) << ( *itEntry )->gtNrEventsReject() << " "
01210                             << std::right << std::setw(10) << ( *itEntry )->gtNrEventsError() << std::setw(2) << " "
01211                             << "\n";
01212                     }
01213                 }
01214 
01215             }
01216         }
01217 
01218             break;
01219         default: {
01220             myCout
01221                 << "\n\nL1GtTrigReport: Error - no print verbosity level = " << m_printVerbosity
01222                 << " defined! \nCheck available values in the cfi file." << "\n";
01223         }
01224 
01225             break;
01226     }
01227 
01228     // TODO for other verbosity levels
01229     // print the trigger menu, the prescale factors and the trigger mask, etc
01230 
01231 
01232     myCout << std::endl;
01233     myCout << "L1T-Report end!" << std::endl;
01234     myCout << std::endl;
01235 
01236     switch (m_printOutput) {
01237         case 0: {
01238 
01239             std::cout << myCout.str() << std::endl;
01240 
01241         }
01242 
01243             break;
01244         case 1: {
01245 
01246             LogTrace("L1GtTrigReport") << myCout.str() << std::endl;
01247 
01248         }
01249             break;
01250 
01251         case 2: {
01252 
01253             edm::LogVerbatim("L1GtTrigReport") << myCout.str() << std::endl;
01254 
01255         }
01256 
01257             break;
01258         case 3: {
01259 
01260             edm::LogInfo("L1GtTrigReport") << myCout.str();
01261 
01262         }
01263 
01264             break;
01265         default: {
01266             std::cout
01267                 << "\n\n  L1GtTrigReport: Error - no print output = " << m_printOutput
01268                 << " defined! \n  Check available values in the cfi file." << "\n" << std::endl;
01269 
01270         }
01271             break;
01272     }
01273 
01274 }
01275