CMS 3D CMS Logo

CaloSimHitAnalysis.cc
Go to the documentation of this file.
3 
13 
16 
22 
27 
32 
33 #include <TH1F.h>
34 #include <TH2F.h>
35 
36 #include <memory>
37 #include <iostream>
38 #include <fstream>
39 #include <vector>
40 #include <map>
41 #include <string>
42 
43 //#define EDM_ML_DEBUG
44 
45 class CaloSimHitAnalysis : public edm::one::EDAnalyzer<edm::one::WatchRuns, edm::one::SharedResources> {
46 public:
48  ~CaloSimHitAnalysis() override {}
49  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
50 
51 protected:
52  void beginJob() override {}
53  void analyze(edm::Event const&, edm::EventSetup const&) override;
54  void beginRun(edm::Run const&, edm::EventSetup const&) override {}
55  void endRun(edm::Run const&, edm::EventSetup const&) override {}
56 
57  void analyzeHits(std::vector<PCaloHit>&, int);
58  void analyzePassiveHits(std::vector<PassiveHit>& hits);
59 
60 private:
62  const std::vector<std::string> hitLab_;
63  const std::vector<double> timeSliceUnit_;
65  const bool testNumber_, passive_;
66  const int allSteps_;
67  const std::vector<std::string> detNames_;
69  const std::vector<edm::EDGetTokenT<edm::PCaloHitContainer> > toks_calo_;
71 
74 
75  static constexpr int nCalo_ = 6;
80  TH2F *h_rz_, *h_rz1_, *h_etaphi_;
82  std::vector<TH1F*> h_edepTk_, h_timeTk_;
83  std::vector<TH2F*> h_rzH_;
84  std::map<int, unsigned int> etaDepth_;
85 };
86 
88  : g4Label_(ps.getUntrackedParameter<std::string>("moduleLabel", "g4SimHits")),
89  hitLab_(ps.getParameter<std::vector<std::string> >("hitCollection")),
90  timeSliceUnit_(ps.getParameter<std::vector<double> >("timeSliceUnit")),
91  maxEnergy_(ps.getUntrackedParameter<double>("maxEnergy", 250.0)),
92  maxTime_(ps.getUntrackedParameter<double>("maxTime", 1000.0)),
93  tMax_(ps.getUntrackedParameter<double>("timeCut", 100.0)),
94  tScale_(ps.getUntrackedParameter<double>("timeScale", 1.0)),
95  tCut_(ps.getUntrackedParameter<double>("timeThreshold", 15.0)),
96  testNumber_(ps.getUntrackedParameter<bool>("testNumbering", false)),
97  passive_(ps.getUntrackedParameter<bool>("passiveHits", false)),
98  allSteps_(ps.getUntrackedParameter<int>("allSteps", 100)),
99  detNames_(ps.getUntrackedParameter<std::vector<std::string> >("detNames")),
101  toks_calo_{edm::vector_transform(
102  hitLab_,
103  [this](const std::string& name) { return consumes<edm::PCaloHitContainer>(edm::InputTag{g4Label_, name}); })},
104  tok_passive_(consumes<edm::PassiveHitContainer>(edm::InputTag(g4Label_, "AllPassiveHits"))) {
105  usesResource(TFileService::kSharedResource);
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 
281  caloGeometry_ = &set.getData(tokGeom_);
283 
284  for (unsigned int i = 0; i < toks_calo_.size(); i++) {
285  const edm::Handle<edm::PCaloHitContainer>& hitsCalo = e.getHandle(toks_calo_[i]);
286  bool getHits = (hitsCalo.isValid());
287 #ifdef EDM_ML_DEBUG
288  edm::LogVerbatim("HitStudy") << "CaloSimHitAnalysis: Input flags Hits[" << i << "]: " << getHits;
289 #endif
290  if (getHits) {
291  std::vector<PCaloHit> caloHits;
292  caloHits.insert(caloHits.end(), hitsCalo->begin(), hitsCalo->end());
293 #ifdef EDM_ML_DEBUG
294  edm::LogVerbatim("HitStudy") << "CaloSimHitAnalysis: Hit buffer [" << i << "] " << caloHits.size();
295 #endif
296  analyzeHits(caloHits, i);
297  }
298  }
299 
300  if (passive_) {
301  const edm::Handle<edm::PassiveHitContainer>& hitsPassive = e.getHandle(tok_passive_);
302  bool getHits = (hitsPassive.isValid());
303 #ifdef EDM_ML_DEBUG
304  edm::LogVerbatim("HitStudy") << "CaloSimHitAnalysis: Passive: " << getHits;
305 #endif
306  if (getHits) {
307  std::vector<PassiveHit> passiveHits;
308  passiveHits.insert(passiveHits.end(), hitsPassive->begin(), hitsPassive->end());
309 #ifdef EDM_ML_DEBUG
310  edm::LogVerbatim("HitStudy") << "CaloSimHitAnalysis: Passive Hit buffer " << passiveHits.size();
311 #endif
312  analyzePassiveHits(passiveHits);
313  }
314  }
315 }
316 
317 void CaloSimHitAnalysis::analyzeHits(std::vector<PCaloHit>& hits, int indx) {
318  int nHit = hits.size();
319  int nHB = 0, nHE = 0, nHO = 0, nHF = 0, nEB = 0, nEE = 0, nBad = 0;
320 #ifdef EDM_ML_DEBUG
321  int iHit = 0;
322 #endif
323 
324  std::map<CaloHitID, std::pair<double, double> > hitMap;
325  double etot[nCalo_], etotG[nCalo_];
326  for (unsigned int k = 0; k < nCalo_; ++k) {
327  etot[k] = etotG[k] = 0;
328  }
329  for (const auto& hit : hits) {
330  double edep = hit.energy();
331  double time = tScale_ * hit.time();
332  uint32_t id = hit.id();
333  int itra = hit.geantTrackId();
334  double edepEM = hit.energyEM();
335  double edepHad = hit.energyHad();
336  int idx(-1);
337  if (indx != 2) {
338  idx = indx;
339  if (indx == 0)
340  ++nEB;
341  else
342  ++nEE;
343  } else {
344  int subdet(0);
345  if (testNumber_) {
346  int ieta(0), phi(0), z(0), lay(0), depth(0);
347  HcalTestNumbering::unpackHcalIndex(id, subdet, z, depth, ieta, phi, lay);
348  id = HcalDetId(static_cast<HcalSubdetector>(subdet), z * ieta, phi, depth).rawId();
349  } else {
350  subdet = HcalDetId(id).subdet();
351  }
352  if (subdet == static_cast<int>(HcalBarrel)) {
353  idx = indx;
354  nHB++;
355  } else if (subdet == static_cast<int>(HcalEndcap)) {
356  idx = indx + 1;
357  nHE++;
358  } else if (subdet == static_cast<int>(HcalOuter)) {
359  idx = indx + 2;
360  nHO++;
361  } else if (subdet == static_cast<int>(HcalForward)) {
362  idx = indx + 3;
363  nHF++;
364  }
365  }
366 #ifdef EDM_ML_DEBUG
367  edm::LogVerbatim("HitStudy") << "Hit[" << iHit << ":" << nHit << ":" << idx << "] E " << edep << ":" << edepEM
368  << ":" << edepHad << " T " << time << " itra " << itra << " ID " << std::hex << id
369  << std::dec;
370  ++iHit;
371 #endif
372  if (idx >= 0) {
373  CaloHitID hid(id, time, itra, 0, timeSliceUnit_[indx]);
374  auto itr = hitMap.find(hid);
375  if (itr == hitMap.end())
376  hitMap[hid] = std::make_pair(time, edep);
377  else
378  ((itr->second).second) += edep;
379  h_edepT_[idx]->Fill(edep);
380  h_timeT_[idx]->Fill(time);
381  if (edepEM > 0)
382  h_edepEM_[idx]->Fill(edepEM);
383  if (edepHad > 0)
384  h_edepHad_[idx]->Fill(edepHad);
385  if (time > tCut_)
386  h_edepT1_[idx]->Fill(edep);
387  } else {
388  ++nBad;
389  }
390  }
391 
392  //Now make plots
393  for (auto itr = hitMap.begin(); itr != hitMap.end(); ++itr) {
394  int idx = -1;
396  DetId id((itr->first).unitID());
397 #ifdef EDM_ML_DEBUG
398  edm::LogVerbatim("HitStudy") << "Index " << indx << " Geom " << caloGeometry_ << ":" << hcalGeom_ << " "
399  << std::hex << id.rawId() << std::dec;
400 #endif
401  if (indx != 2) {
402  idx = indx;
404  } else {
405  int subdet = id.subdetId();
406  if (subdet == static_cast<int>(HcalBarrel)) {
407  idx = indx;
408  } else if (subdet == static_cast<int>(HcalEndcap)) {
409  idx = indx + 1;
410  } else if (subdet == static_cast<int>(HcalOuter)) {
411  idx = indx + 2;
412  } else if (subdet == static_cast<int>(HcalForward)) {
413  idx = indx + 3;
414  }
415  point = hcalGeom_->getPosition(id);
416  }
417  double edep = (itr->second).second;
418  double time = (itr->second).first;
419 #ifdef EDM_ML_DEBUG
420  edm::LogVerbatim("HitStudy") << "Index " << idx << ":" << nCalo_ << " Point " << point << " E " << edep << " T "
421  << time;
422 #endif
423  if (idx >= 0) {
424  h_time_[idx]->Fill(time);
425  h_edep_[idx]->Fill(edep);
426  h_rr_[idx]->Fill(point.perp());
427  h_zz_[idx]->Fill(point.z());
428  h_eta_[idx]->Fill(point.eta());
429  h_phi_[idx]->Fill(point.phi());
430  h_rz_->Fill(std::abs(point.z()), point.perp());
431  h_etaphi_->Fill(std::abs(point.eta()), std::abs(point.phi()));
432  etot[idx] += edep;
433  if (time < tMax_)
434  etotG[idx] += edep;
435  if (time > tCut_) {
436  h_edep1_[idx]->Fill(edep);
437  h_rz1_->Fill(std::abs(point.z()), point.perp());
438  }
439  }
440  }
441 
442  if (indx < 2) {
443  h_etot_[indx]->Fill(etot[indx]);
444  h_etotg_[indx]->Fill(etotG[indx]);
445  if (indx == 0)
446  h_hit_[indx]->Fill(double(nEB));
447  else
448  h_hit_[indx]->Fill(double(nEE));
449  } else {
450  h_hit_[2]->Fill(double(nHB));
451  h_hit_[3]->Fill(double(nHE));
452  h_hit_[4]->Fill(double(nHO));
453  h_hit_[5]->Fill(double(nHF));
454  for (int idx = 2; idx < nCalo_; idx++) {
455  h_etot_[idx]->Fill(etot[idx]);
456  h_etotg_[idx]->Fill(etotG[idx]);
457  }
458  }
459 
460  edm::LogVerbatim("HitStudy") << "CaloSimHitAnalysis::analyzeHits: EB " << nEB << " EE " << nEE << " HB " << nHB
461  << " HE " << nHE << " HO " << nHO << " HF " << nHF << " Bad " << nBad << " All " << nHit
462  << " Reduced " << hitMap.size();
463 }
464 
465 void CaloSimHitAnalysis::analyzePassiveHits(std::vector<PassiveHit>& hits) {
466  const std::string active = "Active";
467  const std::string sensor = "Sensor";
468  std::map<std::pair<std::string, uint32_t>, int> hitx;
469  std::map<int, int> tracks;
470  unsigned int passive1(0), passive2(0);
471  for (auto& hit : hits) {
472  std::string name = hit.vname();
473  std::pair<std::string, uint32_t> volume = std::make_pair(name, (hit.id() % 1000000));
474  auto itr = hitx.find(volume);
475  if (itr == hitx.end())
476  hitx[volume] = 1;
477  else
478  ++(itr->second);
479  auto ktr = tracks.find(hit.trackId());
480  if (ktr == tracks.end())
481  tracks[hit.trackId()] = 1;
482  else
483  ++(ktr->second);
484  h_edepp_->Fill(hit.energy());
485  h_timep_->Fill(hit.time());
486  h_stepp_->Fill(hit.stepLength());
487  if ((name.find(active) != std::string::npos) || (name.find(sensor) != std::string::npos)) {
488  unsigned idet = detNames_.size();
489  for (unsigned int k = 0; k < detNames_.size(); ++k) {
490  if (name.find(detNames_[k]) != std::string::npos) {
491  idet = k;
492  break;
493  }
494  }
495  if (idet < detNames_.size()) {
496  h_edepTk_[idet]->Fill(hit.energy());
497  h_timeTk_[idet]->Fill(hit.time());
498  }
499  }
500 
501  if ((allSteps_ / 100) % 10 > 0) {
502  uint32_t id = hit.id();
503  if (DetId(id).det() == DetId::Hcal) {
504  HcalDetId hid = HcalDetId(id);
505  int indx = (100 * hid.ietaAbs() + hid.depth());
506  auto itr = etaDepth_.find(indx);
507 #ifdef EDM_ML_DEBUG
508  edm::LogVerbatim("HitStudy") << "CaloSimHitAnalysis::ID: " << hid << " Index " << indx << " Iterator "
509  << (itr != etaDepth_.end());
510 #endif
511  ++passive1;
512  if (itr != etaDepth_.end()) {
513  uint32_t ipos = itr->second;
514  double rr = std::sqrt(hit.x() * hit.x() + hit.y() * hit.y());
515  if (ipos < h_rzH_.size()) {
516  h_rzH_[ipos]->Fill(hit.z(), rr);
517  ++passive2;
518  }
519  }
520  }
521  }
522  }
523  h_hitp_->Fill(hitx.size());
524  h_trackp_->Fill(tracks.size());
525  edm::LogVerbatim("HitStudy") << "CaloSimHitAnalysis::analyzPassiveHits: Total " << hits.size() << " Cells "
526  << hitx.size() << " Tracks " << tracks.size() << " Passive " << passive1 << ":"
527  << passive2;
528 }
529 
530 //define this as a plug-in
void analyzeHits(std::vector< PCaloHit > &, int)
static const std::string kSharedResource
Definition: TFileService.h:76
Log< level::Info, true > LogVerbatim
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
std::map< int, unsigned int > etaDepth_
void analyzePassiveHits(std::vector< PassiveHit > &hits)
TH1F * h_edepT1_[nCalo_]
constexpr int ietaAbs() const
get the absolute value of the cell ieta
Definition: HcalDetId.h:148
const std::vector< std::string > hitLab_
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_]
std::vector< TH1F * > h_timeTk_
auto vector_transform(std::vector< InputType > const &input, Function predicate) -> std::vector< typename std::remove_cv< typename std::remove_reference< decltype(predicate(input.front()))>::type >::type >
Definition: transform.h:11
const std::string names[nVars_]
const std::string g4Label_
const std::vector< std::string > detNames_
U second(std::pair< T, U > const &p)
constexpr HcalSubdetector subdet() const
get the subdetector
Definition: HcalDetId.h:138
const edm::ESGetToken< CaloGeometry, CaloGeometryRecord > tokGeom_
void endRun(edm::Run const &, edm::EventSetup const &) override
GlobalPoint getPosition(const DetId &id) const
Get the position of a given detector id.
Definition: CaloGeometry.cc:50
const edm::EDGetTokenT< edm::PassiveHitContainer > tok_passive_
TH1F * h_edep_[nCalo_]
T sqrt(T t)
Definition: SSEVec.h:23
TH1F * h_edepHad_[nCalo_]
Definition: tfile.py:1
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_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
TH1F * h_edep1_[nCalo_]
#define M_PI
unsigned int id
const CaloGeometry * caloGeometry_
Definition: DetId.h:17
const std::vector< edm::EDGetTokenT< edm::PCaloHitContainer > > toks_calo_
TH1F * h_etot_[nCalo_]
void beginRun(edm::Run const &, edm::EventSetup const &) override
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
const std::vector< double > timeSliceUnit_
void add(std::string const &label, ParameterSetDescription const &psetDescription)
std::vector< TH2F * > h_rzH_
bool isValid() const
Definition: HandleBase.h:70
void analyze(edm::Event const &, edm::EventSetup const &) override
GlobalPoint getPosition(const DetId &id) const
void beginJob() override
const CaloSubdetectorGeometry * getSubdetectorGeometry(const DetId &id) const
access the subdetector geometry for the given subdetector directly
Definition: CaloGeometry.cc:34
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_]
constexpr int depth() const
get the tower depth
Definition: HcalDetId.h:164