CMS 3D CMS Logo

PFCandidateMonitor.cc
Go to the documentation of this file.
6 
7 #include <TFile.h>
8 #include <TH1.h>
9 #include <TH2.h>
10 #include <TROOT.h>
11 
12 //
13 // -- Constructor
14 //
16  : Benchmark(mode), candBench_(mode), matchCandBench_(mode), dRMax_(dRMax), matchCharge_(matchCharge) {
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 }
32 
33 //
34 // -- Destructor
35 //
37 
38 //
39 // -- Set Parameters accessing them from ParameterSet
40 //
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 }
58 
59 //
60 // -- Set Parameters
61 //
63  bool matchCharge,
65  float ptmin,
66  float ptmax,
67  float etamin,
68  float etamax,
69  float phimin,
70  float phimax,
71  bool refHistoFlag) {
72  dRMax_ = dRMax;
73  matchCharge_ = matchCharge;
74  mode_ = mode;
75  createReferenceHistos_ = refHistoFlag;
76 
77  setRange(ptmin, ptmax, etamin, etamax, phimin, phimax);
78 
81 }
82 
83 //
84 // -- Create histograms accessing parameters from ParameterSet
85 //
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 }
160 
161 //
162 // -- Create histograms using local parameters
163 //
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 }
191 
192 //
193 // -- Set directory to book histograms using ROOT
194 //
197 
200 }
201 
202 //
203 // -- fill histograms for a single collection
204 //
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 }
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
abstract base class
Definition: Benchmark.h:21
void setDirectory(TDirectory *dir) override
set directory (to use in ROOT)
void setup(DQMStore::IBooker &b)
book histograms
virtual void setDirectory(TDirectory *dir)
Definition: Benchmark.cc:14
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
void setParameters(Mode mode)
Definition: Benchmark.h:39
~PFCandidateMonitor() override
virtual double eta() const =0
momentum pseudorapidity
virtual double pt() const =0
transverse momentum
CandidateBenchmark candBench_
void setRange(float ptMin, float ptMax, float etaMin, float etaMax, float phiMin, float phiMax)
Definition: Benchmark.h:41
void fillOne(const reco::Candidate &cand)
double b
Definition: hdecay.h:120
double ptmin
Definition: HydjetWrapper.h:90
MatchCandidateBenchmark matchCandBench_
Mode mode_
Definition: Benchmark.h:124
void setup(DQMStore::IBooker &b)
book histograms
dbl *** dir
Definition: mlp_gen.cc:35
virtual double phi() const =0
momentum azimuthal angle
ParameterSet const & parameterSet(Provenance const &provenance)
Definition: Provenance.cc:11
PFCandidateMonitor(float dRMax=0.3, bool matchCharge=true, Benchmark::Mode mode=Benchmark::DEFAULT)