#include <PFJetDQMAnalyzer.h>
Public Member Functions | |
PFJetDQMAnalyzer (const edm::ParameterSet ¶meterSet) | |
Private Member Functions | |
void | analyze (edm::Event const &, edm::EventSetup const &) |
void | beginJob () |
void | endJob () |
void | storeBadEvents (edm::Event const &, float &val) |
Private Attributes | |
std::string | benchmarkLabel_ |
edm::InputTag | inputLabel_ |
edm::InputTag | matchLabel_ |
int | nBadEvents_ |
PFJetMonitor | pfJetMonitor_ |
edm::ParameterSet | pSet_ |
Definition at line 13 of file PFJetDQMAnalyzer.h.
PFJetDQMAnalyzer::PFJetDQMAnalyzer | ( | const edm::ParameterSet & | parameterSet | ) |
Definition at line 22 of file PFJetDQMAnalyzer.cc.
References benchmarkLabel_, edm::ParameterSet::getParameter(), inputLabel_, matchLabel_, pfJetMonitor_, pSet_, and PFJetMonitor::setParameters().
{ pSet_ = parameterSet; inputLabel_ = pSet_.getParameter<edm::InputTag>("InputCollection"); matchLabel_ = pSet_.getParameter<edm::InputTag>("MatchCollection"); benchmarkLabel_ = pSet_.getParameter<std::string>("BenchmarkLabel"); pfJetMonitor_.setParameters(parameterSet); }
void PFJetDQMAnalyzer::analyze | ( | edm::Event const & | iEvent, |
edm::EventSetup const & | iSetup | ||
) | [private, virtual] |
Implements edm::EDAnalyzer.
Definition at line 49 of file PFJetDQMAnalyzer.cc.
References PFJetMonitor::fill(), edm::Event::getByLabel(), edm::ParameterSet::getParameter(), inputLabel_, edm::HandleBase::isValid(), matchLabel_, nBadEvents_, pfJetMonitor_, pSet_, and storeBadEvents().
{ edm::Handle< edm::View<reco::Jet> > jetCollection; iEvent.getByLabel(inputLabel_, jetCollection); edm::Handle< edm::View<reco::Jet> > matchedJetCollection; iEvent.getByLabel( matchLabel_, matchedJetCollection); float maxRes = 0.0; float minRes = 99.99; if (jetCollection.isValid() && matchedJetCollection.isValid()) { pfJetMonitor_.fill( *jetCollection, *matchedJetCollection, minRes, maxRes); edm::ParameterSet skimPS = pSet_.getParameter<edm::ParameterSet>("SkimParameter"); if ( (skimPS.getParameter<bool>("switchOn")) && (nBadEvents_ <= skimPS.getParameter<int32_t>("maximumNumberToBeStored")) ){ if ( minRes < skimPS.getParameter<double>("lowerCutOffOnResolution")) { storeBadEvents(iEvent,minRes); nBadEvents_++; } else if (maxRes > skimPS.getParameter<double>("upperCutOffOnResolution")) { storeBadEvents(iEvent,maxRes); nBadEvents_++; } } } }
void PFJetDQMAnalyzer::beginJob | ( | void | ) | [private, virtual] |
Reimplemented from edm::EDAnalyzer.
Definition at line 36 of file PFJetDQMAnalyzer.cc.
References benchmarkLabel_, Benchmark::DQM_, nBadEvents_, cppFunctionSkipper::operator, scaleCards::path, pfJetMonitor_, pSet_, DQMStore::setCurrentFolder(), and PFJetMonitor::setup().
{ Benchmark::DQM_ = edm::Service<DQMStore>().operator->(); // part of the following could be put in the base class std::string path = "ParticleFlow/" + benchmarkLabel_; Benchmark::DQM_->setCurrentFolder(path.c_str()); edm::LogInfo("PFJetDQMAnalyzer") << " PFJetDQMAnalyzer::beginJob " << "Histogram Folder path set to "<< path; pfJetMonitor_.setup(pSet_); nBadEvents_ = 0; }
void PFJetDQMAnalyzer::endJob | ( | void | ) | [private, virtual] |
void PFJetDQMAnalyzer::storeBadEvents | ( | edm::Event const & | iEvent, |
float & | val | ||
) | [private] |
Definition at line 75 of file PFJetDQMAnalyzer.cc.
References benchmarkLabel_, DQMStore::bookFloat(), Benchmark::DQM_, edm::EventID::event(), MonitorElement::Fill(), DQMStore::get(), edm::EventBase::id(), edm::EventID::luminosityBlock(), scaleCards::path, MonitorElement::Reset(), edm::EventID::run(), and DQMStore::setCurrentFolder().
Referenced by analyze().
{ unsigned int runNb = iEvent.id().run(); unsigned int evtNb = iEvent.id().event(); unsigned int lumiNb = iEvent.id().luminosityBlock(); std::string path = "ParticleFlow/" + benchmarkLabel_ + "/BadEvents"; Benchmark::DQM_->setCurrentFolder(path.c_str()); std::ostringstream eventid_str; eventid_str << runNb << "_"<< evtNb << "_" << lumiNb; MonitorElement* me = Benchmark::DQM_->get(path + "/" + eventid_str.str()); if (me) me->Reset(); else me = Benchmark::DQM_->bookFloat(eventid_str.str()); me->Fill(val); }
std::string PFJetDQMAnalyzer::benchmarkLabel_ [private] |
Definition at line 27 of file PFJetDQMAnalyzer.h.
Referenced by beginJob(), PFJetDQMAnalyzer(), and storeBadEvents().
edm::InputTag PFJetDQMAnalyzer::inputLabel_ [private] |
Definition at line 26 of file PFJetDQMAnalyzer.h.
Referenced by analyze(), and PFJetDQMAnalyzer().
edm::InputTag PFJetDQMAnalyzer::matchLabel_ [private] |
Definition at line 25 of file PFJetDQMAnalyzer.h.
Referenced by analyze(), and PFJetDQMAnalyzer().
int PFJetDQMAnalyzer::nBadEvents_ [private] |
Definition at line 32 of file PFJetDQMAnalyzer.h.
Referenced by analyze(), and beginJob().
PFJetMonitor PFJetDQMAnalyzer::pfJetMonitor_ [private] |
Definition at line 29 of file PFJetDQMAnalyzer.h.
Referenced by analyze(), beginJob(), and PFJetDQMAnalyzer().
edm::ParameterSet PFJetDQMAnalyzer::pSet_ [private] |
Definition at line 31 of file PFJetDQMAnalyzer.h.
Referenced by analyze(), beginJob(), and PFJetDQMAnalyzer().