Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:09:13

0001 //
0002 // Revision 1.30  2011/03/30 21:44:03  fwyzard
0003 // make sure HLTConfigProvider is used only if succesfully initialized
0004 //
0005 // Revision 1.29  2011/03/30 21:35:40  fwyzard
0006 // make sure all members are initialized
0007 //
0008 // Revision 1.28  2011/03/29 09:46:03  rekovic
0009 // clean vector pairPDPaths in beginRun and tidy up
0010 //
0011 // Revision 1.27  2011/03/24 18:35:38  rekovic
0012 // Change name for pd histo
0013 //
0014 // Revision 1.26  2011/03/24 18:25:45  rekovic
0015 // Add single 1D plot of streamA content
0016 //
0017 // Revision 1.25  2010/07/20 02:58:27  wmtan
0018 // Add missing #include files
0019 //
0020 // Revision 1.24  2010/03/17 20:54:51  wittich
0021 // add scalers that I manually reset on beginLumi
0022 //
0023 // Revision 1.23  2010/02/25 17:34:01  wdd
0024 // Central migration of TriggerNames class interface
0025 //
0026 // Revision 1.22  2010/02/24 17:43:47  wittich
0027 // - keep trying to get path names if it doesn't work first time
0028 // - move the Bx histograms out of raw to the toplevel directory.
0029 //
0030 // Revision 1.21  2010/02/11 23:54:28  wittich
0031 // modify how the monitoring histo is filled
0032 //
0033 // Revision 1.20  2010/02/11 00:11:08  wmtan
0034 // Adapt to moved framework header
0035 //
0036 // Revision 1.19  2010/02/02 13:53:05  wittich
0037 // fix duplicate histogram name
0038 //
0039 // Revision 1.18  2010/02/02 11:42:53  wittich
0040 // new diagnostic histograms
0041 //
0042 // Revision 1.17  2009/11/20 00:39:12  lorenzo
0043 // fixes
0044 //
0045 // Revision 1.16  2008/09/03 13:59:06  wittich
0046 // make HLT DQM path configurable via python parameter,
0047 // which defaults to HLT/HLTScalers_EvF
0048 //
0049 // Revision 1.15  2008/09/03 02:13:47  wittich
0050 // - bug fix in L1Scalers
0051 // - configurable dqm directory in L1SCalers
0052 // - other minor tweaks in HLTScalers
0053 //
0054 
0055 #include <iostream>
0056 
0057 // FW
0058 #include "FWCore/Framework/interface/Event.h"
0059 #include "FWCore/Framework/interface/Run.h"
0060 
0061 #include "FWCore/ServiceRegistry/interface/Service.h"
0062 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0063 
0064 #include "FWCore/Framework/interface/LuminosityBlock.h"
0065 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0066 
0067 // HLT
0068 #include "DataFormats/Common/interface/TriggerResults.h"
0069 #include "DataFormats/Common/interface/HLTenums.h"
0070 #include "FWCore/Common/interface/TriggerNames.h"
0071 
0072 #include "DQM/TrigXMonitor/interface/HLTScalers.h"
0073 #include "DataFormats/Common/interface/Handle.h"
0074 #include "DQMServices/Core/interface/DQMStore.h"
0075 
0076 using namespace edm;
0077 
0078 HLTScalers::HLTScalers(const edm::ParameterSet& ps)
0079     : folderName_(ps.getUntrackedParameter<std::string>("dqmFolder", "HLT/HLTScalers_EvF")),
0080       processname_(ps.getParameter<std::string>("processname")),
0081       pairPDPaths_(),
0082       trigResultsSource_(consumes<TriggerResults>(ps.getParameter<edm::InputTag>("triggerResults"))),
0083       scalersN_(nullptr),
0084       scalersException_(nullptr),
0085       hltCorrelations_(nullptr),
0086       detailedScalers_(nullptr),
0087       nProc_(nullptr),
0088       nLumiBlock_(nullptr),
0089       hltBx_(nullptr),
0090       hltBxVsPath_(nullptr),
0091       hltOverallScaler_(nullptr),
0092       hltOverallScalerN_(nullptr),
0093       diagnostic_(nullptr),
0094       sentPaths_(false),
0095       monitorDaemon_(ps.getUntrackedParameter<bool>("MonitorDaemon", false)),
0096       nev_(0),
0097       nLumi_(0) {
0098   LogDebug("HLTScalers") << "HLTScalers: constructor....";
0099 }
0100 
0101 void HLTScalers::dqmBeginRun(const edm::Run& run, const edm::EventSetup& c) {
0102   LogDebug("HLTScalers") << "HLTScalers::beginRun, run " << run.id();
0103 
0104   // HLT config does not change within runs!
0105   bool changed = false;
0106 
0107   // clear vector pairPDPaths_
0108   pairPDPaths_.clear();
0109 
0110   if (not hltConfig_.init(run, c, processname_, changed)) {
0111     edm::LogError("TrigXMonitor") << "HLTConfigProvider failed to initialize.";
0112   } else {
0113     // check if trigger name in (new) config
0114     //  cout << "Available TriggerNames are: " << endl;
0115     //  hltConfig_.dump("Triggers");
0116 
0117     if (hltConfig_.streamIndex("A") < hltConfig_.streamNames().size()) {
0118       // get hold of PD names and constituent path names
0119       const std::vector<std::string>& PD = hltConfig_.streamContent("A");
0120 
0121       for (unsigned int i = 0; i < PD.size(); i++) {
0122         const std::vector<std::string>& datasetPaths = hltConfig_.datasetContent(PD[i]);
0123         pairPDPaths_.push_back(make_pair(PD[i], datasetPaths));
0124       }
0125 
0126       // push stream A and its PDs
0127       pairPDPaths_.push_back(make_pair("A", PD));
0128 
0129     } else {
0130       LogDebug("HLTScalers") << "HLTScalers::beginRun, steamm A not in the HLT menu ";
0131     }
0132   }
0133 }
0134 
0135 void HLTScalers::bookHistograms(DQMStore::IBooker& iBooker, edm::Run const&, edm::EventSetup const&) {
0136   std::string rawdir(folderName_ + "/raw");
0137   iBooker.setCurrentFolder(rawdir);
0138 
0139   nProc_ = iBooker.bookInt("nProcessed");
0140   nLumiBlock_ = iBooker.bookInt("nLumiBlock");
0141   diagnostic_ = iBooker.book1D("hltMerge", "HLT merging diagnostic", 1, 0.5, 1.5);
0142   // fill for ever accepted event
0143   hltOverallScaler_ = iBooker.book1D("hltOverallScaler", "HLT Overall Scaler", 1, 0.5, 1.5);
0144   hltOverallScalerN_ = iBooker.book1D("hltOverallScalerN", "Reset HLT Overall Scaler", 1, 0.5, 1.5);
0145 
0146   // DQM: Previously the number of trigger paths was determined on the first
0147   // event, by taking the size of the htlResults to book the histogram.
0148   // Now we use the size of the hltConfig instead.
0149   int npath = hltConfig_.size();
0150   unsigned int nPD = pairPDPaths_.size();
0151 
0152   // need to get maxModules dynamically
0153   int maxModules = 200;
0154 
0155   scalersPD_ = iBooker.book1D("pdScalers", "PD scalers (stream A)", nPD, -0.5, nPD - 0.5);
0156   detailedScalers_ =
0157       iBooker.book2D("detailedHltScalers", "HLT Scalers", npath, -0.5, npath - 0.5, maxModules, 0, maxModules - 1);
0158   scalers_ = iBooker.book1D("hltScalers", "HLT scalers", npath, -0.5, npath - 0.5);
0159   scalersN_ = iBooker.book1D("hltScalersN", "Reset HLT scalers", npath, -0.5, npath - 0.5);
0160   scalersException_ = iBooker.book1D("hltExceptions", "HLT Exception scalers", npath, -0.5, npath - 0.5);
0161   hltCorrelations_ =
0162       iBooker.book2D("hltCorrelations", "HLT Scalers", npath, -0.5, npath - 0.5, npath, -0.5, npath - 0.5);
0163 
0164   // these two belong in top-level
0165   iBooker.setCurrentFolder(folderName_);
0166   hltBxVsPath_ =
0167       iBooker.book2D("hltBxVsPath", "HLT Accept vs Bunch Number", 3600, -0.5, 3599.5, npath, -0.5, npath - 0.5);
0168   hltBx_ = iBooker.book1D("hltBx", "Bx of HLT Accepted Events ", 3600, -0.5, 3599.5);
0169 }
0170 
0171 void HLTScalers::beginLuminosityBlock(const edm::LuminosityBlock& lumiSeg, const edm::EventSetup& c) {
0172   LogDebug("HLTScalers") << "Start of luminosity block.";
0173   // reset the N guys
0174   if (scalersN_)
0175     scalersN_->Reset();
0176   if (hltOverallScalerN_)
0177     hltOverallScalerN_->Reset();
0178 }
0179 
0180 void HLTScalers::analyze(const edm::Event& e, const edm::EventSetup& c) {
0181   nProc_->Fill(++nev_);
0182   diagnostic_->setBinContent(1, 1);  // this ME is never touched -
0183   // it just tells you how the merging is doing.
0184 
0185   edm::Handle<TriggerResults> hltResults;
0186   bool b = e.getByToken(trigResultsSource_, hltResults);
0187   if (!b) {
0188     Labels l;
0189     labelsForToken(trigResultsSource_, l);
0190 
0191     edm::LogInfo("HLTScalers") << "getByLabel for TriggerResults failed"
0192                                << " with label " << l.module;
0193     return;
0194   }
0195   int npath = hltResults->size();
0196   unsigned int nPD = pairPDPaths_.size();
0197   ///////////////////////////////////////////////////////////////////////////
0198 
0199   const edm::TriggerNames& trigNames = e.triggerNames(*hltResults);
0200   // for some reason this doesn't appear to work on the first event sometimes
0201   if (!sentPaths_) {
0202     const edm::TriggerNames& names = e.triggerNames(*hltResults);
0203 
0204     // save path names in DQM-accessible format
0205     int q = 0;
0206     for (TriggerNames::Strings::const_iterator j = names.triggerNames().begin(); j != names.triggerNames().end(); ++j) {
0207       LogDebug("HLTScalers") << q << ": " << *j;
0208       ++q;
0209       scalers_->getTH1()->GetXaxis()->SetBinLabel(q, j->c_str());
0210     }
0211 
0212     for (unsigned int i = 0; i < nPD; i++) {
0213       LogDebug("HLTScalers") << i << ": " << pairPDPaths_[i].first << std::endl;
0214       scalersPD_->getTH1()->GetXaxis()->SetBinLabel(i + 1, pairPDPaths_[i].first.c_str());
0215     }
0216 
0217     sentPaths_ = true;
0218   }
0219 
0220   bool accept = false;
0221   int bx = e.bunchCrossing();
0222   for (int i = 0; i < npath; ++i) {
0223     // state returns 0 on ready, 1 on accept, 2 on fail, 3 on exception.
0224     // these are defined in HLTEnums.h
0225     for (unsigned int j = 0; j < hltResults->index(i); ++j) {
0226       detailedScalers_->Fill(i, j);
0227     }
0228     if (hltResults->state(i) == hlt::Pass) {
0229       scalers_->Fill(i);
0230       scalersN_->Fill(i);
0231       hltBxVsPath_->Fill(bx, i);
0232       accept = true;
0233       for (int j = i + 1; j < npath; ++j) {
0234         if (hltResults->state(j) == hlt::Pass) {
0235           hltCorrelations_->Fill(i, j);  // fill
0236           hltCorrelations_->Fill(j, i);
0237         }
0238       }
0239     } else if (hltResults->state(i) == hlt::Exception) {
0240       scalersException_->Fill(i);
0241     }
0242   }
0243   if (accept) {
0244     hltOverallScaler_->Fill(1.0);
0245     hltOverallScalerN_->Fill(1.0);
0246     hltBx_->Fill(int(bx));
0247   }
0248 
0249   bool anyGroupPassed = false;
0250   for (unsigned int mi = 0; mi < pairPDPaths_.size(); mi++) {
0251     bool groupPassed = false;
0252 
0253     for (unsigned int i = 0; i < pairPDPaths_[mi].second.size(); i++) {
0254       // string hltPathName =  hist_2d->GetXaxis()->GetBinLabel(i);
0255       std::string hltPathName = pairPDPaths_[mi].second[i];
0256 
0257       // check if this is hlt path name
0258       // unsigned int pathByIndex = triggerNames.triggerIndex(hltPathName);
0259       unsigned int pathByIndex = trigNames.triggerIndex(pairPDPaths_[mi].second[i]);
0260       if (pathByIndex >= hltResults->size())
0261         continue;
0262 
0263       // check if its L1 passed
0264       // comment out below but set groupL1Passed to true always
0265       // if(hasL1Passed(hltPathName,triggerNames)) groupL1Passed = true;
0266       // groupL1Passed = true;
0267 
0268       // Fill HLTPassed Matrix and HLTPassFail Matrix
0269       // --------------------------------------------------------
0270 
0271       if (hltResults->accept(pathByIndex)) {
0272         groupPassed = true;
0273         break;
0274       }
0275     }
0276 
0277     if (groupPassed) {
0278       scalersPD_->Fill(mi);
0279       anyGroupPassed = true;
0280     }
0281   }
0282 
0283   if (anyGroupPassed)
0284     scalersPD_->Fill(pairPDPaths_.size() - 1);
0285 }
0286 
0287 void HLTScalers::endLuminosityBlock(const edm::LuminosityBlock& lumiSeg, const edm::EventSetup& c) {
0288   // put this in as a first-pass for figuring out the rate
0289   // each lumi block is 23 seconds in length
0290   nLumiBlock_->Fill(lumiSeg.id().luminosityBlock());
0291 
0292   LogDebug("HLTScalers") << "End of luminosity block.";
0293 }
0294 
0295 void HLTScalers::dqmEndRun(const edm::Run& run, const edm::EventSetup& c) {
0296   LogDebug("HLTScalers") << "HLTScalers::endRun , run " << run.id();
0297 }