CMS 3D CMS Logo

PFMuonDQMAnalyzer.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  createEfficiencyHistos_ = pSet_.getParameter<bool>( "CreateEfficiencyHistos" );
28 
29  pfCandidateMonitor_.setParameters(parameterSet);
30 
31  myCand_ = consumes< edm::View<reco::Muon> >(inputLabel_);
32  myMatchedCand_ = consumes< edm::View<reco::Muon> >(matchLabel_);
33 
34 
36 
37  subsystemname_ = "ParticleFlow" ;
39 
40  nBadEvents_ = 0;
41 
42 }
43 
44 
45 //
46 // -- BookHistograms
47 //
49  edm::Run const & /* iRun */,
50  edm::EventSetup const & /* iSetup */ )
51 {
53 
54  edm::LogInfo("PFMuonDQMAnalyzer") << " PFMuonDQMAnalyzer::bookHistograms " << "Histogram Folder path set to " << eventInfoFolder_;
55 
57 }
58 
59 
60 //
61 // -- Analyze
62 //
64  edm::EventSetup const& iSetup) {
65 
66  edm::Handle< edm::View<reco::Muon> > candCollection;
67  edm::Handle< edm::View<reco::Muon> > matchedCandCollection;
68  if ( !createEfficiencyHistos_ ) {
69  iEvent.getByToken( myCand_, candCollection);
70  iEvent.getByToken( myMatchedCand_, matchedCandCollection);
71  } else {
72  iEvent.getByToken( myMatchedCand_, candCollection);
73  iEvent.getByToken( myCand_, matchedCandCollection);
74  }
75 
76  float maxRes = 0.0;
77  float minRes = 99.99;
78  if (candCollection.isValid() && matchedCandCollection.isValid()) {
79  pfCandidateMonitor_.fill( *candCollection, *matchedCandCollection, minRes, maxRes, pSet_, *matchedCandCollection) ;
80 
81  }
82 }
83 
84 
T getParameter(std::string const &) const
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:579
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
edm::InputTag inputLabel_
edm::ParameterSet pSet_
void setup(DQMStore::IBooker &b)
book histograms
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
PFCandidateMonitor pfCandidateMonitor_
std::string subsystemname_
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
edm::EDGetTokenT< edm::View< reco::Muon > > myMatchedCand_
void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:268
std::string eventInfoFolder_
PFMuonDQMAnalyzer(const edm::ParameterSet &parameterSet)
bool isValid() const
Definition: HandleBase.h:74
void analyze(edm::Event const &, edm::EventSetup const &) override
edm::InputTag matchLabel_
std::string benchmarkLabel_
ParameterSet const & parameterSet(Provenance const &provenance)
Definition: Provenance.cc:11
edm::EDGetTokenT< edm::View< reco::Muon > > myCand_
Definition: Run.h:44
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override