CMS 3D CMS Logo

PFClient.cc
Go to the documentation of this file.
2 
7 
11 
12 //
13 // -- Constructor
14 //
16  folderNames_ = parameterSet.getParameter<std::vector<std::string>>("FolderNames");
17  histogramNames_ = parameterSet.getParameter<std::vector<std::string>>("HistogramNames");
18  efficiencyFlag_ = parameterSet.getParameter<bool>("CreateEfficiencyPlots");
19  effHistogramNames_ = parameterSet.getParameter<std::vector<std::string>>("HistogramNamesForEfficiencyPlots");
20  projectionHistogramNames_ = parameterSet.getParameter<std::vector<std::string>>("HistogramNamesForProjectionPlots");
21  profileFlag_ = parameterSet.getParameter<bool>("CreateProfilePlots");
22  profileHistogramNames_ = parameterSet.getParameter<std::vector<std::string>>("HistogramNamesForProfilePlots");
23 }
24 
25 //
26 // -- EndJobBegin Run
27 //
29  doSummaries(ibooker, igetter);
30  doProjection(ibooker, igetter);
31  if (efficiencyFlag_)
32  doEfficiency(ibooker, igetter);
33  if (profileFlag_)
34  doProfiles(ibooker, igetter);
35 }
36 
37 //
38 // -- Create Summaries
39 //
41  for (std::vector<std::string>::const_iterator ifolder = folderNames_.begin(); ifolder != folderNames_.end();
42  ifolder++) {
43  std::string path = "ParticleFlow/" + (*ifolder);
44 
45  for (std::vector<std::string>::const_iterator ihist = histogramNames_.begin(); ihist != histogramNames_.end();
46  ihist++) {
47  std::string hname = (*ihist);
48  createResolutionPlots(ibooker, igetter, path, hname);
49  }
50  }
51 }
52 
53 //
54 // -- Create Projection
55 //
57  for (std::vector<std::string>::const_iterator ifolder = folderNames_.begin(); ifolder != folderNames_.end();
58  ifolder++) {
59  std::string path = "ParticleFlow/" + (*ifolder);
60 
61  for (std::vector<std::string>::const_iterator ihist = projectionHistogramNames_.begin();
62  ihist != projectionHistogramNames_.end();
63  ihist++) {
64  std::string hname = (*ihist);
65  createProjectionPlots(ibooker, igetter, path, hname);
66  }
67  }
68 }
69 
70 //
71 // -- Create Profile
72 //
74  for (std::vector<std::string>::const_iterator ifolder = folderNames_.begin(); ifolder != folderNames_.end();
75  ifolder++) {
76  std::string path = "ParticleFlow/" + (*ifolder);
77 
78  for (std::vector<std::string>::const_iterator ihist = profileHistogramNames_.begin();
79  ihist != profileHistogramNames_.end();
80  ihist++) {
81  std::string hname = (*ihist);
82  createProfilePlots(ibooker, igetter, path, hname);
83  }
84  }
85 }
86 
87 //
88 // -- Create Efficiency
89 //
91  for (std::vector<std::string>::const_iterator ifolder = folderNames_.begin(); ifolder != folderNames_.end();
92  ifolder++) {
93  std::string path = "ParticleFlow/" + (*ifolder);
94 
95  for (std::vector<std::string>::const_iterator ihist = effHistogramNames_.begin(); ihist != effHistogramNames_.end();
96  ihist++) {
97  std::string hname = (*ihist);
98  createEfficiencyPlots(ibooker, igetter, path, hname);
99  }
100  }
101 }
102 
103 //
104 // -- Create Resolution Plots
105 //
107  DQMStore::IGetter &igetter,
109  std::string &name) {
110  MonitorElement *me = igetter.get(folder + "/" + name);
111 
112  if (!me)
113  return;
114 
115  MonitorElement *me_average;
116  MonitorElement *me_rms;
117  MonitorElement *me_mean;
118  MonitorElement *me_sigma;
119 
121  (me->kind() == MonitorElement::DQM_KIND_TH2D)) {
122  TH2 *th = me->getTH2F();
123  size_t nbinx = me->getNbinsX();
124  size_t nbiny = me->getNbinsY();
125 
126  float ymin = th->GetYaxis()->GetXmin();
127  float ymax = th->GetYaxis()->GetXmax();
128  std::string xtit = th->GetXaxis()->GetTitle();
129  std::string ytit = th->GetYaxis()->GetTitle();
130  float *xbins = new float[nbinx + 1];
131  for (size_t ix = 1; ix < nbinx + 1; ++ix) {
132  xbins[ix - 1] = th->GetXaxis()->GetBinLowEdge(ix);
133  if (ix == nbinx)
134  xbins[ix] = th->GetXaxis()->GetBinUpEdge(ix);
135  }
136  std::string tit_new = ";" + xtit + ";" + ytit;
137  ibooker.setCurrentFolder(folder);
138  MonitorElement *me_slice = ibooker.book1D("PFlowSlice", "PFlowSlice", nbiny, ymin, ymax);
139 
140  tit_new = ";" + xtit + ";Average_" + ytit;
141  me_average = ibooker.book1D("average_" + name, tit_new, nbinx, xbins);
142  me_average->setEfficiencyFlag();
143  tit_new = ";" + xtit + ";RMS_" + ytit;
144  me_rms = ibooker.book1D("rms_" + name, tit_new, nbinx, xbins);
145  me_rms->setEfficiencyFlag();
146  tit_new = ";" + xtit + ";Mean_" + ytit;
147  me_mean = ibooker.book1D("mean_" + name, tit_new, nbinx, xbins);
148  me_mean->setEfficiencyFlag();
149  tit_new = ";" + xtit + ";Sigma_" + ytit;
150  me_sigma = ibooker.book1D("sigma_" + name, tit_new, nbinx, xbins);
151  me_sigma->setEfficiencyFlag();
152 
153  double average, rms, mean, sigma;
154 
155  for (size_t ix = 1; ix < nbinx + 1; ++ix) {
156  me_slice->Reset();
157  for (size_t iy = 1; iy < nbiny + 1; ++iy) {
158  me_slice->setBinContent(iy, th->GetBinContent(ix, iy));
159  }
160  getHistogramParameters(me_slice, average, rms, mean, sigma);
161  me_average->setBinContent(ix, average);
162  me_rms->setBinContent(ix, rms);
163  me_mean->setBinContent(ix, mean);
164  me_sigma->setBinContent(ix, sigma);
165  }
166  if (me_slice)
167  igetter.removeElement(me_slice->getName()); //?
168  delete[] xbins;
169  }
170 }
171 
172 //
173 // -- Create Projection Plots
174 //
176  DQMStore::IGetter &igetter,
178  std::string &name) {
179  MonitorElement *me = igetter.get(folder + "/" + name);
180  if (!me)
181  return;
182 
183  MonitorElement *projection = nullptr;
184 
186  (me->kind() == MonitorElement::DQM_KIND_TH2D)) {
187  TH2 *th = me->getTH2F();
188  size_t nbinx = me->getNbinsX();
189  size_t nbiny = me->getNbinsY();
190 
191  float ymin = th->GetYaxis()->GetXmin();
192  float ymax = th->GetYaxis()->GetXmax();
193  std::string xtit = th->GetXaxis()->GetTitle();
194  std::string ytit = th->GetYaxis()->GetTitle();
195  float *xbins = new float[nbinx + 1];
196  for (size_t ix = 1; ix < nbinx + 1; ++ix) {
197  xbins[ix - 1] = th->GetXaxis()->GetBinLowEdge(ix);
198  if (ix == nbinx)
199  xbins[ix] = th->GetXaxis()->GetBinUpEdge(ix);
200  }
201 
202  std::string tit_new;
203  ibooker.setCurrentFolder(folder);
204 
205  if (folder == "ParticleFlow/PFElectronValidation/CompWithGenElectron") {
206  if (name == "delta_et_Over_et_VS_et_")
207  projection = ibooker.book1D("delta_et_Over_et", "E_{T} resolution;#DeltaE_{T}/E_{T}", nbiny, ymin, ymax);
208  if (name == "delta_et_VS_et_")
209  projection = ibooker.book1D("delta_et_", "#DeltaE_{T};#DeltaE_{T}", nbiny, ymin, ymax);
210  if (name == "delta_eta_VS_et_")
211  projection = ibooker.book1D("delta_eta_", "#Delta#eta;#Delta#eta", nbiny, ymin, ymax);
212  if (name == "delta_phi_VS_et_")
213  projection = ibooker.book1D("delta_phi_", "#Delta#phi;#Delta#phi", nbiny, ymin, ymax);
214  }
215 
216  if (projection) {
217  for (size_t iy = 1; iy < nbiny + 1; ++iy) {
218  projection->setBinContent(iy, th->ProjectionY("e")->GetBinContent(iy));
219  }
220  projection->setEntries(me->getEntries());
221  }
222 
223  delete[] xbins;
224  }
225 }
226 
227 //
228 // -- Create Profile Plots
229 //
231  DQMStore::IGetter &igetter,
233  std::string &name) {
234  MonitorElement *me = igetter.get(folder + "/" + name);
235  if (!me)
236  return;
237 
239  (me->kind() == MonitorElement::DQM_KIND_TH2D)) {
240  TH2 *th = me->getTH2F();
241  size_t nbinx = me->getNbinsX();
242 
243  float ymin = th->GetYaxis()->GetXmin();
244  float ymax = th->GetYaxis()->GetXmax();
245  std::string xtit = th->GetXaxis()->GetTitle();
246  std::string ytit = th->GetYaxis()->GetTitle();
247  double *xbins = new double[nbinx + 1];
248  for (size_t ix = 1; ix < nbinx + 1; ++ix) {
249  xbins[ix - 1] = th->GetXaxis()->GetBinLowEdge(ix);
250  if (ix == nbinx)
251  xbins[ix] = th->GetXaxis()->GetBinUpEdge(ix);
252  }
253 
254  std::string tit_new;
255  ibooker.setCurrentFolder(folder);
256  // TProfiles
257  MonitorElement *me_profile[2];
258  me_profile[0] = ibooker.bookProfile("profile_" + name, tit_new, nbinx, xbins, ymin, ymax, "");
259  me_profile[1] = ibooker.bookProfile("profileRMS_" + name, tit_new, nbinx, xbins, ymin, ymax, "s");
260  TProfile *profileX = th->ProfileX();
261  // size_t nbiny = me->getNbinsY();
262  // TProfile* profileX = th->ProfileX("",0,nbiny+1); add underflow and
263  // overflow
264  static const Int_t NUM_STAT = 7;
265  Double_t stats[NUM_STAT] = {0};
266  th->GetStats(stats);
267 
268  for (Int_t i = 0; i < 2; i++) {
269  if (me_profile[i]) {
270  for (size_t ix = 0; ix <= nbinx + 1; ++ix) {
271  me_profile[i]->setBinContent(ix, profileX->GetBinContent(ix) * profileX->GetBinEntries(ix));
272  // me_profile[i]->Fill( profileX->GetBinCenter(ix),
273  // profileX->GetBinContent(ix)*profileX->GetBinEntries(ix) ) ;
274  me_profile[i]->setBinEntries(ix, profileX->GetBinEntries(ix));
275  me_profile[i]->getTProfile()->GetSumw2()->fArray[ix] = profileX->GetSumw2()->fArray[ix];
276  // me_profile[i]->getTProfile()->GetBinSumw2()->fArray[ix] =
277  // profileX->GetBinSumw2()->fArray[ix]; // segmentation violation
278  }
279  }
280  me_profile[i]->getTProfile()->PutStats(stats);
281  me_profile[i]->setEntries(profileX->GetEntries());
282  }
283 
284  delete[] xbins;
285  }
286 }
287 
288 //
289 // -- Get Histogram Parameters
290 //
292  MonitorElement *me_slice, double &average, double &rms, double &mean, double &sigma) {
293  average = 0.0;
294  rms = 0.0;
295  mean = 0.0;
296  sigma = 0.0;
297 
298  if (!me_slice)
299  return;
300  if (me_slice->kind() == MonitorElement::DQM_KIND_TH1F) {
301  average = me_slice->getMean();
302  rms = me_slice->getRMS();
303  TH1F *th_slice = me_slice->getTH1F();
304  if (th_slice && th_slice->GetEntries() > 0) {
305  // need our own copy for thread safety
306  TF1 gaus("mygaus", "gaus");
307  th_slice->Fit(&gaus, "Q0 SERIAL");
308  sigma = gaus.GetParameter(2);
309  mean = gaus.GetParameter(1);
310  }
311  }
312 }
313 
314 //
315 // -- Create Efficiency Plots
316 //
318  DQMStore::IGetter &igetter,
320  std::string &name) {
321  MonitorElement *me1 = igetter.get(folder + "/" + name + "ref_");
322  MonitorElement *me2 = igetter.get(folder + "/" + name + "gen_");
323  if (!me1 || !me2)
324  return;
325 
326  TH1F *me1_forEff = (TH1F *)me1->getTH1F()->Rebin(2, "me1_forEff");
327  TH1F *me2_forEff = (TH1F *)me2->getTH1F()->Rebin(2, "me2_forEff");
328 
329  MonitorElement *me_eff;
330  if ((me1->kind() == MonitorElement::DQM_KIND_TH1F) && (me1->kind() == MonitorElement::DQM_KIND_TH1F)) {
331  // TH1* th1 = me1->getTH1F();
332  // size_t nbinx = me1->getNbinsX();
333  size_t nbinx = me1_forEff->GetNbinsX();
334 
335  // float xmin = th1->GetXaxis()->GetXmin();
336  // float xmax = th1->GetXaxis()->GetXmax();
337  float xmin = me1_forEff->GetXaxis()->GetXmin();
338  float xmax = me1_forEff->GetXaxis()->GetXmax();
339  std::string xtit = me1->getAxisTitle(1);
340  std::string tit_new;
341  tit_new = ";" + xtit + ";" + xtit + " efficiency";
342 
343  ibooker.setCurrentFolder(folder);
344  me_eff = ibooker.book1D("efficiency_" + name, tit_new, nbinx, xmin, xmax);
345 
346  me_eff->Reset();
347  me_eff->setEfficiencyFlag();
348  /*
349  double efficiency;
350  for (size_t ix = 1; ix < nbinx+1; ++ix) {
351  float val1 = me1->getBinContent(ix);
352  float val2 = me2->getBinContent(ix);
353  if (val2 > 0.0) efficiency = val1/val2;
354  else efficiency = 0;
355  me_eff->setBinContent(ix,efficiency);
356  }
357  */
358  // Binomial errors "B" asked by Florian
359  /*me1_forEff->Sumw2(); me2_forEff->Sumw2();*/ me_eff->getTH1F()->Sumw2();
360  me_eff->getTH1F()->Divide(me1_forEff, me2_forEff, 1, 1, "B");
361  }
362 }
363 
TProfile * getTProfile() const
T getParameter(std::string const &) const
void setBinContent(int binx, double content)
set content of bin (1-D)
const double xbins[]
void createResolutionPlots(DQMStore::IBooker &, DQMStore::IGetter &, std::string &folder, std::string &name)
Definition: PFClient.cc:106
MonitorElement * bookProfile(Args &&...args)
Definition: DQMStore.h:113
TH1F * getTH1F() const
std::string getAxisTitle(int axis=1) const
get x-, y- or z-axis title (axis=1, 2, 3 respectively)
void doProjection(DQMStore::IBooker &, DQMStore::IGetter &)
Definition: PFClient.cc:56
const std::string & getName() const
get name of ME
std::vector< std::string > projectionHistogramNames_
Definition: PFClient.h:34
std::vector< std::string > histogramNames_
Definition: PFClient.h:32
void dqmEndJob(DQMStore::IBooker &, DQMStore::IGetter &) override
Definition: PFClient.cc:28
double getMean(int axis=1) const
get mean value of histogram along x, y or z axis (axis=1, 2, 3 respectively)
void doSummaries(DQMStore::IBooker &, DQMStore::IGetter &)
Definition: PFClient.cc:40
void setEfficiencyFlag()
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:268
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:106
void setEntries(double nentries)
set # of entries
void Reset()
reset ME (ie. contents, errors, etc)
TH2F * getTH2F() const
MonitorElement * get(std::string const &path)
Definition: DQMStore.cc:303
void createEfficiencyPlots(DQMStore::IBooker &, DQMStore::IGetter &, std::string &folder, std::string &name)
Definition: PFClient.cc:317
std::vector< std::string > profileHistogramNames_
Definition: PFClient.h:35
void createProfilePlots(DQMStore::IBooker &, DQMStore::IGetter &, std::string &folder, std::string &name)
Definition: PFClient.cc:230
double getEntries() const
get # of entries
PFClient(const edm::ParameterSet &parameterSet)
Definition: PFClient.cc:15
std::vector< std::string > folderNames_
Definition: PFClient.h:31
int getNbinsY() const
get # of bins in Y-axis
void doEfficiency(DQMStore::IBooker &, DQMStore::IGetter &)
Definition: PFClient.cc:90
bool profileFlag_
Definition: PFClient.h:37
double getRMS(int axis=1) const
get RMS of histogram along x, y or z axis (axis=1, 2, 3 respectively)
std::vector< std::string > effHistogramNames_
Definition: PFClient.h:33
void doProfiles(DQMStore::IBooker &, DQMStore::IGetter &)
Definition: PFClient.cc:73
void createProjectionPlots(DQMStore::IBooker &, DQMStore::IGetter &, std::string &folder, std::string &name)
Definition: PFClient.cc:175
int getNbinsX() const
get # of bins in X-axis
void removeElement(Args &&...args)
Definition: DQMStore.h:199
void setBinEntries(int bin, double nentries)
set # of bin entries (to be used for profiles)
Kind kind() const
Get the type of the monitor element.
ParameterSet const & parameterSet(Provenance const &provenance)
Definition: Provenance.cc:11
bool efficiencyFlag_
Definition: PFClient.h:36
void getHistogramParameters(MonitorElement *me_slice, double &avarage, double &rms, double &mean, double &sigma)
Definition: PFClient.cc:291