CMS 3D CMS Logo

HcalSimHitStudy.cc
Go to the documentation of this file.
4 
6 
8  g4Label = ps.getUntrackedParameter<std::string>("moduleLabel", "g4SimHits");
9  hcalHits = ps.getUntrackedParameter<std::string>("HitCollection", "HcalHits");
10  outFile_ = ps.getUntrackedParameter<std::string>("outputFile", "hcHit.root");
11  verbose_ = ps.getUntrackedParameter<bool>("Verbose", false);
12  testNumber_ = ps.getParameter<bool>("TestNumber");
13  hep17_ = ps.getParameter<bool>("hep17");
14  checkHit_ = true;
15 
16  tok_hits_ = consumes<edm::PCaloHitContainer>(edm::InputTag(g4Label, hcalHits));
17 
18  edm::LogInfo("HcalSim") << "Module Label: " << g4Label << " Hits: " << hcalHits << " / " << checkHit_
19  << " Output: " << outFile_;
20 }
21 
23 
26  es.get<HcalRecNumberingRecord>().get(pHRNDC);
27  hcons = &(*pHRNDC);
35 
36  // Get Phi segmentation from geometry, use the max phi number so that all iphi
37  // values are included.
38 
39  int NphiMax = hcons->getNPhi(0);
40 
41  NphiMax = (hcons->getNPhi(1) > NphiMax ? hcons->getNPhi(1) : NphiMax);
42  NphiMax = (hcons->getNPhi(2) > NphiMax ? hcons->getNPhi(2) : NphiMax);
43  NphiMax = (hcons->getNPhi(3) > NphiMax ? hcons->getNPhi(3) : NphiMax);
44 
45  // Center the iphi bins on the integers
46  iphi_min = 0.5;
47  iphi_max = NphiMax + 0.5;
49 
50  int iEtaHBMax = hcons->getEtaRange(0).second;
51  int iEtaHEMax = hcons->getEtaRange(1).second;
52  int iEtaHFMax = hcons->getEtaRange(2).second;
53  int iEtaHOMax = hcons->getEtaRange(3).second;
54 
55  // Retain classic behavior, all plots have same ieta range.
56  // Comment out code to allow each subdetector to have its on range
57 
58  int iEtaMax = (iEtaHBMax > iEtaHEMax ? iEtaHBMax : iEtaHEMax);
59  iEtaMax = (iEtaMax > iEtaHFMax ? iEtaMax : iEtaHFMax);
60  iEtaMax = (iEtaMax > iEtaHOMax ? iEtaMax : iEtaHOMax);
61 
62  iEtaHBMax = iEtaMax;
63  iEtaHEMax = iEtaMax;
64  iEtaHFMax = iEtaMax;
65  iEtaHOMax = iEtaMax;
66 
67  // Give an empty bin around the subdet ieta range to make it clear that all
68  // ieta rings have been included
69  ieta_min_HB = -iEtaHBMax - 1.5;
70  ieta_max_HB = iEtaHBMax + 1.5;
72 
73  ieta_min_HE = -iEtaHEMax - 1.5;
74  ieta_max_HE = iEtaHEMax + 1.5;
76 
77  ieta_min_HF = -iEtaHFMax - 1.5;
78  ieta_max_HF = iEtaHFMax + 1.5;
80 
81  ieta_min_HO = -iEtaHOMax - 1.5;
82  ieta_max_HO = iEtaHOMax + 1.5;
84 
85  Char_t hname[100];
86  Char_t htitle[100];
87 
88  ib.setCurrentFolder("HcalHitsV/HcalSimHitsTask");
89 
90  // Histograms for Hits
91  if (checkHit_) {
92  meAllNHit_ = ib.book1D("Hit01", "Number of Hits in HCal", 20000, 0., 20000.);
93  meBadDetHit_ = ib.book1D("Hit02", "Hits with wrong Det", 100, 0., 100.);
94  meBadSubHit_ = ib.book1D("Hit03", "Hits with wrong Subdet", 100, 0., 100.);
95  meBadIdHit_ = ib.book1D("Hit04", "Hits with wrong ID", 100, 0., 100.);
96  meHBNHit_ = ib.book1D("Hit05", "Number of Hits in HB", 20000, 0., 20000.);
97  meHENHit_ = ib.book1D("Hit06", "Number of Hits in HE", 10000, 0., 10000.);
98  meHONHit_ = ib.book1D("Hit07", "Number of Hits in HO", 10000, 0., 10000.);
99  meHFNHit_ = ib.book1D("Hit08", "Number of Hits in HF", 10000, 0., 10000.);
100  meDetectHit_ = ib.book1D("Hit09", "Detector ID", 50, 0., 50.);
101  meSubdetHit_ = ib.book1D("Hit10", "Subdetectors in HCal", 50, 0., 50.);
102  meDepthHit_ = ib.book1D("Hit11", "Depths in HCal", 20, 0., 20.);
103  meEtaHit_ = ib.book1D("Hit12", "Eta in HCal", ieta_bins_HF, ieta_min_HF, ieta_max_HF);
104  meEtaPhiHit_ =
105  ib.book2D("Hit12b", "Eta-phi in HCal", ieta_bins_HF, ieta_min_HF, ieta_max_HF, iphi_bins, iphi_min, iphi_max);
106  for (int depth = 1; depth <= maxDepth_; depth++) {
107  sprintf(hname, "Hit12bd%d", depth);
108  sprintf(htitle, "Eta-phi in HCal d%d", depth);
109  meEtaPhiHitDepth_.push_back(
111  }
112  // KC: There are different phi segmentation schemes, this plot uses wider
113  // bins to represent the most sparse segmentation
114  mePhiHit_ = ib.book1D("Hit13", "Phi in HCal (HB,HO)", iphi_bins, iphi_min, iphi_max);
115  mePhiHitb_ = ib.book1D("Hit13b", "Phi in HCal (HE,HF)", iphi_bins, iphi_min, iphi_max);
116  meEnergyHit_ = ib.book1D("Hit14", "Energy in HCal", 2000, 0., 20.);
117  meTimeHit_ = ib.book1D("Hit15", "Time in HCal", 528, 0., 528.);
118  meTimeWHit_ = ib.book1D("Hit16", "Time in HCal (E wtd)", 528, 0., 528.);
119  meHBDepHit_ = ib.book1D("Hit17", "Depths in HB", 20, 0., 20.);
120  meHEDepHit_ = ib.book1D("Hit18", "Depths in HE", 20, 0., 20.);
121  meHODepHit_ = ib.book1D("Hit19", "Depths in HO", 20, 0., 20.);
122  meHFDepHit_ = ib.book1D("Hit20", "Depths in HF", 20, 0., 20.);
123  meHFDepHitw_ = ib.book1D("Hit20b", "Depths in HF (p.e. weighted)", 20, 0., 20.);
124  meHBEtaHit_ = ib.book1D("Hit21", "Eta in HB", ieta_bins_HB, ieta_min_HB, ieta_max_HB);
125  meHEEtaHit_ = ib.book1D("Hit22", "Eta in HE", ieta_bins_HE, ieta_min_HE, ieta_max_HE);
126  meHOEtaHit_ = ib.book1D("Hit23", "Eta in HO", ieta_bins_HO, ieta_min_HO, ieta_max_HO);
127  meHFEtaHit_ = ib.book1D("Hit24", "Eta in HF", ieta_bins_HF, ieta_min_HF, ieta_max_HF);
128  meHBPhiHit_ = ib.book1D("Hit25", "Phi in HB", iphi_bins, iphi_min, iphi_max);
129  meHEPhiHit_ = ib.book1D("Hit26", "Phi in HE", iphi_bins, iphi_min, iphi_max);
130  meHOPhiHit_ = ib.book1D("Hit27", "Phi in HO", iphi_bins, iphi_min, iphi_max);
131  meHFPhiHit_ = ib.book1D("Hit28", "Phi in HF", iphi_bins, iphi_min, iphi_max);
132  meHBEneHit_ = ib.book1D("Hit29", "Energy in HB", 2000, 0., 20.);
133  meHEEneHit_ = ib.book1D("Hit30", "Energy in HE", 500, 0., 5.);
134  meHEP17EneHit_ = ib.book1D("Hit30b", "Energy in HEP17", 500, 0., 5.);
135  meHOEneHit_ = ib.book1D("Hit31", "Energy in HO", 500, 0., 5.);
136  meHFEneHit_ = ib.book1D("Hit32", "Energy in HF", 1001, -0.5, 1000.5);
137 
138  // HxEneMap, HxEneSum, HxEneSum_vs_ieta plot the sum of the simhits energy
139  // within a single ieta-iphi tower.
140 
141  meHBEneMap_ =
142  ib.book2D("HBEneMap", "HBEneMap", ieta_bins_HB, ieta_min_HB, ieta_max_HB, iphi_bins, iphi_min, iphi_max);
143  meHEEneMap_ =
144  ib.book2D("HEEneMap", "HEEneMap", ieta_bins_HE, ieta_min_HE, ieta_max_HE, iphi_bins, iphi_min, iphi_max);
145  meHOEneMap_ =
146  ib.book2D("HOEneMap", "HOEneMap", ieta_bins_HO, ieta_min_HO, ieta_max_HO, iphi_bins, iphi_min, iphi_max);
147  meHFEneMap_ =
148  ib.book2D("HFEneMap", "HFEneMap", ieta_bins_HF, ieta_min_HF, ieta_max_HF, iphi_bins, iphi_min, iphi_max);
149 
150  meHBEneSum_ = ib.book1D("HBEneSum", "HBEneSum", 2000, 0., 20.);
151  meHEEneSum_ = ib.book1D("HEEneSum", "HEEneSum", 500, 0., 5.);
152  meHOEneSum_ = ib.book1D("HOEneSum", "HOEneSum", 500, 0., 5.);
153  meHFEneSum_ = ib.book1D("HFEneSum", "HFEneSum", 1001, -0.5, 1000.5);
154 
156  "HBEneSum_vs_ieta", "HBEneSum_vs_ieta", ieta_bins_HB, ieta_min_HB, ieta_max_HB, 2011, -10.5, 2000.5, " ");
158  "HEEneSum_vs_ieta", "HEEneSum_vs_ieta", ieta_bins_HE, ieta_min_HE, ieta_max_HE, 2011, -10.5, 2000.5, " ");
160  "HOEneSum_vs_ieta", "HOEneSum_vs_ieta", ieta_bins_HO, ieta_min_HO, ieta_max_HO, 2011, -10.5, 2000.5, " ");
162  "HFEneSum_vs_ieta", "HFEneSum_vs_ieta", ieta_bins_HF, ieta_min_HF, ieta_max_HF, 2011, -10.5, 2000.5, " ");
163 
164  meHBTimHit_ = ib.book1D("Hit33", "Time in HB", 528, 0., 528.);
165  meHETimHit_ = ib.book1D("Hit34", "Time in HE", 528, 0., 528.);
166  meHOTimHit_ = ib.book1D("Hit35", "Time in HO", 528, 0., 528.);
167  meHFTimHit_ = ib.book1D("Hit36", "Time in HF", 528, 0., 528.);
168  // These are the zoomed in energy ranges
169  meHBEneHit2_ = ib.book1D("Hit37", "Energy in HB 2", 100, 0., 0.0001);
170  meHEEneHit2_ = ib.book1D("Hit38", "Energy in HE 2", 100, 0., 0.0001);
171  meHEP17EneHit2_ = ib.book1D("Hit38b", "Energy in HEP17 2", 100, 0., 0.0001);
172  meHOEneHit2_ = ib.book1D("Hit39", "Energy in HO 2", 100, 0., 0.0001);
173  meHFEneHit2_ = ib.book1D("Hit40", "Energy in HF 2", 100, 0.5, 100.5);
174  meHBL10Ene_ = ib.book1D("Hit41", "Log10Energy in HB", 140, -10., 4.);
175  meHEL10Ene_ = ib.book1D("Hit42", "Log10Energy in HE", 140, -10., 4.);
176  meHFL10Ene_ = ib.book1D("Hit43", "Log10Energy in HF", 50, -1., 4.);
177  meHOL10Ene_ = ib.book1D("Hit44", "Log10Energy in HO", 140, -10., 4.);
178  meHBL10EneP_ = ib.bookProfile("Hit45", "Log10Energy in HB vs Hit contribution", 140, -10., 4., 100, 0., 1.);
179  meHEL10EneP_ = ib.bookProfile("Hit46", "Log10Energy in HE vs Hit contribution", 140, -10., 4., 100, 0., 1.);
180  meHFL10EneP_ = ib.bookProfile("Hit47", "Log10Energy in HF vs Hit contribution", 140, -10., 4., 100, 0., 1.);
181  meHOL10EneP_ = ib.bookProfile("Hit48", "Log10Energy in HO vs Hit contribution", 140, -10., 4., 100, 0., 1.);
182  }
183 }
184 
185 /*void HcalSimHitStudy::endJob() {
186  if (dbe_ && outFile_.size() > 0) dbe_->save(outFile_);
187 }*/
188 
190  LogDebug("HcalSim") << "Run = " << e.id().run() << " Event = " << e.id().event();
191 
192  std::vector<PCaloHit> caloHits;
194 
195  bool getHits = false;
196  if (checkHit_) {
197  e.getByToken(tok_hits_, hitsHcal);
198  if (hitsHcal.isValid())
199  getHits = true;
200  }
201 
202  LogDebug("HcalSim") << "HcalValidation: Input flags Hits " << getHits;
203 
204  if (getHits) {
205  caloHits.insert(caloHits.end(), hitsHcal->begin(), hitsHcal->end());
206  LogDebug("HcalSim") << "HcalValidation: Hit buffer " << caloHits.size();
207  analyzeHits(caloHits);
208  }
209 }
210 
211 void HcalSimHitStudy::analyzeHits(std::vector<PCaloHit> &hits) {
212  int nHit = hits.size();
213  int nHB = 0, nHE = 0, nHO = 0, nHF = 0, nBad1 = 0, nBad2 = 0, nBad = 0;
214  std::vector<double> encontHB(140, 0.);
215  std::vector<double> encontHE(140, 0.);
216  std::vector<double> encontHF(140, 0.);
217  std::vector<double> encontHO(140, 0.);
218  double entotHB = 0, entotHE = 0, entotHF = 0, entotHO = 0;
219 
220  double HBEneMap[ieta_bins_HB][iphi_bins];
221  double HEEneMap[ieta_bins_HE][iphi_bins];
222  double HOEneMap[ieta_bins_HO][iphi_bins];
223  double HFEneMap[ieta_bins_HF][iphi_bins];
224 
225  // Works in ieta_min_Hx is < 0
226  int eta_offset_HB = -(int)ieta_min_HB;
227  int eta_offset_HE = -(int)ieta_min_HE;
228  int eta_offset_HO = -(int)ieta_min_HO;
229  int eta_offset_HF = -(int)ieta_min_HF;
230 
231  for (int i = 0; i < ieta_bins_HB; i++) {
232  for (int j = 0; j < iphi_bins; j++) {
233  HBEneMap[i][j] = 0.;
234  }
235  }
236 
237  for (int i = 0; i < ieta_bins_HE; i++) {
238  for (int j = 0; j < iphi_bins; j++) {
239  HEEneMap[i][j] = 0.;
240  }
241  }
242 
243  for (int i = 0; i < ieta_bins_HO; i++) {
244  for (int j = 0; j < iphi_bins; j++) {
245  HOEneMap[i][j] = 0.;
246  }
247  }
248 
249  for (int i = 0; i < ieta_bins_HF; i++) {
250  for (int j = 0; j < iphi_bins; j++) {
251  HFEneMap[i][j] = 0.;
252  }
253  }
254 
255  for (int i = 0; i < nHit; i++) {
256  double energy = hits[i].energy();
257  double log10en = log10(energy);
258  int log10i = int((log10en + 10.) * 10.);
259  double time = hits[i].time();
260  unsigned int id_ = hits[i].id();
261  int det, subdet, depth, eta, phi;
262  HcalDetId hid;
263  if (testNumber_)
264  hid = HcalHitRelabeller::relabel(id_, hcons);
265  else
266  hid = HcalDetId(id_);
267  det = hid.det();
268  subdet = hid.subdet();
269  depth = hid.depth();
270  eta = hid.ieta();
271  phi = hid.iphi();
272 
273  LogDebug("HcalSim") << "Hit[" << i << "] ID " << std::hex << id_ << std::dec << " Det " << det << " Sub " << subdet
274  << " depth " << depth << " Eta " << eta << " Phi " << phi << " E " << energy << " time "
275  << time;
276  if (det == 4) { // Check DetId.h
277  if (subdet == static_cast<int>(HcalBarrel))
278  nHB++;
279  else if (subdet == static_cast<int>(HcalEndcap))
280  nHE++;
281  else if (subdet == static_cast<int>(HcalOuter))
282  nHO++;
283  else if (subdet == static_cast<int>(HcalForward))
284  nHF++;
285  else {
286  nBad++;
287  nBad2++;
288  }
289  } else {
290  nBad++;
291  nBad1++;
292  }
293 
294  meDetectHit_->Fill(double(det));
295  if (det == 4) {
296  meSubdetHit_->Fill(double(subdet));
297  meDepthHit_->Fill(double(depth));
298  meEtaHit_->Fill(double(eta));
299  meEtaPhiHit_->Fill(double(eta), double(phi));
300  meEtaPhiHitDepth_[depth - 1]->Fill(double(eta), double(phi));
301 
302  // We will group the phi plots by HB,HO and HE,HF since these groups share
303  // similar segmentation schemes
304  if (subdet == static_cast<int>(HcalBarrel))
305  mePhiHit_->Fill(double(phi));
306  else if (subdet == static_cast<int>(HcalEndcap))
307  mePhiHitb_->Fill(double(phi));
308  else if (subdet == static_cast<int>(HcalOuter))
309  mePhiHit_->Fill(double(phi));
310  else if (subdet == static_cast<int>(HcalForward))
311  mePhiHitb_->Fill(double(phi));
312 
313  // KC: HF energy is in photoelectrons rather than eV, so it will not be
314  // included in total HCal energy
315  if (subdet != static_cast<int>(HcalForward)) {
316  meEnergyHit_->Fill(energy);
317 
318  // Since the HF energy is a different scale it does not make sense to
319  // include it in the Energy Weighted Plot
320  meTimeWHit_->Fill(double(time), energy);
321  }
322  meTimeHit_->Fill(time);
323 
324  if (subdet == static_cast<int>(HcalBarrel)) {
325  meHBDepHit_->Fill(double(depth));
326  meHBEtaHit_->Fill(double(eta));
327  meHBPhiHit_->Fill(double(phi));
328  meHBEneHit_->Fill(energy);
329  meHBEneHit2_->Fill(energy);
330  meHBTimHit_->Fill(time);
331  meHBL10Ene_->Fill(log10en);
332  if (log10i >= 0 && log10i < 140)
333  encontHB[log10i] += energy;
334  entotHB += energy;
335 
336  HBEneMap[eta + eta_offset_HB][phi - 1] += energy;
337 
338  } else if (subdet == static_cast<int>(HcalEndcap)) {
339  meHEDepHit_->Fill(double(depth));
340  meHEEtaHit_->Fill(double(eta));
341  meHEPhiHit_->Fill(double(phi));
342 
343  bool isHEP17 = (phi >= 63) && (phi <= 66) && (eta > 0);
344  if (hep17_) {
345  if (!isHEP17) {
346  meHEEneHit_->Fill(energy);
347  meHEEneHit2_->Fill(energy);
348  } else {
349  meHEP17EneHit_->Fill(energy);
350  meHEP17EneHit2_->Fill(energy);
351  }
352  } else {
353  meHEEneHit_->Fill(energy);
354  meHEEneHit2_->Fill(energy);
355  }
356 
357  meHETimHit_->Fill(time);
358  meHEL10Ene_->Fill(log10en);
359  if (log10i >= 0 && log10i < 140)
360  encontHE[log10i] += energy;
361  entotHE += energy;
362 
363  HEEneMap[eta + eta_offset_HE][phi - 1] += energy;
364 
365  } else if (subdet == static_cast<int>(HcalOuter)) {
366  meHODepHit_->Fill(double(depth));
367  meHOEtaHit_->Fill(double(eta));
368  meHOPhiHit_->Fill(double(phi));
369  meHOEneHit_->Fill(energy);
370  meHOEneHit2_->Fill(energy);
371  meHOTimHit_->Fill(time);
372  meHOL10Ene_->Fill(log10en);
373  if (log10i >= 0 && log10i < 140)
374  encontHO[log10i] += energy;
375  entotHO += energy;
376 
377  HOEneMap[eta + eta_offset_HO][phi - 1] += energy;
378 
379  } else if (subdet == static_cast<int>(HcalForward)) {
380  meHFDepHit_->Fill(double(depth));
381  meHFDepHitw_->Fill(double(depth), energy);
382  meHFEtaHit_->Fill(double(eta));
383  meHFPhiHit_->Fill(double(phi));
384  meHFEneHit_->Fill(energy);
385  meHFEneHit2_->Fill(energy);
386  meHFTimHit_->Fill(time);
387  meHFL10Ene_->Fill(log10en);
388  if (log10i >= 0 && log10i < 140)
389  encontHF[log10i] += energy;
390  entotHF += energy;
391 
392  HFEneMap[eta + eta_offset_HF][phi - 1] += energy;
393  }
394  }
395  }
396  if (entotHB != 0)
397  for (int i = 0; i < 140; i++)
398  meHBL10EneP_->Fill(-10. + (float(i) + 0.5) / 10., encontHB[i] / entotHB);
399  if (entotHE != 0)
400  for (int i = 0; i < 140; i++)
401  meHEL10EneP_->Fill(-10. + (float(i) + 0.5) / 10., encontHE[i] / entotHE);
402  if (entotHF != 0)
403  for (int i = 0; i < 140; i++)
404  meHFL10EneP_->Fill(-10. + (float(i) + 0.5) / 10., encontHF[i] / entotHF);
405  if (entotHO != 0)
406  for (int i = 0; i < 140; i++)
407  meHOL10EneP_->Fill(-10. + (float(i) + 0.5) / 10., encontHO[i] / entotHO);
408 
409  meAllNHit_->Fill(double(nHit));
410  meBadDetHit_->Fill(double(nBad1));
411  meBadSubHit_->Fill(double(nBad2));
412  meBadIdHit_->Fill(double(nBad));
413  meHBNHit_->Fill(double(nHB));
414  meHENHit_->Fill(double(nHE));
415  meHONHit_->Fill(double(nHO));
416  meHFNHit_->Fill(double(nHF));
417 
418  for (int i = 0; i < ieta_bins_HB; i++) {
419  for (int j = 0; j < iphi_bins; j++) {
420  if (HBEneMap[i][j] != 0) {
421  meHBEneSum_->Fill(HBEneMap[i][j]);
422  meHBEneSum_vs_ieta_->Fill((i - eta_offset_HB), HBEneMap[i][j]);
423  meHBEneMap_->Fill((i - eta_offset_HB), j + 1, HBEneMap[i][j]);
424  }
425  }
426  }
427 
428  for (int i = 0; i < ieta_bins_HE; i++) {
429  for (int j = 0; j < iphi_bins; j++) {
430  if (HEEneMap[i][j] != 0) {
431  meHEEneSum_->Fill(HEEneMap[i][j]);
432  meHEEneSum_vs_ieta_->Fill((i - eta_offset_HE), HEEneMap[i][j]);
433  meHEEneMap_->Fill((i - eta_offset_HE), j + 1, HEEneMap[i][j]);
434  }
435  }
436  }
437 
438  for (int i = 0; i < ieta_bins_HO; i++) {
439  for (int j = 0; j < iphi_bins; j++) {
440  if (HOEneMap[i][j] != 0) {
441  meHOEneSum_->Fill(HOEneMap[i][j]);
442  meHOEneSum_vs_ieta_->Fill((i - eta_offset_HO), HOEneMap[i][j]);
443  meHOEneMap_->Fill((i - eta_offset_HO), j + 1, HOEneMap[i][j]);
444  }
445  }
446  }
447 
448  for (int i = 0; i < ieta_bins_HF; i++) {
449  for (int j = 0; j < iphi_bins; j++) {
450  if (HFEneMap[i][j] != 0) {
451  meHFEneSum_->Fill(HFEneMap[i][j]);
452  meHFEneSum_vs_ieta_->Fill((i - eta_offset_HF), HFEneMap[i][j]);
453  meHFEneMap_->Fill((i - eta_offset_HF), j + 1, HFEneMap[i][j]);
454  }
455  }
456  }
457 
458  LogDebug("HcalSim") << "HcalSimHitStudy::analyzeHits: HB " << nHB << " HE " << nHE << " HO " << nHO << " HF " << nHF
459  << " Bad " << nBad << " All " << nHit;
460 }
#define LogDebug(id)
RunNumber_t run() const
Definition: EventID.h:39
MonitorElement * mePhiHit_
MonitorElement * meHBEneMap_
edm::EDGetTokenT< edm::PCaloHitContainer > tok_hits_
T getParameter(std::string const &) const
EventNumber_t event() const
Definition: EventID.h:41
T getUntrackedParameter(std::string const &, T const &) const
MonitorElement * meHEDepHit_
MonitorElement * meHOL10Ene_
MonitorElement * meHOEneSum_vs_ieta_
MonitorElement * meHOEneSum_
HcalSubdetector subdet() const
get the subdetector
Definition: HcalDetId.h:146
MonitorElement * meHOPhiHit_
MonitorElement * bookProfile(Args &&...args)
Definition: DQMStore.h:113
MonitorElement * meDetectHit_
MonitorElement * meHFNHit_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
void analyze(const edm::Event &e, const edm::EventSetup &c) override
MonitorElement * meHFEtaHit_
MonitorElement * meHFL10EneP_
MonitorElement * meHBEneSum_vs_ieta_
MonitorElement * meHFEneHit2_
MonitorElement * meHOL10EneP_
MonitorElement * meHFDepHitw_
MonitorElement * meHBNHit_
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
MonitorElement * meHFTimHit_
MonitorElement * meHEEneSum_vs_ieta_
std::pair< int, int > getEtaRange(const int &i) const
MonitorElement * meHBEneHit2_
MonitorElement * meHONHit_
void Fill(long long x)
MonitorElement * meAllNHit_
const HcalDDDRecConstants * hcons
MonitorElement * meHBL10Ene_
MonitorElement * meHEEneSum_
MonitorElement * meBadDetHit_
MonitorElement * meHEP17EneHit2_
int depth() const
get the tower depth
Definition: HcalDetId.h:166
MonitorElement * meHFL10Ene_
MonitorElement * mePhiHitb_
void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:268
MonitorElement * meHETimHit_
MonitorElement * meTimeHit_
void analyzeHits(std::vector< PCaloHit > &)
MonitorElement * meHFDepHit_
MonitorElement * meHBEneSum_
int ieta() const
get the cell ieta
Definition: HcalDetId.h:159
MonitorElement * meHOEneMap_
MonitorElement * meHBTimHit_
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:106
MonitorElement * meHFEneSum_vs_ieta_
MonitorElement * meTimeWHit_
bool isValid() const
Definition: HandleBase.h:74
MonitorElement * meHBDepHit_
MonitorElement * meHFEneSum_
MonitorElement * meHFEneMap_
MonitorElement * meBadSubHit_
int iphi() const
get the cell iphi
Definition: HcalDetId.h:161
MonitorElement * meHENHit_
MonitorElement * book2D(Args &&...args)
Definition: DQMStore.h:109
MonitorElement * meBadIdHit_
MonitorElement * meHEL10EneP_
MonitorElement * meEnergyHit_
std::string hcalHits
int getMaxDepth(const int &type) const
MonitorElement * meHOEtaHit_
MonitorElement * meEtaPhiHit_
MonitorElement * meHEEneMap_
MonitorElement * meEtaHit_
MonitorElement * meHBEneHit_
MonitorElement * meHOTimHit_
MonitorElement * meHBL10EneP_
edm::EventID id() const
Definition: EventBase.h:59
MonitorElement * meHEEneHit2_
MonitorElement * meHFPhiHit_
T get() const
Definition: EventSetup.h:71
std::string outFile_
MonitorElement * meDepthHit_
MonitorElement * meHBEtaHit_
MonitorElement * meHODepHit_
std::vector< MonitorElement * > meEtaPhiHitDepth_
MonitorElement * meHEL10Ene_
DetId relabel(const uint32_t testId) const
MonitorElement * meHEEneHit_
std::string g4Label
MonitorElement * meHEPhiHit_
HcalSimHitStudy(const edm::ParameterSet &ps)
MonitorElement * meSubdetHit_
MonitorElement * meHBPhiHit_
MonitorElement * meHEP17EneHit_
Definition: Run.h:45
MonitorElement * meHFEneHit_
~HcalSimHitStudy() override
int getNPhi(const int &type) const
MonitorElement * meHEEtaHit_
ib
Definition: cuy.py:662
MonitorElement * meHOEneHit_
constexpr Detector det() const
get the detector field from this detid
Definition: DetId.h:39
MonitorElement * meHOEneHit2_