File indexing completed on 2024-04-06 12:07:36
0001 #include "FWCore/Framework/interface/Run.h"
0002 #include "FWCore/Framework/interface/MakerMacros.h"
0003 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0004 #include "DataFormats/Common/interface/Handle.h"
0005 #include "DataFormats/HLTReco/interface/TriggerObject.h"
0006 #include "DataFormats/HLTReco/interface/TriggerEvent.h"
0007 #include "DQMServices/Core/interface/DQMStore.h"
0008
0009
0010 #include "FourVectorHLT.h"
0011
0012 using namespace edm;
0013
0014 FourVectorHLT::FourVectorHLT(const edm::ParameterSet& iConfig) {
0015 LogDebug("FourVectorHLT") << "constructor....";
0016
0017 usesResource("DQMStore");
0018 dbe_ = Service<DQMStore>().operator->();
0019 if (!dbe_) {
0020 LogWarning("Status") << "unable to get DQMStore service?";
0021 }
0022
0023 dirname_ = "HLT/FourVectorHLT";
0024
0025 if (dbe_ != nullptr) {
0026 LogDebug("Status") << "Setting current directory to " << dirname_;
0027 dbe_->setCurrentFolder(dirname_);
0028 }
0029
0030
0031 ptMin_ = iConfig.getUntrackedParameter<double>("ptMin", 0.);
0032 ptMax_ = iConfig.getUntrackedParameter<double>("ptMax", 200.);
0033 nBins_ = iConfig.getUntrackedParameter<unsigned int>("Nbins", 50);
0034
0035 plotAll_ = iConfig.getUntrackedParameter<bool>("plotAll", false);
0036
0037
0038 std::vector<edm::ParameterSet> filters = iConfig.getParameter<std::vector<edm::ParameterSet> >("filters");
0039 for (std::vector<edm::ParameterSet>::iterator filterconf = filters.begin(); filterconf != filters.end();
0040 filterconf++) {
0041 std::string me = filterconf->getParameter<std::string>("name");
0042 int objectType = filterconf->getParameter<unsigned int>("type");
0043 float ptMin = filterconf->getUntrackedParameter<double>("ptMin");
0044 float ptMax = filterconf->getUntrackedParameter<double>("ptMax");
0045 hltPaths_.push_back(PathInfo(me, objectType, ptMin, ptMax));
0046 }
0047 if (!hltPaths_.empty() && plotAll_) {
0048
0049 LogWarning("Configuration") << "Using both plotAll and a list. "
0050 "list will be ignored.";
0051 hltPaths_.clear();
0052 }
0053 triggerSummaryLabel_ = iConfig.getParameter<edm::InputTag>("triggerSummaryLabel");
0054
0055
0056 triggerSummaryToken_ = consumes<trigger::TriggerEvent>(iConfig.getParameter<edm::InputTag>("triggerSummaryLabel"));
0057 }
0058
0059 FourVectorHLT::~FourVectorHLT() {
0060
0061
0062 }
0063
0064
0065
0066
0067
0068
0069 void FourVectorHLT::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0070 using namespace edm;
0071 using namespace trigger;
0072 ++nev_;
0073 LogDebug("Status") << "analyze";
0074
0075 edm::Handle<TriggerEvent> triggerObj;
0076 iEvent.getByToken(triggerSummaryToken_, triggerObj);
0077 if (!triggerObj.isValid()) {
0078 edm::LogInfo("Status") << "Summary HLT object (TriggerEvent) not found, "
0079 "skipping event";
0080 return;
0081 }
0082
0083 const trigger::TriggerObjectCollection& toc(triggerObj->getObjects());
0084
0085 if (plotAll_) {
0086 for (size_t ia = 0; ia < triggerObj->sizeFilters(); ++ia) {
0087 std::string fullname = triggerObj->filterTag(ia).encode();
0088
0089
0090 std::string name;
0091 size_t p = fullname.find_first_of(':');
0092 if (p != std::string::npos) {
0093 name = fullname.substr(0, p);
0094 } else {
0095 name = fullname;
0096 }
0097
0098 LogDebug("Parameter") << "filter " << ia << ", full name = " << fullname << ", p = " << p
0099 << ", abbreviated = " << name;
0100
0101 PathInfoCollection::iterator pic = hltPaths_.find(name);
0102 if (pic == hltPaths_.end()) {
0103
0104 MonitorElement *et(nullptr), *eta(nullptr), *phi(nullptr), *etavsphi(nullptr);
0105
0106 std::string histoname(name + "_et");
0107 LogDebug("Status") << "new histo with name " << histoname;
0108 dbe_->setCurrentFolder(dirname_);
0109 std::string title(name + " E_{T}");
0110 et = dbe_->book1D(histoname.c_str(), title.c_str(), nBins_, 0, 100);
0111
0112 histoname = name + "_eta";
0113 title = name + " #eta";
0114 eta = dbe_->book1D(histoname.c_str(), title.c_str(), nBins_, -2.7, 2.7);
0115
0116 histoname = name + "_phi";
0117 title = name + " #phi";
0118 phi = dbe_->book1D(histoname.c_str(), title.c_str(), nBins_, -3.14, 3.14);
0119
0120 histoname = name + "_etaphi";
0121 title = name + " #eta vs #phi";
0122 etavsphi = dbe_->book2D(histoname.c_str(), title.c_str(), nBins_, -2.7, 2.7, nBins_, -3.14, 3.14);
0123
0124
0125 PathInfo e(name, 0, et, eta, phi, etavsphi, ptMin_, ptMax_);
0126 hltPaths_.push_back(e);
0127 pic = hltPaths_.begin() + hltPaths_.size() - 1;
0128 }
0129 const trigger::Keys& k = triggerObj->filterKeys(ia);
0130 for (trigger::Keys::const_iterator ki = k.begin(); ki != k.end(); ++ki) {
0131 LogDebug("Parameters") << "pt, eta, phi = " << toc[*ki].pt() << ", " << toc[*ki].eta() << ", "
0132 << toc[*ki].phi();
0133 pic->getEtHisto()->Fill(toc[*ki].pt());
0134 pic->getEtaHisto()->Fill(toc[*ki].eta());
0135 pic->getPhiHisto()->Fill(toc[*ki].phi());
0136 pic->getEtaVsPhiHisto()->Fill(toc[*ki].eta(), toc[*ki].phi());
0137 }
0138 }
0139
0140 } else {
0141 for (PathInfoCollection::iterator v = hltPaths_.begin(); v != hltPaths_.end(); ++v) {
0142 const int index = triggerObj->filterIndex(v->getName());
0143 if (index >= triggerObj->sizeFilters()) {
0144 continue;
0145 }
0146 LogDebug("Status") << "filling ... ";
0147 const trigger::Keys& k = triggerObj->filterKeys(index);
0148 for (trigger::Keys::const_iterator ki = k.begin(); ki != k.end(); ++ki) {
0149 v->getEtHisto()->Fill(toc[*ki].pt());
0150 v->getEtaHisto()->Fill(toc[*ki].eta());
0151 v->getPhiHisto()->Fill(toc[*ki].phi());
0152 v->getEtaVsPhiHisto()->Fill(toc[*ki].eta(), toc[*ki].phi());
0153 }
0154 }
0155 }
0156 }
0157
0158
0159 void FourVectorHLT::beginJob() {
0160 nev_ = 0;
0161 DQMStore* dbe = nullptr;
0162 dbe = Service<DQMStore>().operator->();
0163
0164 if (dbe) {
0165 dbe->setCurrentFolder(dirname_);
0166
0167 if (!plotAll_) {
0168 for (PathInfoCollection::iterator v = hltPaths_.begin(); v != hltPaths_.end(); ++v) {
0169 MonitorElement *et, *eta, *phi, *etavsphi = nullptr;
0170 std::string histoname(v->getName() + "_et");
0171 std::string title(v->getName() + " E_t");
0172 et = dbe->book1D(histoname.c_str(), title.c_str(), nBins_, v->getPtMin(), v->getPtMax());
0173
0174 histoname = v->getName() + "_eta";
0175 title = v->getName() + " #eta";
0176 eta = dbe->book1D(histoname.c_str(), title.c_str(), nBins_, -2.7, 2.7);
0177
0178 histoname = v->getName() + "_phi";
0179 title = v->getName() + " #phi";
0180 phi = dbe->book1D(histoname.c_str(), histoname.c_str(), nBins_, -3.14, 3.14);
0181
0182 histoname = v->getName() + "_etaphi";
0183 title = v->getName() + " #eta vs #phi";
0184 etavsphi = dbe->book2D(histoname.c_str(), title.c_str(), nBins_, -2.7, 2.7, nBins_, -3.14, 3.14);
0185
0186 v->setHistos(et, eta, phi, etavsphi);
0187 }
0188 }
0189 }
0190 }
0191
0192
0193 void FourVectorHLT::endJob() {
0194 LogInfo("Status") << "endJob: analyzed " << nev_ << " events";
0195 return;
0196 }
0197
0198
0199 void FourVectorHLT::beginRun(const edm::Run& run, const edm::EventSetup& c) {
0200 LogDebug("Status") << "beginRun, run " << run.id();
0201 }
0202
0203
0204 void FourVectorHLT::endRun(const edm::Run& run, const edm::EventSetup& c) {
0205 LogDebug("Status") << "endRun, run " << run.id();
0206 }