File indexing completed on 2024-04-06 12:18:57
0001 #include "DataFormats/HLTReco/interface/TriggerObject.h"
0002 #include "FWCore/Common/interface/TriggerNames.h"
0003 #include "FWCore/Framework/interface/Frameworkfwd.h"
0004 #include "FWCore/Framework/interface/MakerMacros.h"
0005 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0006 #include "HLTriggerOffline/SUSYBSM/interface/SUSY_HLT_ElecFakes.h"
0007 #include <iostream>
0008
0009 SUSY_HLT_ElecFakes::SUSY_HLT_ElecFakes(const edm::ParameterSet &ps) {
0010 edm::LogInfo("SUSY_HLT_ElecFakes") << "Constructor SUSY_HLT_ElecFakes::SUSY_HLT_ElecFakes " << std::endl;
0011
0012 theTrigSummary_ = consumes<trigger::TriggerEvent>(ps.getParameter<edm::InputTag>("trigSummary"));
0013 triggerResults_ = consumes<edm::TriggerResults>(ps.getParameter<edm::InputTag>("TriggerResults"));
0014 HLTProcess_ = ps.getParameter<std::string>("HLTProcess");
0015 triggerPath_ = ps.getParameter<std::string>("TriggerPath");
0016 triggerFilter_ = ps.getParameter<edm::InputTag>("TriggerFilter");
0017 triggerJetFilter_ = ps.getParameter<edm::InputTag>("TriggerJetFilter");
0018 }
0019
0020 SUSY_HLT_ElecFakes::~SUSY_HLT_ElecFakes() {
0021 edm::LogInfo("SUSY_HLT_ElecFakes") << "Destructor SUSY_HLT_ElecFakes::~SUSY_HLT_ElecFakes " << std::endl;
0022 }
0023
0024 void SUSY_HLT_ElecFakes::dqmBeginRun(edm::Run const &run, edm::EventSetup const &e) {
0025 bool changed;
0026
0027 if (!fHltConfig.init(run, e, HLTProcess_, changed)) {
0028 edm::LogError("SUSY_HLT_ElecFakes") << "Initialization of HLTConfigProvider failed!!";
0029 return;
0030 }
0031
0032 bool pathFound = false;
0033 const std::vector<std::string> allTrigNames = fHltConfig.triggerNames();
0034 for (size_t j = 0; j < allTrigNames.size(); ++j) {
0035 if (allTrigNames[j].find(triggerPath_) != std::string::npos) {
0036 pathFound = true;
0037 }
0038 }
0039
0040 if (!pathFound) {
0041 edm::LogInfo("SUSY_HLT_ElecFakes") << "Path not found"
0042 << "\n";
0043 return;
0044 }
0045
0046 edm::LogInfo("SUSY_HLT_ElecFakes") << "SUSY_HLT_ElecFakes::beginRun" << std::endl;
0047 }
0048
0049 void SUSY_HLT_ElecFakes::bookHistograms(DQMStore::IBooker &ibooker_, edm::Run const &, edm::EventSetup const &) {
0050 edm::LogInfo("SUSY_HLT_ElecFakes") << "SUSY_HLT_ElecFakes::bookHistograms" << std::endl;
0051
0052 bookHistos(ibooker_);
0053 }
0054
0055 void SUSY_HLT_ElecFakes::analyze(edm::Event const &e, edm::EventSetup const &eSetup) {
0056 edm::LogInfo("SUSY_HLT_ElecFakes") << "SUSY_HLT_ElecFakes::analyze" << std::endl;
0057
0058
0059
0060
0061 edm::Handle<edm::TriggerResults> hltresults;
0062 e.getByToken(triggerResults_, hltresults);
0063 if (!hltresults.isValid()) {
0064 edm::LogError("SUSY_HLT_ElecFakes") << "invalid collection: TriggerResults"
0065 << "\n";
0066 return;
0067 }
0068 edm::Handle<trigger::TriggerEvent> triggerSummary;
0069 e.getByToken(theTrigSummary_, triggerSummary);
0070 if (!triggerSummary.isValid()) {
0071 edm::LogError("SUSY_HLT_ElecFakes") << "invalid collection: TriggerSummary"
0072 << "\n";
0073 return;
0074 }
0075
0076
0077 size_t filterIndex = triggerSummary->filterIndex(triggerFilter_);
0078 trigger::TriggerObjectCollection triggerObjects = triggerSummary->getObjects();
0079 if (!(filterIndex >= triggerSummary->sizeFilters())) {
0080 const trigger::Keys &keys = triggerSummary->filterKeys(filterIndex);
0081 for (size_t j = 0; j < keys.size(); ++j) {
0082 trigger::TriggerObject foundObject = triggerObjects[keys[j]];
0083
0084 h_triggerElPt->Fill(foundObject.pt());
0085 h_triggerElEta->Fill(foundObject.eta());
0086 h_triggerElPhi->Fill(foundObject.phi());
0087
0088 }
0089 }
0090
0091 filterIndex = triggerSummary->filterIndex(triggerJetFilter_);
0092
0093 if (!(filterIndex >= triggerSummary->sizeFilters())) {
0094 const trigger::Keys &keys = triggerSummary->filterKeys(filterIndex);
0095 for (size_t j = 0; j < keys.size(); ++j) {
0096 trigger::TriggerObject foundObject = triggerObjects[keys[j]];
0097 h_triggerJetPt->Fill(foundObject.pt());
0098 h_triggerJetEta->Fill(foundObject.eta());
0099 h_triggerJetPhi->Fill(foundObject.phi());
0100 }
0101 }
0102 }
0103
0104 void SUSY_HLT_ElecFakes::bookHistos(DQMStore::IBooker &ibooker_) {
0105 ibooker_.cd();
0106 ibooker_.setCurrentFolder("HLT/SUSYBSM/" + triggerPath_);
0107
0108
0109 h_triggerElPt = ibooker_.book1D("triggerElPt", "Trigger El Pt; GeV", 50, 0.0, 100.0);
0110 h_triggerElEta = ibooker_.book1D("triggerElEta", "Trigger El Eta", 20, -2.5, 2.5);
0111 h_triggerElPhi = ibooker_.book1D("triggerElPhi", "Trigger El Phi", 20, -3.5, 3.5);
0112
0113 h_triggerJetPt = ibooker_.book1D("triggerJetPt", "Trigger Jet Pt; GeV", 20, 0.0, 200.0);
0114 h_triggerJetEta = ibooker_.book1D("triggerJetEta", "Trigger Jet Eta", 20, -3.0, 3.0);
0115 h_triggerJetPhi = ibooker_.book1D("triggerJetPhi", "Trigger Jet Phi", 20, -3.5, 3.5);
0116
0117
0118
0119
0120
0121 ibooker_.cd();
0122 }
0123
0124
0125 DEFINE_FWK_MODULE(SUSY_HLT_ElecFakes);