CMS 3D CMS Logo

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