CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Member Functions | Private Attributes
PFMETDQMAnalyzer Class Reference

#include <PFMETDQMAnalyzer.h>

Inheritance diagram for PFMETDQMAnalyzer:
edm::EDAnalyzer

Public Member Functions

 PFMETDQMAnalyzer (const edm::ParameterSet &parameterSet)
 
- Public Member Functions inherited from edm::EDAnalyzer
 EDAnalyzer ()
 
std::string workerType () const
 
virtual ~EDAnalyzer ()
 

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_
 
PFMETMonitor pfMETMonitor_
 
edm::ParameterSet pSet_
 

Additional Inherited Members

- Public Types inherited from edm::EDAnalyzer
typedef EDAnalyzer ModuleType
 
typedef WorkerT< EDAnalyzerWorkerType
 
- Static Public Member Functions inherited from edm::EDAnalyzer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &)
 
- Protected Member Functions inherited from edm::EDAnalyzer
CurrentProcessingContext const * currentContext () const
 

Detailed Description

Definition at line 13 of file PFMETDQMAnalyzer.h.

Constructor & Destructor Documentation

PFMETDQMAnalyzer::PFMETDQMAnalyzer ( const edm::ParameterSet parameterSet)

Definition at line 20 of file PFMETDQMAnalyzer.cc.

References benchmarkLabel_, edm::ParameterSet::getParameter(), inputLabel_, matchLabel_, pfMETMonitor_, pSet_, and PFMETMonitor::setParameters().

22 {
23  pSet_ = parameterSet;
24  inputLabel_ = pSet_.getParameter<edm::InputTag>("InputCollection");
25  matchLabel_ = pSet_.getParameter<edm::InputTag>("MatchCollection");
26  benchmarkLabel_ = pSet_.getParameter<std::string>("BenchmarkLabel");
27 
28  pfMETMonitor_.setParameters(parameterSet);
29 
30 }
T getParameter(std::string const &) const
void setParameters(Benchmark::Mode mode, float ptmin, float ptmax, float etamin, float etamax, float phimin, float phimax, bool metSpHistos)
set the parameters locally
Definition: PFMETMonitor.cc:58
edm::InputTag inputLabel_
edm::InputTag matchLabel_
PFMETMonitor pfMETMonitor_
std::string benchmarkLabel_
edm::ParameterSet pSet_

Member Function Documentation

void PFMETDQMAnalyzer::analyze ( edm::Event const &  iEvent,
edm::EventSetup const &  iSetup 
)
privatevirtual

Implements edm::EDAnalyzer.

Definition at line 47 of file PFMETDQMAnalyzer.cc.

References PFMETMonitor::fillOne(), edm::Event::getByLabel(), edm::ParameterSet::getParameter(), inputLabel_, edm::HandleBase::isValid(), matchLabel_, nBadEvents_, pfMETMonitor_, pSet_, and storeBadEvents().

48  {
49  edm::Handle< edm::View<reco::MET> > metCollection;
50  iEvent.getByLabel(inputLabel_, metCollection);
51 
52  edm::Handle< edm::View<reco::MET> > matchedMetCollection;
53  iEvent.getByLabel( matchLabel_, matchedMetCollection);
54 
55  if (metCollection.isValid() && matchedMetCollection.isValid()) {
56  float maxRes = 0.0;
57  float minRes = 99.99;
58  pfMETMonitor_.fillOne( (*metCollection)[0], (*matchedMetCollection)[0], minRes, maxRes);
59  edm::ParameterSet skimPS = pSet_.getParameter<edm::ParameterSet>("SkimParameter");
60  if ( (skimPS.getParameter<bool>("switchOn")) &&
61  (nBadEvents_ <= skimPS.getParameter<int32_t>("maximumNumberToBeStored")) ) {
62  if ( minRes < skimPS.getParameter<double>("lowerCutOffOnResolution")) {
63  storeBadEvents(iEvent,minRes);
64  nBadEvents_++;
65  } else if (maxRes > skimPS.getParameter<double>("upperCutOffOnResolution")) {
66  nBadEvents_++;
67  storeBadEvents(iEvent,maxRes);
68  }
69  }
70  }
71 }
T getParameter(std::string const &) const
edm::InputTag inputLabel_
edm::InputTag matchLabel_
PFMETMonitor pfMETMonitor_
int iEvent
Definition: GenABIO.cc:243
bool isValid() const
Definition: HandleBase.h:76
void fillOne(const reco::MET &met, const reco::MET &matchedMet, float &minVal, float &maxVal)
void storeBadEvents(edm::Event const &, float &val)
edm::ParameterSet pSet_
void PFMETDQMAnalyzer::beginJob ( void  )
privatevirtual

Reimplemented from edm::EDAnalyzer.

Definition at line 34 of file PFMETDQMAnalyzer.cc.

References benchmarkLabel_, Benchmark::DQM_, nBadEvents_, cppFunctionSkipper::operator, scaleCards::path, pfMETMonitor_, pSet_, DQMStore::setCurrentFolder(), and PFMETMonitor::setup().

34  {
35 
37  // part of the following could be put in the base class
38  std::string path = "ParticleFlow/" + benchmarkLabel_;
40  edm::LogInfo("PFJMETDQMAnalyzer") << " PFMETDQMAnalyzer::beginJob " <<"Histogram Folder path set to "<< path;
42  nBadEvents_ = 0;
43 }
PFMETMonitor pfMETMonitor_
static DQMStore * DQM_
Definition: Benchmark.h:39
list path
Definition: scaleCards.py:51
std::string benchmarkLabel_
edm::ParameterSet pSet_
void setup()
book histograms
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:429
void PFMETDQMAnalyzer::endJob ( void  )
privatevirtual

Reimplemented from edm::EDAnalyzer.

Definition at line 90 of file PFMETDQMAnalyzer.cc.

90  {
91 }
void PFMETDQMAnalyzer::storeBadEvents ( edm::Event const &  iEvent,
float &  val 
)
private

Definition at line 72 of file PFMETDQMAnalyzer.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().

72  {
73  unsigned int runNb = iEvent.id().run();
74  unsigned int evtNb = iEvent.id().event();
75  unsigned int lumiNb = iEvent.id().luminosityBlock();
76 
77  std::string path = "ParticleFlow/" + benchmarkLabel_ + "/BadEvents";
78  Benchmark::DQM_->setCurrentFolder(path.c_str());
79  std::ostringstream eventid_str;
80  eventid_str << runNb << "_"<< evtNb << "_" << lumiNb;
81  MonitorElement* me = Benchmark::DQM_->get(path + "/" + eventid_str.str());
82  if (me) me->Reset();
83  else me = Benchmark::DQM_->bookFloat(eventid_str.str());
84  me->Fill(val);
85 }
static DQMStore * DQM_
Definition: Benchmark.h:39
MonitorElement * bookFloat(const char *name)
Book float.
Definition: DQMStore.cc:654
void Fill(long long x)
list path
Definition: scaleCards.py:51
int iEvent
Definition: GenABIO.cc:243
MonitorElement * get(const std::string &path) const
get ME from full pathname (e.g. &quot;my/long/dir/my_histo&quot;)
Definition: DQMStore.cc:1468
std::string benchmarkLabel_
void Reset(void)
reset ME (ie. contents, errors, etc)
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:429

Member Data Documentation

std::string PFMETDQMAnalyzer::benchmarkLabel_
private

Definition at line 27 of file PFMETDQMAnalyzer.h.

Referenced by beginJob(), PFMETDQMAnalyzer(), and storeBadEvents().

edm::InputTag PFMETDQMAnalyzer::inputLabel_
private

Definition at line 26 of file PFMETDQMAnalyzer.h.

Referenced by analyze(), and PFMETDQMAnalyzer().

edm::InputTag PFMETDQMAnalyzer::matchLabel_
private

Definition at line 25 of file PFMETDQMAnalyzer.h.

Referenced by analyze(), and PFMETDQMAnalyzer().

int PFMETDQMAnalyzer::nBadEvents_
private

Definition at line 33 of file PFMETDQMAnalyzer.h.

Referenced by analyze(), and beginJob().

PFMETMonitor PFMETDQMAnalyzer::pfMETMonitor_
private

Definition at line 29 of file PFMETDQMAnalyzer.h.

Referenced by analyze(), beginJob(), and PFMETDQMAnalyzer().

edm::ParameterSet PFMETDQMAnalyzer::pSet_
private

Definition at line 31 of file PFMETDQMAnalyzer.h.

Referenced by analyze(), beginJob(), and PFMETDQMAnalyzer().