CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PFMETDQMAnalyzer.cc
Go to the documentation of this file.
2 
4 
9 
13 
17 //
18 // -- Constructor
19 //
21 
22 {
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  myMET_ = consumes< edm::View<reco::MET> >(inputLabel_);
31  myMatchedMET_ = consumes< edm::View<reco::MET> >(matchLabel_);
32 }
33 //
34 // -- BeginJob
35 //
37 
39  // part of the following could be put in the base class
40  std::string path = "ParticleFlow/" + benchmarkLabel_;
41  Benchmark::DQM_->setCurrentFolder(path.c_str());
42  edm::LogInfo("PFJMETDQMAnalyzer") << " PFMETDQMAnalyzer::beginJob " <<"Histogram Folder path set to "<< path;
44  nBadEvents_ = 0;
45 }
46 //
47 // -- Analyze
48 //
50  edm::EventSetup const& iSetup) {
51  edm::Handle< edm::View<reco::MET> > metCollection;
52  iEvent.getByToken(myMET_, metCollection);
53 
54  edm::Handle< edm::View<reco::MET> > matchedMetCollection;
55  iEvent.getByToken(myMatchedMET_, matchedMetCollection);
56 
57  if (metCollection.isValid() && matchedMetCollection.isValid()) {
58  float maxRes = 0.0;
59  float minRes = 99.99;
60  pfMETMonitor_.fillOne( (*metCollection)[0], (*matchedMetCollection)[0], minRes, maxRes);
61  //pfMETMonitor_.fillOne( (*metCollection)[0], (*matchedMetCollection)[0], minRes, maxRes, pSet_);
62  edm::ParameterSet skimPS = pSet_.getParameter<edm::ParameterSet>("SkimParameter");
63  if ( (skimPS.getParameter<bool>("switchOn")) &&
64  (nBadEvents_ <= skimPS.getParameter<int32_t>("maximumNumberToBeStored")) ) {
65  if ( minRes < skimPS.getParameter<double>("lowerCutOffOnResolution")) {
66  storeBadEvents(iEvent,minRes);
67  nBadEvents_++;
68  } else if (maxRes > skimPS.getParameter<double>("upperCutOffOnResolution")) {
69  nBadEvents_++;
70  storeBadEvents(iEvent,maxRes);
71  }
72  }
73  }
74 }
76  unsigned int runNb = iEvent.id().run();
77  unsigned int evtNb = iEvent.id().event();
78  unsigned int lumiNb = iEvent.id().luminosityBlock();
79 
80  std::string path = "ParticleFlow/" + benchmarkLabel_ + "/BadEvents";
81  Benchmark::DQM_->setCurrentFolder(path.c_str());
82  std::ostringstream eventid_str;
83  eventid_str << runNb << "_"<< evtNb << "_" << lumiNb;
84  MonitorElement* me = Benchmark::DQM_->get(path + "/" + eventid_str.str());
85  if (me) me->Reset();
86  else me = Benchmark::DQM_->bookFloat(eventid_str.str());
87  me->Fill(val);
88 }
89 
90 //
91 // -- EndJob
92 //
94 }
RunNumber_t run() const
Definition: EventID.h:42
T getParameter(std::string const &) const
EventNumber_t event() const
Definition: EventID.h:44
void analyze(edm::Event const &, edm::EventSetup 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:67
edm::EDGetTokenT< edm::View< reco::MET > > myMET_
edm::EDGetTokenT< edm::View< reco::MET > > myMatchedMET_
edm::InputTag inputLabel_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
edm::InputTag matchLabel_
PFMETMonitor pfMETMonitor_
LuminosityBlockNumber_t luminosityBlock() const
Definition: EventID.h:43
static DQMStore * DQM_
Definition: Benchmark.h:40
MonitorElement * bookFloat(const char *name)
Book float.
Definition: DQMStore.cc:809
void Fill(long long x)
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:1623
bool isValid() const
Definition: HandleBase.h:76
void fillOne(const reco::MET &met, const reco::MET &matchedMet, float &minVal, float &maxVal)
std::string benchmarkLabel_
void storeBadEvents(edm::Event const &, float &val)
edm::EventID id() const
Definition: EventBase.h:56
edm::ParameterSet pSet_
PFMETDQMAnalyzer(const edm::ParameterSet &parameterSet)
void setup()
book histograms
void Reset(void)
reset ME (ie. contents, errors, etc)
ParameterSet const & parameterSet(Provenance const &provenance)
Definition: Provenance.cc:11
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:584