Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0002 #include "FWCore/Framework/interface/Frameworkfwd.h"
0003 #include "FWCore/Framework/interface/MakerMacros.h"
0004 #include "FWCore/Common/interface/TriggerNames.h"
0005 #include "DataFormats/HLTReco/interface/TriggerObject.h"
0006 #include "DQMOffline/Trigger/interface/EventShapeDQM.h"
0007 
0008 EventShapeDQM::EventShapeDQM(const edm::ParameterSet& ps) {
0009   triggerResults_ = consumes<edm::TriggerResults>(ps.getParameter<edm::InputTag>("triggerResults"));
0010   theEPCollection_ = consumes<reco::EvtPlaneCollection>(ps.getParameter<edm::InputTag>("EPlabel"));
0011   triggerPath_ = ps.getParameter<std::string>("triggerPath");
0012 
0013   order_ = ps.getParameter<int>("order");
0014   EPidx_ = ps.getParameter<int>("EPidx");
0015   EPlvl_ = ps.getParameter<int>("EPlvl");
0016 }
0017 
0018 EventShapeDQM::~EventShapeDQM() = default;
0019 
0020 void EventShapeDQM::bookHistograms(DQMStore::IBooker& ibooker_, edm::Run const&, edm::EventSetup const&) {
0021   ibooker_.cd();
0022   ;
0023   ibooker_.setCurrentFolder("HLT/HI/" + triggerPath_);
0024 
0025   h_Q = ibooker_.book1D("hQn", Form("Q%i;Q%i", order_, order_), 500, 0, 0.5);
0026 
0027   ibooker_.cd();
0028 }
0029 
0030 void EventShapeDQM::analyze(edm::Event const& e, edm::EventSetup const& eSetup) {
0031   edm::Handle<edm::TriggerResults> hltresults;
0032   e.getByToken(triggerResults_, hltresults);
0033   if (!hltresults.isValid()) {
0034     return;
0035   }
0036 
0037   bool hasFired = false;
0038   const edm::TriggerNames& trigNames = e.triggerNames(*hltresults);
0039   unsigned int numTriggers = trigNames.size();
0040   for (unsigned int hltIndex = 0; hltIndex < numTriggers; ++hltIndex) {
0041     if (trigNames.triggerName(hltIndex).find(triggerPath_) != std::string::npos && hltresults->wasrun(hltIndex) &&
0042         hltresults->accept(hltIndex)) {
0043       hasFired = true;
0044     }
0045   }
0046 
0047   edm::Handle<reco::EvtPlaneCollection> ep_;
0048   e.getByToken(theEPCollection_, ep_);
0049   if (!ep_.isValid()) {
0050     return;
0051   }
0052 
0053   if (hasFired) {
0054     h_Q->Fill((*ep_)[EPidx_].vn(EPlvl_));
0055   }
0056 }
0057 
0058 DEFINE_FWK_MODULE(EventShapeDQM);