CMS 3D CMS Logo

List of all members | Public Member Functions | Protected Attributes
PFCandidateMonitor Class Reference

#include <PFCandidateMonitor.h>

Inheritance diagram for PFCandidateMonitor:
Benchmark

Public Member Functions

template<class T , class C >
void fill (const T &candidateCollection, const C &matchedCandCollection, float &minVal, float &maxVal, const edm::ParameterSet &parameterSet)
 fill histograms with all particle More...
 
template<class T , class C , class M >
void fill (const T &candidateCollection, const C &matchedCandCollection, float &minVal, float &maxVal, const edm::ParameterSet &parameterSet, const M &muonMatchedCandCollection)
 
void fillOne (const reco::Candidate &cand)
 
 PFCandidateMonitor (float dRMax=0.3, bool matchCharge=true, Benchmark::Mode mode=Benchmark::DEFAULT)
 
void setDirectory (TDirectory *dir) override
 set directory (to use in ROOT) More...
 
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 More...
 
void setParameters (const edm::ParameterSet &parameterSet)
 set the parameters accessing them from ParameterSet More...
 
void setup (DQMStore::IBooker &b)
 book histograms More...
 
void setup (DQMStore::IBooker &b, const edm::ParameterSet &parameterSet)
 
 ~PFCandidateMonitor () override
 
- Public Member Functions inherited from Benchmark
 Benchmark (Mode mode=DEFAULT)
 
bool isInRange (float pt, float eta, float phi) const
 
void setParameters (Mode mode)
 
void setRange (float ptMin, float ptMax, float etaMin, float etaMax, float phiMin, float phiMax)
 
void write ()
 
virtual ~Benchmark ()(false)
 

Protected Attributes

CandidateBenchmark candBench_
 
bool createEfficiencyHistos_
 
bool createReferenceHistos_
 
TH1F * deltaR_
 
float dRMax_
 
TH1F * eta_gen_
 
TH1F * eta_ref_
 
bool histogramBooked_
 
MatchCandidateBenchmark matchCandBench_
 
bool matchCharge_
 
bool matching_done_
 
TH1F * phi_gen_
 
TH1F * phi_ref_
 
TH1F * pt_gen_
 
TH1F * pt_ref_
 
- Protected Attributes inherited from Benchmark
TDirectory * dir_
 
float etaMax_
 
float etaMin_
 
Mode mode_
 
float phiMax_
 
float phiMin_
 
float ptMax_
 
float ptMin_
 

Additional Inherited Members

- Public Types inherited from Benchmark
enum  Mode { DEFAULT, DQMOFFLINE, VALIDATION }
 
- Protected Member Functions inherited from Benchmark
TH1F * book1D (DQMStore::IBooker &b, const char *histname, const char *title, int nbins, float xmin, float xmax)
 book a 1D histogram, either through IBooker or plain root More...
 
TH2F * book2D (DQMStore::IBooker &b, const char *histname, const char *title, int nbinsx, float xmin, float xmax, int nbinsy, float ymin, float ymax)
 book a 2D histogram, either through IBooker or plain root More...
 
TH2F * book2D (DQMStore::IBooker &b, const char *histname, const char *title, int nbinsx, float *xbins, int nbinsy, float ymin, float ymax)
 book a 2D histogram, either through IBooker or plain root More...
 
TProfile * bookProfile (DQMStore::IBooker &b, const char *histname, const char *title, int nbinsx, float xmin, float xmax, float ymin, float ymax, const char *option)
 book a TProfile, either through IBooker or plain root More...
 
TProfile * bookProfile (DQMStore::IBooker &b, const char *histname, const char *title, int nbinsx, float *xbins, float ymin, float ymax, const char *option)
 book a TProfile, either through IBooker or plain root More...
 

Detailed Description

Definition at line 16 of file PFCandidateMonitor.h.

Constructor & Destructor Documentation

PFCandidateMonitor::PFCandidateMonitor ( float  dRMax = 0.3,
bool  matchCharge = true,
Benchmark::Mode  mode = Benchmark::DEFAULT 
)

Definition at line 15 of file PFCandidateMonitor.cc.

References createReferenceHistos_, deltaR_, eta_gen_, eta_ref_, histogramBooked_, phi_gen_, phi_ref_, pt_gen_, pt_ref_, and Benchmark::setRange().

17  setRange(0.0, 10e10, -10.0, 10.0, -3.14, 3.14);
18 
19  pt_gen_ = nullptr;
20  eta_gen_ = nullptr;
21  phi_gen_ = nullptr;
22 
23  pt_ref_ = nullptr;
24  eta_ref_ = nullptr;
25  phi_ref_ = nullptr;
26 
27  deltaR_ = nullptr;
28 
29  createReferenceHistos_ = false;
30  histogramBooked_ = false;
31 }
Benchmark(Mode mode=DEFAULT)
Definition: Benchmark.h:34
CandidateBenchmark candBench_
void setRange(float ptMin, float ptMax, float etaMin, float etaMax, float phiMin, float phiMax)
Definition: Benchmark.h:41
MatchCandidateBenchmark matchCandBench_
PFCandidateMonitor::~PFCandidateMonitor ( )
override

Definition at line 36 of file PFCandidateMonitor.cc.

36 {}

Member Function Documentation

template<class T , class C >
void PFCandidateMonitor::fill ( const T candidateCollection,
const C &  matchedCandCollection,
float &  minVal,
float &  maxVal,
const edm::ParameterSet parameterSet 
)

fill histograms with all particle

Definition at line 87 of file PFCandidateMonitor.h.

References candBench_, createEfficiencyHistos_, createReferenceHistos_, reco::deltaR(), deltaR_, dRMax_, PVValHelper::eta, reco::Candidate::eta(), CandidateBenchmark::fillOne(), MatchCandidateBenchmark::fillOne(), fillOne(), mps_fire::i, Benchmark::isInRange(), PFB::match(), matchCandBench_, matchCharge_, matching_done_, phi, reco::Candidate::phi(), EnergyCorrector::pt, reco::Candidate::pt(), and findQualityFiles::size.

Referenced by PFCandidateDQMAnalyzer::analyze(), and PFMuonDQMAnalyzer::analyze().

91  {
92  matching_done_ = false;
94  for (unsigned i = 0; i < candCollection.size(); ++i) {
95  if (!isInRange(candCollection[i].pt(), candCollection[i].eta(), candCollection[i].phi()))
96  continue;
97  fillOne(candCollection[i]); // fill pt_gen, eta_gen and phi_gen histos for
98  // UNMATCHED generated candidate
99 
100  for (unsigned j = 0; j < matchedCandCollection.size(); ++j) // for DeltaR spectrum
101  if (deltaR_)
102  deltaR_->Fill(reco::deltaR(candCollection[i], matchedCandCollection[j]));
103  }
104  }
105 
106  std::vector<int> matchIndices;
107  PFB::match(candCollection, matchedCandCollection, matchIndices, matchCharge_, dRMax_);
108  // PFB::match( candCollection, matchedCandCollection, matchIndices,
109  // parameterSet, matchCharge_, dRMax_ );
110  // now matchIndices[i] stores the j-th closest matched jet
111  matching_done_ = true;
112 
113  for (unsigned int i = 0; i < (candCollection).size(); i++) {
114  const reco::Candidate &cand = candCollection[i];
115 
116  if (!isInRange(cand.pt(), cand.eta(), cand.phi()))
117  continue;
118 
119  int iMatch = matchIndices[i];
120  assert(iMatch < static_cast<int>(matchedCandCollection.size()));
121 
122  if (iMatch != -1) {
123  const reco::Candidate &matchedCand = matchedCandCollection[iMatch];
124  if (!isInRange(matchedCand.pt(), matchedCand.eta(), matchedCand.phi()))
125  continue;
126  // std::cout <<"PFJet pT " <<cand.pt() <<" eta " <<cand.eta() <<" phi "
127  // <<cand.phi() ; std::cout <<"\nmatched genJet pT " <<matchedCand.pt()
128  // <<" eta " <<matchedCand.eta() <<" phi " <<matchedCand.phi() <<"\n"
129  // <<std::endl ;
130  float ptRes = (cand.pt() - matchedCand.pt()) / matchedCand.pt();
131 
132  if (ptRes > maxVal)
133  maxVal = ptRes;
134  if (ptRes < minVal)
135  minVal = ptRes;
136 
138  candBench_.fillOne(cand); // fill pt, eta phi and charge histos for MATCHED candidate
139  // matchCandBench_.fillOne(cand, matchedCand); // fill delta_x_VS_y
140  // histos for matched couple
141  matchCandBench_.fillOne(cand, matchedCand,
142  parameterSet); // fill delta_x_VS_y histos for matched couple
144  fillOne(matchedCand); // fill pt_ref, eta_ref and phi_ref histos for
145  // MATCHED reference candidate
146  } else {
147  candBench_.fillOne(matchedCand); // fill pt, eta phi and charge histos
148  // for MATCHED candidate
149  // matchCandBench_.fillOne(matchedCand, cand); // fill delta_x_VS_y
150  // histos for matched couple
151  matchCandBench_.fillOne(cand, matchedCand,
152  parameterSet); // fill delta_x_VS_y histos for matched couple
154  fillOne(cand); // fill pt_ref, eta_ref and phi_ref histos for MATCHED
155  // reference candidate
156  }
157  }
158  }
159 }
size
Write out results.
void fillOne(const reco::Candidate &candidate, const reco::Candidate &matchedCandidate)
fill histograms with a given particle
void match(const C &candCollection, const M &matchedCandCollection, std::vector< int > &matchIndices, bool matchCharge=false, float dRMax=-1)
Definition: Matchers.h:17
void fillOne(const reco::Candidate &candidate)
fill histograms with a given particle
constexpr auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:28
virtual double eta() const =0
momentum pseudorapidity
virtual double pt() const =0
transverse momentum
CandidateBenchmark candBench_
void fillOne(const reco::Candidate &cand)
MatchCandidateBenchmark matchCandBench_
bool isInRange(float pt, float eta, float phi) const
Definition: Benchmark.h:50
virtual double phi() const =0
momentum azimuthal angle
template<class T , class C , class M >
void PFCandidateMonitor::fill ( const T candidateCollection,
const C &  matchedCandCollection,
float &  minVal,
float &  maxVal,
const edm::ParameterSet parameterSet,
const M &  muonMatchedCandCollection 
)

Definition at line 164 of file PFCandidateMonitor.h.

References candBench_, createEfficiencyHistos_, createReferenceHistos_, reco::deltaR(), deltaR_, dRMax_, PVValHelper::eta, reco::Candidate::eta(), CandidateBenchmark::fillOne(), MatchCandidateBenchmark::fillOne(), fillOne(), mps_fire::i, Benchmark::isInRange(), PFB::match(), matchCandBench_, matchCharge_, matching_done_, phi, reco::Candidate::phi(), EnergyCorrector::pt, reco::Candidate::pt(), and findQualityFiles::size.

169  {
170  matching_done_ = false;
172  for (unsigned i = 0; i < candCollection.size(); ++i) {
173  if (!isInRange(candCollection[i].pt(), candCollection[i].eta(), candCollection[i].phi()))
174  continue;
175  fillOne(candCollection[i]); // fill pt_gen, eta_gen and phi_gen histos for
176  // UNMATCHED generated candidate
177 
178  for (unsigned j = 0; j < matchedCandCollection.size(); ++j) // for DeltaR spectrum
179  if (deltaR_)
180  deltaR_->Fill(reco::deltaR(candCollection[i], matchedCandCollection[j]));
181  }
182  }
183 
184  std::vector<int> matchIndices;
185  // PFB::match( candCollection, matchedCandCollection, matchIndices,
186  // matchCharge_, dRMax_ ); PFB::match( candCollection, matchedCandCollection,
187  // matchIndices, parameterSet, matchCharge_, dRMax_ );
188  PFB::match(candCollection,
189  matchedCandCollection,
190  matchIndices,
191  parameterSet,
192  muonMatchedCandCollection,
193  matchCharge_,
194  dRMax_);
195  // now matchIndices[i] stores the j-th closest matched jet
196  matching_done_ = true;
197 
198  for (unsigned int i = 0; i < (candCollection).size(); i++) {
199  const reco::Candidate &cand = candCollection[i];
200 
201  if (!isInRange(cand.pt(), cand.eta(), cand.phi()))
202  continue;
203 
204  int iMatch = matchIndices[i];
205  assert(iMatch < static_cast<int>(matchedCandCollection.size()));
206 
207  if (iMatch != -1) {
208  const reco::Candidate &matchedCand = matchedCandCollection[iMatch];
209  if (!isInRange(matchedCand.pt(), matchedCand.eta(), matchedCand.phi()))
210  continue;
211  // std::cout <<"PFJet pT " <<cand.pt() <<" eta " <<cand.eta() <<" phi "
212  // <<cand.phi() ; std::cout <<"\nmatched genJet pT " <<matchedCand.pt()
213  // <<" eta " <<matchedCand.eta() <<" phi " <<matchedCand.phi() <<"\n"
214  // <<std::endl ;
215  float ptRes = (cand.pt() - matchedCand.pt()) / matchedCand.pt();
216 
217  if (ptRes > maxVal)
218  maxVal = ptRes;
219  if (ptRes < minVal)
220  minVal = ptRes;
221 
223  candBench_.fillOne(cand); // fill pt, eta phi and charge histos for MATCHED candidate
224  matchCandBench_.fillOne(cand, matchedCand,
225  parameterSet); // fill delta_x_VS_y histos for matched couple
227  fillOne(matchedCand); // fill pt_ref, eta_ref and phi_ref histos for
228  // MATCHED reference candidate
229  } else {
230  candBench_.fillOne(matchedCand); // fill pt, eta phi and charge histos
231  // for MATCHED candidate
232  matchCandBench_.fillOne(cand, matchedCand,
233  parameterSet); // fill delta_x_VS_y histos for matched couple
235  fillOne(cand); // fill pt_ref, eta_ref and phi_ref histos for MATCHED
236  // reference candidate
237  }
238  }
239  }
240 }
size
Write out results.
void fillOne(const reco::Candidate &candidate, const reco::Candidate &matchedCandidate)
fill histograms with a given particle
void match(const C &candCollection, const M &matchedCandCollection, std::vector< int > &matchIndices, bool matchCharge=false, float dRMax=-1)
Definition: Matchers.h:17
void fillOne(const reco::Candidate &candidate)
fill histograms with a given particle
constexpr auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:28
virtual double eta() const =0
momentum pseudorapidity
virtual double pt() const =0
transverse momentum
CandidateBenchmark candBench_
void fillOne(const reco::Candidate &cand)
MatchCandidateBenchmark matchCandBench_
bool isInRange(float pt, float eta, float phi) const
Definition: Benchmark.h:50
virtual double phi() const =0
momentum azimuthal angle
void PFCandidateMonitor::fillOne ( const reco::Candidate cand)

Definition at line 205 of file PFCandidateMonitor.cc.

References createEfficiencyHistos_, createReferenceHistos_, reco::Candidate::eta(), eta_gen_, eta_ref_, histogramBooked_, matching_done_, reco::Candidate::phi(), phi_gen_, phi_ref_, reco::Candidate::pt(), pt_gen_, and pt_ref_.

Referenced by fill().

205  {
206  if (matching_done_) {
208  if (pt_ref_)
209  pt_ref_->Fill(cand.pt());
210  if (eta_ref_)
211  eta_ref_->Fill(cand.eta());
212  if (phi_ref_)
213  phi_ref_->Fill(cand.phi());
214  }
216  if (pt_gen_)
217  pt_gen_->Fill(cand.pt());
218  if (eta_gen_)
219  eta_gen_->Fill(cand.eta());
220  if (phi_gen_)
221  phi_gen_->Fill(cand.phi());
222  }
223 }
virtual double eta() const =0
momentum pseudorapidity
virtual double pt() const =0
transverse momentum
virtual double phi() const =0
momentum azimuthal angle
void PFCandidateMonitor::setDirectory ( TDirectory *  dir)
overridevirtual

set directory (to use in ROOT)

Reimplemented from Benchmark.

Definition at line 195 of file PFCandidateMonitor.cc.

References candBench_, matchCandBench_, and Benchmark::setDirectory().

195  {
197 
200 }
virtual void setDirectory(TDirectory *dir)
Definition: Benchmark.cc:14
CandidateBenchmark candBench_
MatchCandidateBenchmark matchCandBench_
dbl *** dir
Definition: mlp_gen.cc:35
void PFCandidateMonitor::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

Definition at line 62 of file PFCandidateMonitor.cc.

References candBench_, createReferenceHistos_, allElectronIsolations_cfi::dRMax, dRMax_, matchCandBench_, matchCharge_, ALCARECOPromptCalibProdSiPixelAli0T_cff::mode, Benchmark::mode_, Benchmark::setParameters(), and Benchmark::setRange().

Referenced by PFCandidateDQMAnalyzer::PFCandidateDQMAnalyzer(), and PFMuonDQMAnalyzer::PFMuonDQMAnalyzer().

71  {
72  dRMax_ = dRMax;
73  matchCharge_ = matchCharge;
74  mode_ = mode;
75  createReferenceHistos_ = refHistoFlag;
76 
77  setRange(ptmin, ptmax, etamin, etamax, phimin, phimax);
78 
81 }
void setParameters(Mode mode)
Definition: Benchmark.h:39
CandidateBenchmark candBench_
void setRange(float ptMin, float ptMax, float etaMin, float etaMax, float phiMin, float phiMax)
Definition: Benchmark.h:41
double ptmin
Definition: HydjetWrapper.h:90
MatchCandidateBenchmark matchCandBench_
Mode mode_
Definition: Benchmark.h:124
void PFCandidateMonitor::setParameters ( const edm::ParameterSet parameterSet)

set the parameters accessing them from ParameterSet

Definition at line 41 of file PFCandidateMonitor.cc.

References candBench_, createEfficiencyHistos_, createReferenceHistos_, dRMax_, edm::ParameterSet::getParameter(), matchCandBench_, matchCharge_, Benchmark::mode_, Benchmark::setParameters(), and Benchmark::setRange().

41  {
42  dRMax_ = parameterSet.getParameter<double>("deltaRMax");
43  matchCharge_ = parameterSet.getParameter<bool>("matchCharge");
44  mode_ = (Benchmark::Mode)parameterSet.getParameter<int>("mode");
45  createReferenceHistos_ = parameterSet.getParameter<bool>("CreateReferenceHistos");
46  createEfficiencyHistos_ = parameterSet.getParameter<bool>("CreateEfficiencyHistos");
47 
48  setRange(parameterSet.getParameter<double>("ptMin"),
49  parameterSet.getParameter<double>("ptMax"),
50  parameterSet.getParameter<double>("etaMin"),
51  parameterSet.getParameter<double>("etaMax"),
52  parameterSet.getParameter<double>("phiMin"),
53  parameterSet.getParameter<double>("phiMax"));
54 
57 }
T getParameter(std::string const &) const
void setParameters(Mode mode)
Definition: Benchmark.h:39
CandidateBenchmark candBench_
void setRange(float ptMin, float ptMax, float etaMin, float etaMax, float phiMin, float phiMax)
Definition: Benchmark.h:41
MatchCandidateBenchmark matchCandBench_
Mode mode_
Definition: Benchmark.h:124
void PFCandidateMonitor::setup ( DQMStore::IBooker b)

book histograms

Definition at line 164 of file PFCandidateMonitor.cc.

References Benchmark::book1D(), candBench_, createEfficiencyHistos_, createReferenceHistos_, eta_gen_, eta_ref_, histogramBooked_, Benchmark::PhaseSpace::m, Benchmark::PhaseSpace::M, matchCandBench_, Benchmark::PhaseSpace::n, phi_gen_, phi_ref_, pt_gen_, pt_ref_, CandidateBenchmark::setup(), and MatchCandidateBenchmark::setup().

Referenced by PFCandidateDQMAnalyzer::bookHistograms(), and PFMuonDQMAnalyzer::bookHistograms().

164  {
165  candBench_.setup(b);
167 
169  PhaseSpace ptPS(100, 0, 100);
170  PhaseSpace phiPS(360, -3.1416, 3.1416);
171  PhaseSpace etaPS(100, -5, 5);
172 
173  pt_ref_ = book1D(b, "pt_ref_", "p_{T}_ref;p_{T} (GeV)", ptPS.n, ptPS.m, ptPS.M);
175  pt_gen_ = book1D(b, "pt_gen_", "p_{T}_gen;p_{T} (GeV)", ptPS.n, ptPS.m, ptPS.M);
176  }
177 
178  eta_ref_ = book1D(b, "eta_ref_", "#eta_ref;#eta", etaPS.n, etaPS.m, etaPS.M);
180  eta_gen_ = book1D(b, "eta_gen_", "#eta_gen;#eta", etaPS.n, etaPS.m, etaPS.M);
181  }
182 
183  phi_ref_ = book1D(b, "phi_ref_", "#phi_ref;#phi", phiPS.n, phiPS.m, phiPS.M);
185  phi_gen_ = book1D(b, "phi_gen_", "#phi_gen;#phi", phiPS.n, phiPS.m, phiPS.M);
186  }
187 
188  histogramBooked_ = true;
189  }
190 }
TH1F * book1D(DQMStore::IBooker &b, const char *histname, const char *title, int nbins, float xmin, float xmax)
book a 1D histogram, either through IBooker or plain root
Definition: Benchmark.cc:16
void setup(DQMStore::IBooker &b)
book histograms
CandidateBenchmark candBench_
MatchCandidateBenchmark matchCandBench_
void setup(DQMStore::IBooker &b)
book histograms
void PFCandidateMonitor::setup ( DQMStore::IBooker b,
const edm::ParameterSet parameterSet 
)

Definition at line 86 of file PFCandidateMonitor.cc.

References Benchmark::book1D(), candBench_, createEfficiencyHistos_, createReferenceHistos_, deltaR_, eta_gen_, eta_ref_, edm::ParameterSet::getParameter(), histogramBooked_, matchCandBench_, phi_gen_, phi_ref_, pt_gen_, pt_ref_, CandidateBenchmark::setup(), and MatchCandidateBenchmark::setup().

86  {
87  candBench_.setup(b, parameterSet);
88  matchCandBench_.setup(b, parameterSet);
89 
91  edm::ParameterSet ptPS = parameterSet.getParameter<edm::ParameterSet>("PtHistoParameter");
92  edm::ParameterSet etaPS = parameterSet.getParameter<edm::ParameterSet>("EtaHistoParameter");
93  edm::ParameterSet phiPS = parameterSet.getParameter<edm::ParameterSet>("PhiHistoParameter");
94 
95  edm::ParameterSet dR = parameterSet.getParameter<edm::ParameterSet>("DeltaRHistoParameter");
96 
97  if (ptPS.getParameter<bool>("switchOn")) {
98  pt_ref_ = book1D(b,
99  "pt_ref_",
100  "p_{T}_ref;p_{T} (GeV)",
101  ptPS.getParameter<int32_t>("nBin"),
102  ptPS.getParameter<double>("xMin"),
103  ptPS.getParameter<double>("xMax"));
105  pt_gen_ = book1D(b,
106  "pt_gen_",
107  "p_{T}_gen;p_{T} (GeV)",
108  ptPS.getParameter<int32_t>("nBin"),
109  ptPS.getParameter<double>("xMin"),
110  ptPS.getParameter<double>("xMax"));
111  }
112  }
113 
114  if (etaPS.getParameter<bool>("switchOn")) {
115  eta_ref_ = book1D(b,
116  "eta_ref_",
117  "#eta_ref;#eta",
118  etaPS.getParameter<int32_t>("nBin"),
119  etaPS.getParameter<double>("xMin"),
120  etaPS.getParameter<double>("xMax"));
122  eta_gen_ = book1D(b,
123  "eta_gen_",
124  "#eta_gen;#eta",
125  etaPS.getParameter<int32_t>("nBin"),
126  etaPS.getParameter<double>("xMin"),
127  etaPS.getParameter<double>("xMax"));
128  }
129  }
130 
131  if (phiPS.getParameter<bool>("switchOn")) {
132  phi_ref_ = book1D(b,
133  "phi_ref_",
134  "#phi_ref;#phi",
135  phiPS.getParameter<int32_t>("nBin"),
136  phiPS.getParameter<double>("xMin"),
137  phiPS.getParameter<double>("xMax"));
139  phi_gen_ = book1D(b,
140  "phi_gen_",
141  "#phi_gen;#phi",
142  phiPS.getParameter<int32_t>("nBin"),
143  phiPS.getParameter<double>("xMin"),
144  phiPS.getParameter<double>("xMax"));
145  }
146  }
147 
148  if (createEfficiencyHistos_ && dR.getParameter<bool>("switchOn")) {
149  deltaR_ = book1D(b,
150  "deltaR_",
151  "#DeltaR;#DeltaR",
152  dR.getParameter<int32_t>("nBin"),
153  dR.getParameter<double>("xMin"),
154  dR.getParameter<double>("xMax"));
155  }
156 
157  histogramBooked_ = true;
158  }
159 }
T getParameter(std::string const &) const
TH1F * book1D(DQMStore::IBooker &b, const char *histname, const char *title, int nbins, float xmin, float xmax)
book a 1D histogram, either through IBooker or plain root
Definition: Benchmark.cc:16
void setup(DQMStore::IBooker &b)
book histograms
CandidateBenchmark candBench_
MatchCandidateBenchmark matchCandBench_
void setup(DQMStore::IBooker &b)
book histograms

Member Data Documentation

CandidateBenchmark PFCandidateMonitor::candBench_
protected

Definition at line 64 of file PFCandidateMonitor.h.

Referenced by fill(), setDirectory(), setParameters(), and setup().

bool PFCandidateMonitor::createEfficiencyHistos_
protected

Definition at line 82 of file PFCandidateMonitor.h.

Referenced by fill(), fillOne(), setParameters(), and setup().

bool PFCandidateMonitor::createReferenceHistos_
protected

Definition at line 78 of file PFCandidateMonitor.h.

Referenced by fill(), fillOne(), PFCandidateMonitor(), setParameters(), and setup().

TH1F* PFCandidateMonitor::deltaR_
protected

Definition at line 75 of file PFCandidateMonitor.h.

Referenced by fill(), PFCandidateMonitor(), and setup().

float PFCandidateMonitor::dRMax_
protected

Definition at line 76 of file PFCandidateMonitor.h.

Referenced by fill(), and setParameters().

TH1F* PFCandidateMonitor::eta_gen_
protected

Definition at line 68 of file PFCandidateMonitor.h.

Referenced by fillOne(), PFCandidateMonitor(), and setup().

TH1F* PFCandidateMonitor::eta_ref_
protected

Definition at line 72 of file PFCandidateMonitor.h.

Referenced by fillOne(), PFCandidateMonitor(), and setup().

bool PFCandidateMonitor::histogramBooked_
protected

Definition at line 79 of file PFCandidateMonitor.h.

Referenced by fillOne(), PFCandidateMonitor(), and setup().

MatchCandidateBenchmark PFCandidateMonitor::matchCandBench_
protected

Definition at line 65 of file PFCandidateMonitor.h.

Referenced by fill(), setDirectory(), setParameters(), and setup().

bool PFCandidateMonitor::matchCharge_
protected

Definition at line 77 of file PFCandidateMonitor.h.

Referenced by fill(), and setParameters().

bool PFCandidateMonitor::matching_done_
protected

Definition at line 81 of file PFCandidateMonitor.h.

Referenced by fill(), and fillOne().

TH1F* PFCandidateMonitor::phi_gen_
protected

Definition at line 69 of file PFCandidateMonitor.h.

Referenced by fillOne(), PFCandidateMonitor(), and setup().

TH1F* PFCandidateMonitor::phi_ref_
protected

Definition at line 73 of file PFCandidateMonitor.h.

Referenced by fillOne(), PFCandidateMonitor(), and setup().

TH1F* PFCandidateMonitor::pt_gen_
protected

Definition at line 67 of file PFCandidateMonitor.h.

Referenced by fillOne(), PFCandidateMonitor(), and setup().

TH1F* PFCandidateMonitor::pt_ref_
protected

Definition at line 71 of file PFCandidateMonitor.h.

Referenced by fillOne(), PFCandidateMonitor(), and setup().