CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PFJetDQMAnalyzer.cc
Go to the documentation of this file.
2 
4 
9 
15 
19 //
20 // -- Constructor
21 //
23 
24 {
26  inputLabel_ = pSet_.getParameter<edm::InputTag>("InputCollection");
27  matchLabel_ = pSet_.getParameter<edm::InputTag>("MatchCollection");
28  benchmarkLabel_ = pSet_.getParameter<std::string>("BenchmarkLabel");
29 
30  pfJetMonitor_.setParameters(parameterSet); // set parameters for booking histograms and validating jet
31 
32  myJet_ = consumes< edm::View<reco::Jet> >(inputLabel_);
33  myMatchedJet_ = consumes< edm::View<reco::Jet> >(matchLabel_);
34 }
35 //
36 // -- BeginJob
37 //
39 
41  // part of the following could be put in the base class
42  std::string path = "ParticleFlow/" + benchmarkLabel_;
43  Benchmark::DQM_->setCurrentFolder(path.c_str());
44  edm::LogInfo("PFJetDQMAnalyzer") << " PFJetDQMAnalyzer::beginJob " << "Histogram Folder path set to "<< path;
45  pfJetMonitor_.setup(pSet_); // booking histograms of type delta_frac_VS_frac from PFJetMonitor, pt_ eta_ phi_ and charge_ from CandidateBenchmark, delta_x_VS_y from MatchCandidateBenchmark
46  nBadEvents_ = 0;
47 }
48 //
49 // -- Analyze
50 //
52  edm::EventSetup const& iSetup) {
53 
54  edm::Handle< edm::View<reco::Jet> > jetCollection;
55  iEvent.getByToken(myJet_, jetCollection);
56 
57  edm::Handle< edm::View<reco::Jet> > matchedJetCollection;
58  iEvent.getByToken(myMatchedJet_, matchedJetCollection);
59 
60  float maxRes = 0.0;
61  float minRes = 99.99;
62  float jetpT = 0.0;
63  if (jetCollection.isValid() && matchedJetCollection.isValid()) {
64  //pfJetMonitor_.fill( *jetCollection, *matchedJetCollection, minRes, maxRes); // match collections and fill pt eta phi and charge histos for candidate jet, fill delta_x_VS_y histos for matched couples, book and fill delta_frac_VS_frac histos for matched couples
65  pfJetMonitor_.fill( *jetCollection, *matchedJetCollection, minRes, maxRes, jetpT, pSet_); // match collections and fill pt eta phi and charge histos for candidate jet, fill delta_x_VS_y histos for matched couples, book and fill delta_frac_VS_frac histos for matched couples
66  edm::ParameterSet skimPS = pSet_.getParameter<edm::ParameterSet>("SkimParameter");
67  if ( (skimPS.getParameter<bool>("switchOn")) &&
68  (nBadEvents_ <= skimPS.getParameter<int32_t>("maximumNumberToBeStored")) ) {
69  if (jetpT > skimPS.getParameter<double>("minimumJetpT")) {
70  if ( minRes < skimPS.getParameter<double>("lowerCutOffOnResolution")) {
71  storeBadEvents(iEvent,minRes);
72  nBadEvents_++;
73  } else if (maxRes > skimPS.getParameter<double>("upperCutOffOnResolution")) {
74  storeBadEvents(iEvent,maxRes);
75  nBadEvents_++;
76  }
77  } // minimum jet pT check
78  }
79  }
80 }
82  unsigned int runNb = iEvent.id().run();
83  unsigned int evtNb = iEvent.id().event();
84  unsigned int lumiNb = iEvent.id().luminosityBlock();
85 
86  std::string path = "ParticleFlow/" + benchmarkLabel_ + "/BadEvents";
87  Benchmark::DQM_->setCurrentFolder(path.c_str());
88  std::ostringstream eventid_str;
89  eventid_str << runNb << "_"<< evtNb << "_" << lumiNb;
90  MonitorElement* me = Benchmark::DQM_->get(path + "/" + eventid_str.str());
91  if (me) me->Reset();
92  else me = Benchmark::DQM_->bookFloat(eventid_str.str());
93  me->Fill(val);
94 }
95 //
96 // -- EndJob
97 //
99 }
RunNumber_t run() const
Definition: EventID.h:42
T getParameter(std::string const &) const
EventNumber_t event() const
Definition: EventID.h:44
void storeBadEvents(edm::Event const &, float &val)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
edm::InputTag inputLabel_
PFJetMonitor pfJetMonitor_
void fill(const T &jetCollection, const C &matchedJetCollection, float &minVal, float &maxVal)
fill histograms with all particle
Definition: PFJetMonitor.h:81
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:891
void Fill(long long x)
int iEvent
Definition: GenABIO.cc:230
tuple path
else: Piece not in the list, fine.
void setup()
book histograms
void setParameters(float dRMax, bool matchCharge, Benchmark::Mode mode, float ptmin, float ptmax, float etamin, float etamax, float phimin, float phimax, bool fracHistoFlag=true)
set the parameters locally
Definition: PFJetMonitor.cc:64
std::string benchmarkLabel_
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
edm::EDGetTokenT< edm::View< reco::Jet > > myMatchedJet_
edm::ParameterSet pSet_
edm::InputTag matchLabel_
edm::EventID id() const
Definition: EventBase.h:56
edm::EDGetTokenT< edm::View< reco::Jet > > myJet_
void analyze(edm::Event const &, edm::EventSetup const &)
void Reset(void)
reset ME (ie. contents, errors, etc)
PFJetDQMAnalyzer(const edm::ParameterSet &parameterSet)
ParameterSet const & parameterSet(Provenance const &provenance)
Definition: Provenance.cc:11
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:667