CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
CaloSimHitAnalysis.cc
Go to the documentation of this file.
3 
12 
15 
21 
26 
31 
32 #include <TH1F.h>
33 #include <TH2F.h>
34 
35 #include <memory>
36 #include <iostream>
37 #include <fstream>
38 #include <vector>
39 #include <map>
40 #include <string>
41 
42 //#define EDM_ML_DEBUG
43 
44 class CaloSimHitAnalysis : public edm::one::EDAnalyzer<edm::one::WatchRuns, edm::one::SharedResources> {
45 public:
47  ~CaloSimHitAnalysis() override {}
48  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
49 
50 protected:
51  void beginJob() override {}
52  void analyze(edm::Event const&, edm::EventSetup const&) override;
53  void beginRun(edm::Run const&, edm::EventSetup const&) override {}
54  void endRun(edm::Run const&, edm::EventSetup const&) override {}
55 
56  void analyzeHits(std::vector<PCaloHit>&, int);
57  void analyzePassiveHits(std::vector<PassiveHit>& hits);
58 
59 private:
61  const std::vector<std::string> hitLab_;
62  const std::vector<double> timeSliceUnit_;
64  const bool testNumber_, passive_;
65  const int allSteps_;
66  const std::vector<std::string> detNames_;
68  std::vector<edm::EDGetTokenT<edm::PCaloHitContainer> > toks_calo_;
70 
73 
74  static constexpr int nCalo_ = 6;
79  TH2F *h_rz_, *h_rz1_, *h_etaphi_;
81  std::vector<TH1F*> h_edepTk_, h_timeTk_;
82  std::vector<TH2F*> h_rzH_;
83  std::map<int, unsigned int> etaDepth_;
84 };
85 
87  : g4Label_(ps.getUntrackedParameter<std::string>("moduleLabel", "g4SimHits")),
88  hitLab_(ps.getParameter<std::vector<std::string> >("hitCollection")),
89  timeSliceUnit_(ps.getParameter<std::vector<double> >("timeSliceUnit")),
90  maxEnergy_(ps.getUntrackedParameter<double>("maxEnergy", 250.0)),
91  maxTime_(ps.getUntrackedParameter<double>("maxTime", 1000.0)),
92  tMax_(ps.getUntrackedParameter<double>("timeCut", 100.0)),
93  tScale_(ps.getUntrackedParameter<double>("timeScale", 1.0)),
94  tCut_(ps.getUntrackedParameter<double>("timeThreshold", 15.0)),
95  testNumber_(ps.getUntrackedParameter<bool>("testNumbering", false)),
96  passive_(ps.getUntrackedParameter<bool>("passiveHits", false)),
97  allSteps_(ps.getUntrackedParameter<int>("allSteps", 100)),
98  detNames_(ps.getUntrackedParameter<std::vector<std::string> >("detNames")),
100  usesResource(TFileService::kSharedResource);
101 
102  // register for data access
103  for (unsigned int i = 0; i < hitLab_.size(); i++)
104  toks_calo_.emplace_back(consumes<edm::PCaloHitContainer>(edm::InputTag(g4Label_, hitLab_[i])));
105  tok_passive_ = consumes<edm::PassiveHitContainer>(edm::InputTag(g4Label_, "AllPassiveHits"));
106 
107  edm::LogVerbatim("HitStudy") << "Module Label: " << g4Label_ << " Hits|timeSliceUnit:";
108  for (unsigned int i = 0; i < hitLab_.size(); i++)
109  edm::LogVerbatim("HitStudy") << "[" << i << "] " << hitLab_[i] << " " << timeSliceUnit_[i];
110  edm::LogVerbatim("HitStudy") << "Passive Hits " << passive_ << " from AllPassiveHits";
111  edm::LogVerbatim("HitStudy") << "maxEnergy: " << maxEnergy_ << " maxTime: " << maxTime_ << " tMax: " << tMax_;
112  for (unsigned int k = 0; k < detNames_.size(); ++k)
113  edm::LogVerbatim("HitStudy") << "Detector[" << k << "] " << detNames_[k];
114 
116  if (!tfile.isAvailable())
117  throw cms::Exception("BadConfig") << "TFileService unavailable: "
118  << "please add it to config file";
119  char name[29], title[120];
120  std::string dets[nCalo_] = {"EB", "EE", "HB", "HE", "HO", "HF"};
121  for (int i = 0; i < nCalo_; i++) {
122  sprintf(name, "Hit%d", i);
123  sprintf(title, "Number of hits in %s", dets[i].c_str());
124  h_hit_[i] = tfile->make<TH1F>(name, title, 100, 0., 20000.);
125  h_hit_[i]->GetXaxis()->SetTitle(title);
126  h_hit_[i]->GetYaxis()->SetTitle("Events");
127  sprintf(name, "Time%d", i);
128  sprintf(title, "Time of the hit (ns) in %s", dets[i].c_str());
129  h_time_[i] = tfile->make<TH1F>(name, title, 100, 0., 200.);
130  h_time_[i]->GetXaxis()->SetTitle(title);
131  h_time_[i]->GetYaxis()->SetTitle("Hits");
132  sprintf(name, "TimeT%d", i);
133  sprintf(title, "Time of each hit (ns) in %s", dets[i].c_str());
134  h_timeT_[i] = tfile->make<TH1F>(name, title, 100, 0., 200.);
135  h_timeT_[i]->GetXaxis()->SetTitle(title);
136  h_timeT_[i]->GetYaxis()->SetTitle("Hits");
137  double ymax = (i > 1) ? 0.01 : 0.1;
138  double ymx0 = (i > 1) ? 0.0025 : 0.025;
139  sprintf(name, "Edep%d", i);
140  sprintf(title, "Energy deposit (GeV) in %s", dets[i].c_str());
141  h_edep_[i] = tfile->make<TH1F>(name, title, 100, 0., ymax);
142  h_edep_[i]->GetXaxis()->SetTitle(title);
143  h_edep_[i]->GetYaxis()->SetTitle("Hits");
144  sprintf(name, "EdepT%d", i);
145  sprintf(title, "Energy deposit (GeV) of each hit in %s", dets[i].c_str());
146  h_edepT_[i] = tfile->make<TH1F>(name, title, 100, 0., ymx0);
147  h_edepT_[i]->GetXaxis()->SetTitle(title);
148  h_edepT_[i]->GetYaxis()->SetTitle("Hits");
149  sprintf(name, "EdepEM%d", i);
150  sprintf(title, "Energy deposit (GeV) by EM particles in %s", dets[i].c_str());
151  h_edepEM_[i] = tfile->make<TH1F>(name, title, 100, 0., ymx0);
152  h_edepEM_[i]->GetXaxis()->SetTitle(title);
153  h_edepEM_[i]->GetYaxis()->SetTitle("Hits");
154  sprintf(name, "EdepHad%d", i);
155  sprintf(title, "Energy deposit (GeV) by hadrons in %s", dets[i].c_str());
156  h_edepHad_[i] = tfile->make<TH1F>(name, title, 100, 0., ymx0);
157  h_edepHad_[i]->GetXaxis()->SetTitle(title);
158  h_edepHad_[i]->GetYaxis()->SetTitle("Hits");
159  sprintf(name, "Edep15%d", i);
160  sprintf(title, "Energy deposit (GeV) for T > %4.0f ns in %s", tCut_, dets[i].c_str());
161  h_edep1_[i] = tfile->make<TH1F>(name, title, 100, 0., ymax);
162  h_edep1_[i]->GetXaxis()->SetTitle(title);
163  h_edep1_[i]->GetYaxis()->SetTitle("Hits");
164  sprintf(name, "EdepT15%d", i);
165  sprintf(title, "Energy deposit (GeV) of each hit for T > %4.0f in %s", tCut_, dets[i].c_str());
166  h_edepT1_[i] = tfile->make<TH1F>(name, title, 100, 0., ymx0);
167  h_edepT1_[i]->GetXaxis()->SetTitle(title);
168  h_edepT1_[i]->GetYaxis()->SetTitle("Hits");
169  ymax = (i > 1) ? 1.0 : maxEnergy_;
170  sprintf(name, "Etot%d", i);
171  sprintf(title, "Total energy deposit (GeV) in %s", dets[i].c_str());
172  h_etot_[i] = tfile->make<TH1F>(name, title, 50, 0., ymax);
173  h_etot_[i]->GetXaxis()->SetTitle(title);
174  h_etot_[i]->GetYaxis()->SetTitle("Events");
175  sprintf(name, "EtotG%d", i);
176  sprintf(title, "Total energy deposit (GeV) in %s (t < 100 ns)", dets[i].c_str());
177  h_etotg_[i] = tfile->make<TH1F>(name, title, 50, 0., ymax);
178  h_etotg_[i]->GetXaxis()->SetTitle(title);
179  h_etotg_[i]->GetYaxis()->SetTitle("Events");
180  sprintf(name, "rr%d", i);
181  sprintf(title, "R of hit point (cm) in %s", dets[i].c_str());
182  h_rr_[i] = tfile->make<TH1F>(name, title, 100, 0., 250.);
183  h_rr_[i]->GetXaxis()->SetTitle(title);
184  h_rr_[i]->GetYaxis()->SetTitle("Hits");
185  sprintf(name, "zz%d", i);
186  sprintf(title, "z of hit point (cm) in %s", dets[i].c_str());
187  h_zz_[i] = tfile->make<TH1F>(name, title, 240, -600., 600.);
188  h_zz_[i]->GetXaxis()->SetTitle(title);
189  h_zz_[i]->GetYaxis()->SetTitle("Hits");
190  sprintf(name, "eta%d", i);
191  sprintf(title, "#eta of hit point in %s", dets[i].c_str());
192  h_eta_[i] = tfile->make<TH1F>(name, title, 100, -5.0, 5.0);
193  h_eta_[i]->GetXaxis()->SetTitle(title);
194  h_eta_[i]->GetYaxis()->SetTitle("Hits");
195  sprintf(name, "phi%d", i);
196  sprintf(title, "#phi of hit point in %s", dets[i].c_str());
197  h_phi_[i] = tfile->make<TH1F>(name, title, 100, -M_PI, M_PI);
198  h_phi_[i]->GetXaxis()->SetTitle(title);
199  h_phi_[i]->GetYaxis()->SetTitle("Hits");
200  }
201  sprintf(title, "R vs Z of hit point");
202  h_rz_ = tfile->make<TH2F>("rz", title, 120, 0., 600., 100, 0., 250.);
203  h_rz_->GetXaxis()->SetTitle("z (cm)");
204  h_rz_->GetYaxis()->SetTitle("R (cm)");
205  sprintf(title, "R vs Z of hit point for hits with T > %4.0f ns", tCut_);
206  h_rz1_ = tfile->make<TH2F>("rz2", title, 120, 0., 600., 100, 0., 250.);
207  h_rz1_->GetXaxis()->SetTitle("z (cm)");
208  h_rz1_->GetYaxis()->SetTitle("R (cm)");
209  sprintf(title, "#phi vs #eta of hit point");
210  h_etaphi_ = tfile->make<TH2F>("etaphi", title, 100, 0., 5., 100, 0., M_PI);
211  h_etaphi_->GetXaxis()->SetTitle("#eta");
212  h_etaphi_->GetYaxis()->SetTitle("#phi");
213 
214  if (passive_) {
215  h_hitp_ = tfile->make<TH1F>("hitp", "All Steps", 100, 0.0, 20000.0);
216  h_hitp_->GetXaxis()->SetTitle("Hits");
217  h_hitp_->GetYaxis()->SetTitle("Events");
218  h_trackp_ = tfile->make<TH1F>("trackp", "All Steps", 100, 0.0, 200000.0);
219  h_hitp_->GetXaxis()->SetTitle("Tracks");
220  h_hitp_->GetYaxis()->SetTitle("Events");
221  h_edepp_ = tfile->make<TH1F>("edepp", "All Steps", 100, 0.0, 50.0);
222  h_edepp_->GetXaxis()->SetTitle("Energy Deposit (MeV)");
223  h_edepp_->GetYaxis()->SetTitle("Hits");
224  h_timep_ = tfile->make<TH1F>("timep", "All Steps", 100, 0.0, 100.0);
225  h_timep_->GetXaxis()->SetTitle("Hits");
226  h_timep_->GetYaxis()->SetTitle("Hit Time (ns)");
227  h_stepp_ = tfile->make<TH1F>("stepp", "All Steps", 1000, 0.0, 100.0);
228  h_stepp_->GetXaxis()->SetTitle("Hits");
229  h_stepp_->GetYaxis()->SetTitle("Step length (cm)");
230  for (unsigned int k = 0; k < detNames_.size(); ++k) {
231  sprintf(name, "edept%d", k);
232  sprintf(title, "Energy Deposit (MeV) in %s", detNames_[k].c_str());
233  h_edepTk_.emplace_back(tfile->make<TH1F>(name, title, 100, 0.0, 1.0));
234  h_edepTk_.back()->GetYaxis()->SetTitle("Hits");
235  h_edepTk_.back()->GetXaxis()->SetTitle(title);
236  sprintf(name, "timet%d", k);
237  sprintf(title, "Hit Time (ns) in %s", detNames_[k].c_str());
238  h_timeTk_.emplace_back(tfile->make<TH1F>(name, title, 100, 0.0, 100.0));
239  h_timeTk_.back()->GetYaxis()->SetTitle("Hits");
240  h_timeTk_.back()->GetXaxis()->SetTitle(title);
241  }
242  if ((allSteps_ / 100) % 10 > 0) {
243  for (int eta = 1; eta <= 29; ++eta) {
244  int dmax = (eta < 16) ? 4 : 7;
245  for (int depth = 1; depth <= dmax; ++depth) {
246  sprintf(name, "Eta%dDepth%d", eta, depth);
247  sprintf(title, "R vs Z (#eta = %d, depth = %d)", eta, depth);
248  etaDepth_[eta * 100 + depth] = h_rzH_.size();
249  h_rzH_.emplace_back(tfile->make<TH2F>(name, title, 120, 0., 600., 100, 0., 250.));
250  h_rzH_.back()->GetXaxis()->SetTitle("z (cm)");
251  h_rzH_.back()->GetYaxis()->SetTitle("R (cm)");
252  }
253  }
254  }
255  }
256 }
257 
260  std::vector<std::string> labels = {"EcalHitsEB1", "EcalHitsEE1", "HcalHits1"};
261  std::vector<double> times = {1, 1, 1};
262  desc.addUntracked<std::string>("moduleLabel", "g4SimHits");
263  desc.add<std::vector<std::string> >("hitCollection", labels);
264  desc.add<std::vector<double> >("timeSliceUnit", times);
265  desc.addUntracked<double>("maxEnergy", 250.0);
266  desc.addUntracked<double>("maxTime", 1000.0);
267  desc.addUntracked<double>("timeCut", 100.0);
268  desc.addUntracked<double>("timeScale", 1.0);
269  desc.addUntracked<double>("timeThreshold", 15.0);
270  desc.addUntracked<bool>("testNumbering", false);
271  desc.addUntracked<bool>("passiveHits", false);
272  std::vector<std::string> names = {"PixelBarrel", "PixelForward", "TIB", "TID", "TOB", "TEC"};
273  desc.addUntracked<std::vector<std::string> >("detNames", names);
274  desc.addUntracked<int>("allStep", 100);
275  descriptions.add("caloSimHitAnalysis", desc);
276 }
277 
279  edm::LogVerbatim("HitStudy") << "CaloSimHitAnalysis:Run = " << e.id().run() << " Event = " << e.id().event();
280 
283 
284  for (unsigned int i = 0; i < toks_calo_.size(); i++) {
286  e.getByToken(toks_calo_[i], hitsCalo);
287  bool getHits = (hitsCalo.isValid());
288 #ifdef EDM_ML_DEBUG
289  edm::LogVerbatim("HitStudy") << "CaloSimHitAnalysis: Input flags Hits[" << i << "]: " << getHits;
290 #endif
291  if (getHits) {
292  std::vector<PCaloHit> caloHits;
293  caloHits.insert(caloHits.end(), hitsCalo->begin(), hitsCalo->end());
294 #ifdef EDM_ML_DEBUG
295  edm::LogVerbatim("HitStudy") << "CaloSimHitAnalysis: Hit buffer [" << i << "] " << caloHits.size();
296 #endif
297  analyzeHits(caloHits, i);
298  }
299  }
300 
301  if (passive_) {
303  e.getByToken(tok_passive_, hitsPassive);
304  bool getHits = (hitsPassive.isValid());
305 #ifdef EDM_ML_DEBUG
306  edm::LogVerbatim("HitStudy") << "CaloSimHitAnalysis: Passive: " << getHits;
307 #endif
308  if (getHits) {
309  std::vector<PassiveHit> passiveHits;
310  passiveHits.insert(passiveHits.end(), hitsPassive->begin(), hitsPassive->end());
311 #ifdef EDM_ML_DEBUG
312  edm::LogVerbatim("HitStudy") << "CaloSimHitAnalysis: Passive Hit buffer " << passiveHits.size();
313 #endif
314  analyzePassiveHits(passiveHits);
315  }
316  }
317 }
318 
319 void CaloSimHitAnalysis::analyzeHits(std::vector<PCaloHit>& hits, int indx) {
320  int nHit = hits.size();
321  int nHB = 0, nHE = 0, nHO = 0, nHF = 0, nEB = 0, nEE = 0, nBad = 0, iHit = 0;
322  std::map<CaloHitID, std::pair<double, double> > hitMap;
323  double etot[nCalo_], etotG[nCalo_];
324  for (unsigned int k = 0; k < nCalo_; ++k) {
325  etot[k] = etotG[k] = 0;
326  }
327  for (const auto& hit : hits) {
328  double edep = hit.energy();
329  double time = tScale_ * hit.time();
330  uint32_t id = hit.id();
331  int itra = hit.geantTrackId();
332  double edepEM = hit.energyEM();
333  double edepHad = hit.energyHad();
334  int idx(-1);
335  if (indx != 2) {
336  idx = indx;
337  if (indx == 0)
338  ++nEB;
339  else
340  ++nEE;
341  } else {
342  int subdet(0);
343  if (testNumber_) {
344  int ieta(0), phi(0), z(0), lay(0), depth(0);
345  HcalTestNumbering::unpackHcalIndex(id, subdet, z, depth, ieta, phi, lay);
346  id = HcalDetId(static_cast<HcalSubdetector>(subdet), z * ieta, phi, depth).rawId();
347  } else {
348  subdet = HcalDetId(id).subdet();
349  }
350  if (subdet == static_cast<int>(HcalBarrel)) {
351  idx = indx;
352  nHB++;
353  } else if (subdet == static_cast<int>(HcalEndcap)) {
354  idx = indx + 1;
355  nHE++;
356  } else if (subdet == static_cast<int>(HcalOuter)) {
357  idx = indx + 2;
358  nHO++;
359  } else if (subdet == static_cast<int>(HcalForward)) {
360  idx = indx + 3;
361  nHF++;
362  }
363  }
364 #ifdef EDM_ML_DEBUG
365  edm::LogVerbatim("HitStudy") << "Hit[" << iHit << ":" << nHit << ":" << idx << "] E " << edep << ":" << edepEM
366  << ":" << edepHad << " T " << time << " itra " << itra << " ID " << std::hex << id
367  << std::dec;
368 #endif
369  ++iHit;
370  if (idx >= 0) {
371  CaloHitID hid(id, time, itra, 0, timeSliceUnit_[indx]);
372  auto itr = hitMap.find(hid);
373  if (itr == hitMap.end())
374  hitMap[hid] = std::make_pair(time, edep);
375  else
376  ((itr->second).second) += edep;
377  h_edepT_[idx]->Fill(edep);
378  h_timeT_[idx]->Fill(time);
379  if (edepEM > 0)
380  h_edepEM_[idx]->Fill(edepEM);
381  if (edepHad > 0)
382  h_edepHad_[idx]->Fill(edepHad);
383  if (time > tCut_)
384  h_edepT1_[idx]->Fill(edep);
385  } else {
386  ++nBad;
387  }
388  }
389 
390  //Now make plots
391  for (auto itr = hitMap.begin(); itr != hitMap.end(); ++itr) {
392  int idx = -1;
394  DetId id((itr->first).unitID());
395 #ifdef EDM_ML_DEBUG
396  edm::LogVerbatim("HitStudy") << "Index " << indx << " Geom " << caloGeometry_ << ":" << hcalGeom_ << " "
397  << std::hex << id.rawId() << std::dec;
398 #endif
399  if (indx != 2) {
400  idx = indx;
401  point = caloGeometry_->getPosition(id);
402  } else {
403  int subdet = id.subdetId();
404  if (subdet == static_cast<int>(HcalBarrel)) {
405  idx = indx;
406  } else if (subdet == static_cast<int>(HcalEndcap)) {
407  idx = indx + 1;
408  } else if (subdet == static_cast<int>(HcalOuter)) {
409  idx = indx + 2;
410  } else if (subdet == static_cast<int>(HcalForward)) {
411  idx = indx + 3;
412  }
413  point = hcalGeom_->getPosition(id);
414  }
415  double edep = (itr->second).second;
416  double time = (itr->second).first;
417 #ifdef EDM_ML_DEBUG
418  edm::LogVerbatim("HitStudy") << "Index " << idx << ":" << nCalo_ << " Point " << point << " E " << edep << " T "
419  << time;
420 #endif
421  if (idx >= 0) {
422  h_time_[idx]->Fill(time);
423  h_edep_[idx]->Fill(edep);
424  h_rr_[idx]->Fill(point.perp());
425  h_zz_[idx]->Fill(point.z());
426  h_eta_[idx]->Fill(point.eta());
427  h_phi_[idx]->Fill(point.phi());
428  h_rz_->Fill(std::abs(point.z()), point.perp());
429  h_etaphi_->Fill(std::abs(point.eta()), std::abs(point.phi()));
430  etot[idx] += edep;
431  if (time < tMax_)
432  etotG[idx] += edep;
433  if (time > tCut_) {
434  h_edep1_[idx]->Fill(edep);
435  h_rz1_->Fill(std::abs(point.z()), point.perp());
436  }
437  }
438  }
439 
440  if (indx < 2) {
441  h_etot_[indx]->Fill(etot[indx]);
442  h_etotg_[indx]->Fill(etotG[indx]);
443  if (indx == 0)
444  h_hit_[indx]->Fill(double(nEB));
445  else
446  h_hit_[indx]->Fill(double(nEE));
447  } else {
448  h_hit_[2]->Fill(double(nHB));
449  h_hit_[3]->Fill(double(nHE));
450  h_hit_[4]->Fill(double(nHO));
451  h_hit_[5]->Fill(double(nHF));
452  for (int idx = 2; idx < nCalo_; idx++) {
453  h_etot_[idx]->Fill(etot[idx]);
454  h_etotg_[idx]->Fill(etotG[idx]);
455  }
456  }
457 
458  edm::LogVerbatim("HitStudy") << "CaloSimHitAnalysis::analyzeHits: EB " << nEB << " EE " << nEE << " HB " << nHB
459  << " HE " << nHE << " HO " << nHO << " HF " << nHF << " Bad " << nBad << " All " << nHit
460  << " Reduced " << hitMap.size();
461 }
462 
463 void CaloSimHitAnalysis::analyzePassiveHits(std::vector<PassiveHit>& hits) {
464  const std::string active = "Active";
465  const std::string sensor = "Sensor";
466  std::map<std::pair<std::string, uint32_t>, int> hitx;
467  std::map<int, int> tracks;
468  unsigned int passive1(0), passive2(0);
469  for (auto& hit : hits) {
470  std::string name = hit.vname();
471  std::pair<std::string, uint32_t> volume = std::make_pair(name, (hit.id() % 1000000));
472  auto itr = hitx.find(volume);
473  if (itr == hitx.end())
474  hitx[volume] = 1;
475  else
476  ++(itr->second);
477  auto ktr = tracks.find(hit.trackId());
478  if (ktr == tracks.end())
479  tracks[hit.trackId()] = 1;
480  else
481  ++(ktr->second);
482  h_edepp_->Fill(hit.energy());
483  h_timep_->Fill(hit.time());
484  h_stepp_->Fill(hit.stepLength());
485  if ((name.find(active) != std::string::npos) || (name.find(sensor) != std::string::npos)) {
486  unsigned idet = detNames_.size();
487  for (unsigned int k = 0; k < detNames_.size(); ++k) {
488  if (name.find(detNames_[k]) != std::string::npos) {
489  idet = k;
490  break;
491  }
492  }
493  if (idet < detNames_.size()) {
494  h_edepTk_[idet]->Fill(hit.energy());
495  h_timeTk_[idet]->Fill(hit.time());
496  }
497  }
498 
499  if ((allSteps_ / 100) % 10 > 0) {
500  uint32_t id = hit.id();
501  if (DetId(id).det() == DetId::Hcal) {
502  HcalDetId hid = HcalDetId(id);
503  int indx = (100 * hid.ietaAbs() + hid.depth());
504  auto itr = etaDepth_.find(indx);
505 #ifdef EDM_ML_DEBUG
506  edm::LogVerbatim("HitStudy") << "CaloSimHitAnalysis::ID: " << hid << " Index " << indx << " Iterator "
507  << (itr != etaDepth_.end());
508 #endif
509  ++passive1;
510  if (itr != etaDepth_.end()) {
511  uint32_t ipos = itr->second;
512  double rr = std::sqrt(hit.x() * hit.x() + hit.y() * hit.y());
513  if (ipos < h_rzH_.size()) {
514  h_rzH_[ipos]->Fill(hit.z(), rr);
515  ++passive2;
516  }
517  }
518  }
519  }
520  }
521  h_hitp_->Fill(hitx.size());
522  h_trackp_->Fill(tracks.size());
523  edm::LogVerbatim("HitStudy") << "CaloSimHitAnalysis::analyzPassiveHits: Total " << hits.size() << " Cells "
524  << hitx.size() << " Tracks " << tracks.size() << " Passive " << passive1 << ":"
525  << passive2;
526 }
527 
528 //define this as a plug-in
RunNumber_t run() const
Definition: EventID.h:38
void analyzeHits(std::vector< PCaloHit > &, int)
static const std::string kSharedResource
Definition: TFileService.h:76
Log< level::Info, true > LogVerbatim
const CaloSubdetectorGeometry * getSubdetectorGeometry(const DetId &id) const
access the subdetector geometry for the given subdetector directly
Definition: CaloGeometry.cc:34
EventNumber_t event() const
Definition: EventID.h:40
constexpr int ietaAbs() const
get the absolute value of the cell ieta
Definition: HcalDetId.h:148
edm::EDGetTokenT< edm::PassiveHitContainer > tok_passive_
T perp() const
Definition: PV3DBase.h:69
ParameterDescriptionBase * addUntracked(U const &iLabel, T const &value)
uint16_t *__restrict__ id
std::map< int, unsigned int > etaDepth_
void analyzePassiveHits(std::vector< PassiveHit > &hits)
TH1F * h_edepT1_[nCalo_]
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
const std::vector< std::string > hitLab_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
CaloSimHitAnalysis(const edm::ParameterSet &ps)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
TH1F * h_timeT_[nCalo_]
TH1F * h_edepT_[nCalo_]
std::vector< TH1F * > h_edepTk_
TH1F * h_time_[nCalo_]
Geom::Phi< T > phi() const
Definition: PV3DBase.h:66
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
T * make(const Args &...args) const
make new ROOT object
Definition: TFileService.h:64
auto const & tracks
cannot be loose
std::vector< TH1F * > h_timeTk_
const std::string names[nVars_]
const std::string g4Label_
const std::vector< std::string > detNames_
bool getData(T &iHolder) const
Definition: EventSetup.h:128
U second(std::pair< T, U > const &p)
constexpr HcalSubdetector subdet() const
get the subdetector
Definition: HcalDetId.h:138
std::vector< edm::EDGetTokenT< edm::PCaloHitContainer > > toks_calo_
const edm::ESGetToken< CaloGeometry, CaloGeometryRecord > tokGeom_
void endRun(edm::Run const &, edm::EventSetup const &) override
TH1F * h_edep_[nCalo_]
T sqrt(T t)
Definition: SSEVec.h:19
TH1F * h_edepHad_[nCalo_]
T z() const
Definition: PV3DBase.h:61
bool isAvailable() const
Definition: Service.h:40
GlobalPoint getPosition(const DetId &id) const
Get the position of a given detector id.
Definition: CaloGeometry.cc:50
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
static constexpr int nCalo_
static void unpackHcalIndex(const uint32_t &idx, int &det, int &z, int &depth, int &eta, int &phi, int &lay)
const HcalGeometry * hcalGeom_
TH1F * h_edep1_[nCalo_]
ParameterDescriptionBase * add(U const &iLabel, T const &value)
bool isValid() const
Definition: HandleBase.h:70
GlobalPoint getPosition(const DetId &id) const
#define M_PI
unsigned int id
const CaloGeometry * caloGeometry_
Definition: DetId.h:17
TH1F * h_etot_[nCalo_]
void beginRun(edm::Run const &, edm::EventSetup const &) override
const std::vector< double > timeSliceUnit_
void add(std::string const &label, ParameterSetDescription const &psetDescription)
std::vector< TH2F * > h_rzH_
T eta() const
Definition: PV3DBase.h:73
void analyze(edm::Event const &, edm::EventSetup const &) override
edm::EventID id() const
Definition: EventBase.h:59
tuple tfile
Definition: compare.py:324
constexpr int depth() const
get the tower depth
Definition: HcalDetId.h:164
void beginJob() override
ESGetTokenH3DDVariant esConsumes(std::string const &Reccord, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
TH1F * h_etotg_[nCalo_]
*vegas h *****************************************************used in the default bin number in original ***version of VEGAS is ***a higher bin number might help to derive a more precise ***grade subtle point
Definition: invegas.h:5
~CaloSimHitAnalysis() override
Definition: Run.h:45
TH1F * h_edepEM_[nCalo_]