CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PFCandidateDQMAnalyzer.cc
Go to the documentation of this file.
2 
4 
9 
12 
16 //
17 // -- Constructor
18 //
20 
21 {
23  inputLabel_ = pSet_.getParameter<edm::InputTag>("InputCollection");
24  matchLabel_ = pSet_.getParameter<edm::InputTag>("MatchCollection");
25  benchmarkLabel_ = pSet_.getParameter<std::string>("BenchmarkLabel");
26  createEfficiencyHistos_ = pSet_.getParameter<bool>( "CreateEfficiencyHistos" );
27 
28  pfCandidateMonitor_.setParameters(parameterSet);
29 
30  myCand_ = consumes< edm::View<reco::Candidate> >(inputLabel_);
31  myMatchedCand_ = consumes< edm::View<reco::Candidate> >(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("PFCandidateDQMAnalyzer") << " PFCandidateDQMAnalyzer::beginJob " << "Histogram Folder path set to "<< path;
44  nBadEvents_ = 0;
45 }
46 //
47 // -- Analyze
48 //
50  edm::EventSetup const& iSetup) {
51 
53  edm::Handle< edm::View<reco::Candidate> > matchedCandCollection;
54  if ( !createEfficiencyHistos_ ) {
55  iEvent.getByToken( myCand_, candCollection);
56  iEvent.getByToken( myMatchedCand_, matchedCandCollection);
57  } else {
58  iEvent.getByToken( myMatchedCand_, candCollection);
59  iEvent.getByToken( myCand_, matchedCandCollection);
60  }
61 
62  float maxRes = 0.0;
63  float minRes = 99.99;
64  if (candCollection.isValid() && matchedCandCollection.isValid()) {
65  //pfCandidateMonitor_.fill( *candCollection, *matchedCandCollection, minRes, maxRes);
66  pfCandidateMonitor_.fill( *candCollection, *matchedCandCollection, minRes, maxRes, pSet_);
67 
68  edm::ParameterSet skimPS = pSet_.getParameter<edm::ParameterSet>("SkimParameter");
69  if ( (skimPS.getParameter<bool>("switchOn")) &&
70  (nBadEvents_ <= skimPS.getParameter<int32_t>("maximumNumberToBeStored")) ) {
71  if ( minRes < skimPS.getParameter<double>("lowerCutOffOnResolution")) {
72  nBadEvents_++;
73  storeBadEvents(iEvent,minRes);
74  } else if (maxRes > skimPS.getParameter<double>("upperCutOffOnResolution")) {
75  nBadEvents_++;
76  storeBadEvents(iEvent,maxRes);
77  }
78  }
79 
80  }
81 }
83  unsigned int runNb = iEvent.id().run();
84  unsigned int evtNb = iEvent.id().event();
85  unsigned int lumiNb = iEvent.id().luminosityBlock();
86 
87  std::string path = "ParticleFlow/" + benchmarkLabel_ + "/BadEvents";
88  Benchmark::DQM_->setCurrentFolder(path.c_str());
89  std::ostringstream eventid_str;
90  eventid_str << runNb << "_"<< evtNb << "_" << lumiNb;
91  MonitorElement* me = Benchmark::DQM_->get(path + "/" + eventid_str.str());
92  if (me) me->Reset();
93  else me = Benchmark::DQM_->bookFloat(eventid_str.str());
94  me->Fill(val);
95 }
96 //
97 // -- EndJob
98 //
100 }
RunNumber_t run() const
Definition: EventID.h:42
T getParameter(std::string const &) const
EventNumber_t event() const
Definition: EventID.h:44
void setup()
book histograms
edm::EDGetTokenT< edm::View< reco::Candidate > > myMatchedCand_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
PFCandidateDQMAnalyzer(const edm::ParameterSet &parameterSet)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
LuminosityBlockNumber_t luminosityBlock() const
Definition: EventID.h:43
static DQMStore * DQM_
Definition: Benchmark.h:40
void setParameters(float dRMax, bool matchCharge, Benchmark::Mode mode, float ptmin, float ptmax, float etamin, float etamax, float phimin, float phimax, bool refHistoFlag)
set the parameters locally
MonitorElement * bookFloat(const char *name)
Book float.
Definition: DQMStore.cc:891
void Fill(long long x)
void analyze(edm::Event const &, edm::EventSetup const &)
void fill(const T &candidateCollection, const C &matchedCandCollection, float &minVal, float &maxVal, const edm::ParameterSet &parameterSet)
fill histograms with all particle
int iEvent
Definition: GenABIO.cc:230
tuple path
else: Piece not in the list, fine.
edm::EDGetTokenT< edm::View< reco::Candidate > > myCand_
MonitorElement * get(const std::string &path) const
get ME from full pathname (e.g. &quot;my/long/dir/my_histo&quot;)
Definition: DQMStore.cc:1708
bool isValid() const
Definition: HandleBase.h:76
PFCandidateMonitor pfCandidateMonitor_
edm::EventID id() const
Definition: EventBase.h:56
void storeBadEvents(edm::Event const &, float &val)
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:667