CMS 3D CMS Logo

CaloSimHitStudy.cc
Go to the documentation of this file.
3 
7 
14 
17 
23 
30 
35 
36 #include <TH1F.h>
37 
38 #include <memory>
39 #include <iostream>
40 #include <fstream>
41 #include <vector>
42 #include <map>
43 #include <string>
44 
45 //#define EDM_ML_DEBUG
46 
47 class CaloSimHitStudy : public edm::one::EDAnalyzer<edm::one::WatchRuns, edm::one::SharedResources> {
48 public:
50  ~CaloSimHitStudy() override {}
51  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
52 
53 protected:
54  void beginJob() override {}
55  void analyze(edm::Event const&, edm::EventSetup const&) override;
56  void beginRun(edm::Run const&, edm::EventSetup const&) override {}
57  void endRun(edm::Run const&, edm::EventSetup const&) override {}
58 
59  void analyzeHits(std::vector<PCaloHit>&, int);
61  void analyzeHits(std::vector<PCaloHit>&, std::vector<PCaloHit>&, std::vector<PCaloHit>&);
62 
63 private:
65  const std::vector<std::string> hitLab_;
66  const double maxEnergy_, tmax_, eMIP_;
67  const bool storeRL_, testNumber_;
70  const std::vector<edm::EDGetTokenT<edm::PCaloHitContainer>> toks_calo_;
71  std::vector<edm::EDGetTokenT<edm::PSimHitContainer>> toks_track_;
72  std::vector<edm::EDGetTokenT<edm::PSimHitContainer>> toks_tkHigh_;
73  std::vector<edm::EDGetTokenT<edm::PSimHitContainer>> toks_tkLow_;
74 
77 
78  const std::vector<std::string> muonLab_ = {"MuonRPCHits", "MuonCSCHits", "MuonDTHits", "MuonGEMHits"};
79  const std::vector<std::string> tkHighLab_ = {"TrackerHitsPixelBarrelHighTof",
80  "TrackerHitsPixelEndcapHighTof",
81  "TrackerHitsTECHighTof",
82  "TrackerHitsTIBHighTof",
83  "TrackerHitsTIDHighTof",
84  "TrackerHitsTOBHighTof"};
85  const std::vector<std::string> tkLowLab_ = {"TrackerHitsPixelBarrelLowTof",
86  "TrackerHitsPixelEndcapLowTof",
87  "TrackerHitsTECLowTof",
88  "TrackerHitsTIBLowTof",
89  "TrackerHitsTIDLowTof",
90  "TrackerHitsTOBLowTof"};
91 
92  static constexpr int nCalo_ = 9;
93  static constexpr int nTrack_ = 16;
99 };
100 
102  : g4Label_(ps.getUntrackedParameter<std::string>("ModuleLabel")),
103  hitLab_(ps.getUntrackedParameter<std::vector<std::string>>("CaloCollection")),
104  maxEnergy_(ps.getUntrackedParameter<double>("MaxEnergy", 200.0)),
105  tmax_(ps.getUntrackedParameter<double>("TimeCut", 100.0)),
106  eMIP_(ps.getUntrackedParameter<double>("MIPCut", 0.70)),
107  storeRL_(ps.getUntrackedParameter<bool>("StoreRL", false)),
108  testNumber_(ps.getUntrackedParameter<bool>("TestNumbering", true)),
110  tok_evt_(consumes<edm::HepMCProduct>(
111  edm::InputTag(ps.getUntrackedParameter<std::string>("SourceLabel", "VtxSmeared")))),
112  toks_calo_{edm::vector_transform(hitLab_, [this](const std::string& name) {
113  return consumes<edm::PCaloHitContainer>(edm::InputTag{g4Label_, name});
114  })} {
115  usesResource(TFileService::kSharedResource);
116 
117  for (unsigned i = 0; i != muonLab_.size(); i++)
118  toks_track_.emplace_back(consumes<edm::PSimHitContainer>(edm::InputTag(g4Label_, muonLab_[i])));
119  for (unsigned i = 0; i != tkHighLab_.size(); i++)
120  toks_tkHigh_.emplace_back(consumes<edm::PSimHitContainer>(edm::InputTag(g4Label_, tkHighLab_[i])));
121  for (unsigned i = 0; i != tkLowLab_.size(); i++)
122  toks_tkLow_.emplace_back(consumes<edm::PSimHitContainer>(edm::InputTag(g4Label_, tkLowLab_[i])));
123 
124  edm::LogVerbatim("HitStudy") << "Module Label: " << g4Label_ << " Hits: " << hitLab_[0] << ", " << hitLab_[1]
125  << ", " << hitLab_[2] << ", " << hitLab_[3] << " MaxEnergy: " << maxEnergy_
126  << " Tmax: " << tmax_ << " MIP Cut: " << eMIP_;
127 
129 
130  if (!tfile.isAvailable())
131  throw cms::Exception("BadConfig") << "TFileService unavailable: "
132  << "please add it to config file";
133  char name[20], title[120];
134  sprintf(title, "Incident PT (GeV)");
135  ptInc_ = tfile->make<TH1F>("PtInc", title, 1000, 0., maxEnergy_);
136  ptInc_->GetXaxis()->SetTitle(title);
137  ptInc_->GetYaxis()->SetTitle("Events");
138  sprintf(title, "Incident Energy (GeV)");
139  eneInc_ = tfile->make<TH1F>("EneInc", title, 1000, 0., maxEnergy_);
140  eneInc_->GetXaxis()->SetTitle(title);
141  eneInc_->GetYaxis()->SetTitle("Events");
142  sprintf(title, "Incident #eta");
143  etaInc_ = tfile->make<TH1F>("EtaInc", title, 200, -5., 5.);
144  etaInc_->GetXaxis()->SetTitle(title);
145  etaInc_->GetYaxis()->SetTitle("Events");
146  sprintf(title, "Incident #phi");
147  phiInc_ = tfile->make<TH1F>("PhiInc", title, 200, -3.1415926, 3.1415926);
148  phiInc_->GetXaxis()->SetTitle(title);
149  phiInc_->GetYaxis()->SetTitle("Events");
150 #ifdef EDM_ML_DEBUG
151  edm::LogVerbatim("HitStudy") << "CaloSimHitStudy: Completed defining histos for incident particle";
152 #endif
153  std::string dets[nCalo_] = {"EB", "EB(APD)", "EB(ATJ)", "EE", "ES", "HB", "HE", "HO", "HF"};
154  double nhcMax[nCalo_] = {40000., 2000., 2000., 40000., 10000., 10000., 10000., 2000., 10000.};
155  for (int i = 0; i < nCalo_; i++) {
156  sprintf(name, "Hit%d", i);
157  sprintf(title, "Number of hits in %s", dets[i].c_str());
158  hit_[i] = tfile->make<TH1F>(name, title, 1000, 0., nhcMax[i]);
159  hit_[i]->GetXaxis()->SetTitle(title);
160  hit_[i]->GetYaxis()->SetTitle("Events");
161  sprintf(name, "Time%d", i);
162  sprintf(title, "Time of the hit (ns) in %s", dets[i].c_str());
163  time_[i] = tfile->make<TH1F>(name, title, 1000, 0., 1000.);
164  time_[i]->GetXaxis()->SetTitle(title);
165  time_[i]->GetYaxis()->SetTitle("Hits");
166  sprintf(name, "TimeAll%d", i);
167  sprintf(title, "Hit time (ns) in %s with no check on Track ID", dets[i].c_str());
168  timeAll_[i] = tfile->make<TH1F>(name, title, 1000, 0., 1000.);
169  timeAll_[i]->GetXaxis()->SetTitle(title);
170  timeAll_[i]->GetYaxis()->SetTitle("Hits");
171  double ymax = maxEnergy_;
172  if (i == 1 || i == 2 || i == 4)
173  ymax = 1.0;
174  else if (i > 4 && i < 8)
175  ymax = 10.0;
176  sprintf(name, "Edep%d", i);
177  sprintf(title, "Energy deposit (GeV) in %s", dets[i].c_str());
178  edep_[i] = tfile->make<TH1F>(name, title, 5000, 0., ymax);
179  edep_[i]->GetXaxis()->SetTitle(title);
180  edep_[i]->GetYaxis()->SetTitle("Hits");
181  sprintf(name, "EdepEM%d", i);
182  sprintf(title, "Energy deposit (GeV) by EM particles in %s", dets[i].c_str());
183  edepEM_[i] = tfile->make<TH1F>(name, title, 5000, 0., ymax);
184  edepEM_[i]->GetXaxis()->SetTitle(title);
185  edepEM_[i]->GetYaxis()->SetTitle("Hits");
186  sprintf(name, "EdepHad%d", i);
187  sprintf(title, "Energy deposit (GeV) by hadrons in %s", dets[i].c_str());
188  edepHad_[i] = tfile->make<TH1F>(name, title, 5000, 0., ymax);
189  edepHad_[i]->GetXaxis()->SetTitle(title);
190  edepHad_[i]->GetYaxis()->SetTitle("Hits");
191  sprintf(name, "Etot%d", i);
192  sprintf(title, "Total energy deposit (GeV) in %s", dets[i].c_str());
193  etot_[i] = tfile->make<TH1F>(name, title, 5000, 0., ymax);
194  etot_[i]->GetXaxis()->SetTitle(title);
195  etot_[i]->GetYaxis()->SetTitle("Events");
196  sprintf(name, "EtotG%d", i);
197  sprintf(title, "Total energy deposit (GeV) in %s (t < 100 ns)", dets[i].c_str());
198  etotg_[i] = tfile->make<TH1F>(name, title, 5000, 0., ymax);
199  etotg_[i]->GetXaxis()->SetTitle(title);
200  etotg_[i]->GetYaxis()->SetTitle("Events");
201  sprintf(name, "eta%d", i);
202  sprintf(title, "#eta of hit point in %s", dets[i].c_str());
203  eta_[i] = tfile->make<TH1F>(name, title, 100, -5.0, 5.0);
204  eta_[i]->GetXaxis()->SetTitle(title);
205  eta_[i]->GetYaxis()->SetTitle("Hits");
206  sprintf(name, "phi%d", i);
207  sprintf(title, "#phi of hit point in %s", dets[i].c_str());
208  phi_[i] = tfile->make<TH1F>(name, title, 100, -M_PI, M_PI);
209  phi_[i]->GetXaxis()->SetTitle(title);
210  phi_[i]->GetYaxis()->SetTitle("Hits");
211  }
212 #ifdef EDM_ML_DEBUG
213  edm::LogVerbatim("HitStudy") << "CaloSimHitStudy: Completed defining histos for first level of Calorimeter";
214 #endif
215  std::string detx[9] = {"EB/EE (MIP)",
216  "HB/HE (MIP)",
217  "HB/HE/HO (MIP)",
218  "EB/EE (no MIP)",
219  "HB/HE (no MIP)",
220  "HB/HE/HO (no MIP)",
221  "EB/EE (All)",
222  "HB/HE (All)",
223  "HB/HE/HO (All)"};
224  for (int i = 0; i < 9; i++) {
225  double ymax = 1.0;
226  if (i == 0 || i == 3 || i == 6)
227  ymax = maxEnergy_;
228  sprintf(name, "EdepCal%d", i);
229  sprintf(title, "Energy deposit in %s", detx[i].c_str());
230  edepC_[i] = tfile->make<TH1F>(name, title, 5000, 0., ymax);
231  edepC_[i]->GetXaxis()->SetTitle(title);
232  edepC_[i]->GetYaxis()->SetTitle("Events");
233  sprintf(name, "EdepCalT%d", i);
234  sprintf(title, "Energy deposit (t < %f ns) in %s", tmax_, detx[i].c_str());
235  edepT_[i] = tfile->make<TH1F>(name, title, 5000, 0., ymax);
236  edepT_[i]->GetXaxis()->SetTitle(title);
237  edepT_[i]->GetYaxis()->SetTitle("Events");
238  }
239 #ifdef EDM_ML_DEBUG
240  edm::LogVerbatim("HitStudy") << "CaloSimHitStudy: Completed defining histos for second level of Calorimeter";
241 #endif
242  hitLow = tfile->make<TH1F>("HitLow", "Number of hits in Track (Low)", 1000, 0, 20000.);
243  hitLow->GetXaxis()->SetTitle("Number of hits in Track (Low)");
244  hitLow->GetYaxis()->SetTitle("Events");
245  hitHigh = tfile->make<TH1F>("HitHigh", "Number of hits in Track (High)", 1000, 0, 5000.);
246  hitHigh->GetXaxis()->SetTitle("Number of hits in Track (High)");
247  hitHigh->GetYaxis()->SetTitle("Events");
248  hitMu = tfile->make<TH1F>("HitMu", "Number of hits in Track (Muon)", 1000, 0, 2000.);
249  hitMu->GetXaxis()->SetTitle("Number of hits in Muon");
250  hitMu->GetYaxis()->SetTitle("Events");
251 #ifdef EDM_ML_DEBUG
252  edm::LogVerbatim("HitStudy") << "CaloSimHitStudy: Completed defining histos for general tracking hits";
253 #endif
254  std::string dett[nTrack_] = {"Pixel Barrel (High)",
255  "Pixel Endcap (High)",
256  "TEC (High)",
257  "TIB (High)",
258  "TID (High)",
259  "TOB (High)",
260  "Pixel Barrel (Low)",
261  "Pixel Endcap (Low)",
262  "TEC (Low)",
263  "TIB (Low)",
264  "TID (Low)",
265  "TOB (Low)",
266  "RPC",
267  "CSC",
268  "DT",
269  "GEM"};
270  double nhtMax[nTrack_] = {
271  500., 500., 1000., 1000., 500., 1000., 5000., 2000., 10000., 5000., 2000., 5000., 500., 1000., 1000., 500.};
272  for (int i = 0; i < nTrack_; i++) {
273  sprintf(name, "HitTk%d", i);
274  sprintf(title, "Number of hits in %s", dett[i].c_str());
275  hitTk_[i] = tfile->make<TH1F>(name, title, 1000, 0., nhtMax[i]);
276  hitTk_[i]->GetXaxis()->SetTitle(title);
277  hitTk_[i]->GetYaxis()->SetTitle("Events");
278  sprintf(name, "TimeTk%d", i);
279  sprintf(title, "Time of the hit (ns) in %s", dett[i].c_str());
280  tofTk_[i] = tfile->make<TH1F>(name, title, 1000, 0., 200.);
281  tofTk_[i]->GetXaxis()->SetTitle(title);
282  tofTk_[i]->GetYaxis()->SetTitle("Hits");
283  sprintf(name, "EdepTk%d", i);
284  sprintf(title, "Energy deposit (GeV) in %s", dett[i].c_str());
285  double ymax = (i < 12) ? 0.25 : 0.005;
286  edepTk_[i] = tfile->make<TH1F>(name, title, 250, 0., ymax);
287  edepTk_[i]->GetXaxis()->SetTitle(title);
288  edepTk_[i]->GetYaxis()->SetTitle("Hits");
289  }
290 #ifdef EDM_ML_DEBUG
291  edm::LogVerbatim("HitStudy") << "CaloSimHitStudy: Completed defining histos for SimHit objects";
292 #endif
293 }
294 
297  std::vector<std::string> calonames = {"EcalHitsEB", "EcalHitsEE", "EcalHitsES", "HcalHits"};
298  desc.addUntracked<std::string>("SourceLabel", "generatorSmeared");
299  desc.addUntracked<std::string>("ModuleLabel", "g4SimHits");
300  desc.addUntracked<std::vector<std::string>>("CaloCollection", calonames);
301  desc.addUntracked<double>("MaxEnergy", 200.0);
302  desc.addUntracked<double>("TimeCut", 100.0);
303  desc.addUntracked<double>("MIPCut", 0.70);
304  desc.addUntracked<bool>("StoreRL", false);
305  desc.addUntracked<bool>("TestNumbering", true);
306  descriptions.add("CaloSimHitStudy", desc);
307 }
308 
310  edm::LogVerbatim("HitStudy") << "CaloSimHitStudy:Run = " << e.id().run() << " Event = " << e.id().event();
311 
312  caloGeometry_ = &set.getData(tokGeom_);
314 
315  const edm::Handle<edm::HepMCProduct> EvtHandle = e.getHandle(tok_evt_);
316  const HepMC::GenEvent* myGenEvent = EvtHandle->GetEvent();
317 
318  double eInc = 0, etaInc = 0, phiInc = 0;
319  HepMC::GenEvent::particle_const_iterator p = myGenEvent->particles_begin();
320  if (p != myGenEvent->particles_end()) {
321  eInc = (*p)->momentum().e();
322  etaInc = (*p)->momentum().eta();
323  phiInc = (*p)->momentum().phi();
324  }
325  double ptInc = eInc / std::cosh(etaInc);
326  ptInc_->Fill(ptInc);
327  eneInc_->Fill(eInc);
328  etaInc_->Fill(etaInc);
329  phiInc_->Fill(phiInc);
330 
331  std::vector<PCaloHit> ebHits, eeHits, hcHits;
332  for (unsigned int i = 0; i < toks_calo_.size(); i++) {
333  bool getHits = false;
334  const edm::Handle<edm::PCaloHitContainer>& hitsCalo = e.getHandle(toks_calo_[i]);
335  if (hitsCalo.isValid())
336  getHits = true;
337 #ifdef EDM_ML_DEBUG
338  edm::LogVerbatim("HitStudy") << "CaloSimHitStudy: Input flags Hits " << getHits;
339 #endif
340  if (getHits) {
341  std::vector<PCaloHit> caloHits;
342  caloHits.insert(caloHits.end(), hitsCalo->begin(), hitsCalo->end());
343  if (i == 0)
344  ebHits.insert(ebHits.end(), hitsCalo->begin(), hitsCalo->end());
345  else if (i == 1)
346  eeHits.insert(eeHits.end(), hitsCalo->begin(), hitsCalo->end());
347  else if (i == 3)
348  hcHits.insert(hcHits.end(), hitsCalo->begin(), hitsCalo->end());
349 #ifdef EDM_ML_DEBUG
350  edm::LogVerbatim("HitStudy") << "CaloSimHitStudy: Hit buffer " << caloHits.size();
351 #endif
352  analyzeHits(caloHits, i);
353  }
354  }
355  analyzeHits(ebHits, eeHits, hcHits);
356 
357  std::vector<PSimHit> muonHits;
358  for (unsigned int i = 0; i < toks_track_.size(); i++) {
359  const edm::Handle<edm::PSimHitContainer>& hitsTrack = e.getHandle(toks_track_[i]);
360  if (hitsTrack.isValid()) {
361  muonHits.insert(muonHits.end(), hitsTrack->begin(), hitsTrack->end());
362  analyzeHits(hitsTrack, i + 12);
363  }
364  }
365  unsigned int nhmu = muonHits.size();
366  hitMu->Fill(double(nhmu));
367  std::vector<PSimHit> tkHighHits;
368  for (unsigned int i = 0; i < toks_tkHigh_.size(); i++) {
369  const edm::Handle<edm::PSimHitContainer>& hitsTrack = e.getHandle(toks_tkHigh_[i]);
370  if (hitsTrack.isValid()) {
371  tkHighHits.insert(tkHighHits.end(), hitsTrack->begin(), hitsTrack->end());
372  analyzeHits(hitsTrack, i);
373  }
374  }
375  unsigned int nhtkh = tkHighHits.size();
376  hitHigh->Fill(double(nhtkh));
377  std::vector<PSimHit> tkLowHits;
378  for (unsigned int i = 0; i < toks_tkLow_.size(); i++) {
379  const edm::Handle<edm::PSimHitContainer>& hitsTrack = e.getHandle(toks_tkLow_[i]);
380  if (hitsTrack.isValid()) {
381  tkLowHits.insert(tkLowHits.end(), hitsTrack->begin(), hitsTrack->end());
382  analyzeHits(hitsTrack, i + 6);
383  }
384  }
385  unsigned int nhtkl = tkLowHits.size();
386  hitLow->Fill(double(nhtkl));
387 }
388 
389 void CaloSimHitStudy::analyzeHits(std::vector<PCaloHit>& hits, int indx) {
390  int nHit = hits.size();
391  int nHB(0), nHE(0), nHO(0), nHF(0), nEB(0), nEBAPD(0), nEBATJ(0), nEE(0), nES(0);
392 #ifdef EDM_ML_DEBUG
393  int nBad(0);
394 #endif
395  std::map<unsigned int, double> hitMap;
396  std::vector<double> etot(nCalo_, 0), etotG(nCalo_, 0);
397  for (int i = 0; i < nHit; i++) {
398  double edep = hits[i].energy();
399  double time = hits[i].time();
400  unsigned int id = hits[i].id();
401  double edepEM = hits[i].energyEM();
402  double edepHad = hits[i].energyHad();
403  if (indx == 0) {
404  int dep = (hits[i].depth()) & PCaloHit::kEcalDepthIdMask;
405  if (dep == 1)
406  id |= 0x20000;
407  else if (dep == 2)
408  id |= 0x40000;
409  }
410  std::map<unsigned int, double>::const_iterator it = hitMap.find(id);
411  if (it == hitMap.end()) {
412  hitMap.insert(std::pair<unsigned int, double>(id, time));
413  }
414  int idx = -1;
415  if (indx != 3) {
416  if (indx == 0) {
417  if (storeRL_)
418  idx = 0;
419  else
420  idx = ((hits[i].depth()) & PCaloHit::kEcalDepthIdMask);
421  } else
422  idx = indx + 2;
423  time_[idx]->Fill(time);
424  edep_[idx]->Fill(edep);
425  edepEM_[idx]->Fill(edepEM);
426  edepHad_[idx]->Fill(edepHad);
427  if (idx == 0)
428  nEB++;
429  else if (idx == 1)
430  nEBAPD++;
431  else if (idx == 2)
432  nEBATJ++;
433  else if (idx == 3)
434  nEE++;
435  else if (idx == 4)
436  nES++;
437 #ifdef EDM_ML_DEBUG
438  else
439  nBad++;
440 #endif
441  if (indx >= 0 && indx < 3) {
442  etot[idx] += edep;
443  if (time < 100)
444  etotG[idx] += edep;
445  }
446  } else {
447  int subdet(0);
448  if (testNumber_) {
449  int ieta(0), phi(0), z(0), lay(0), depth(0);
450  HcalTestNumbering::unpackHcalIndex(id, subdet, z, depth, ieta, phi, lay);
451  } else {
452  subdet = HcalDetId(id).subdet();
453  }
454  if (subdet == static_cast<int>(HcalBarrel)) {
455  idx = indx + 2;
456  nHB++;
457  } else if (subdet == static_cast<int>(HcalEndcap)) {
458  idx = indx + 3;
459  nHE++;
460  } else if (subdet == static_cast<int>(HcalOuter)) {
461  idx = indx + 4;
462  nHO++;
463  } else if (subdet == static_cast<int>(HcalForward)) {
464  idx = indx + 5;
465  nHF++;
466  }
467  if (idx > 0) {
468  time_[idx]->Fill(time);
469  edep_[idx]->Fill(edep);
470  edepEM_[idx]->Fill(edepEM);
471  edepHad_[idx]->Fill(edepHad);
472  etot[idx] += edep;
473  if (time < 100)
474  etotG[idx] += edep;
475  } else {
476 #ifdef EDM_ML_DEBUG
477  nBad++;
478 #endif
479  }
480  }
481  }
482  if (indx < 3) {
483  etot_[indx + 2]->Fill(etot[indx + 2]);
484  etotg_[indx + 2]->Fill(etotG[indx + 2]);
485  if (indx == 0) {
486  etot_[indx]->Fill(etot[indx]);
487  etotg_[indx]->Fill(etotG[indx]);
488  etot_[indx + 1]->Fill(etot[indx + 1]);
489  etotg_[indx + 1]->Fill(etotG[indx + 1]);
490  hit_[indx]->Fill(double(nEB));
491  hit_[indx + 1]->Fill(double(nEBAPD));
492  hit_[indx + 2]->Fill(double(nEBATJ));
493  } else {
494  hit_[indx + 2]->Fill(double(nHit));
495  }
496  } else if (indx == 3) {
497  hit_[5]->Fill(double(nHB));
498  hit_[6]->Fill(double(nHE));
499  hit_[7]->Fill(double(nHO));
500  hit_[8]->Fill(double(nHF));
501  for (int idx = 5; idx < 9; idx++) {
502  etot_[idx]->Fill(etot[idx]);
503  etotg_[idx]->Fill(etotG[idx]);
504  }
505  }
506 
507 #ifdef EDM_ML_DEBUG
508  edm::LogVerbatim("HitStudy") << "CaloSimHitStudy::analyzeHits: EB " << nEB << ", " << nEBAPD << ", " << nEBATJ
509  << " EE " << nEE << " ES " << nES << " HB " << nHB << " HE " << nHE << " HO " << nHO
510  << " HF " << nHF << " Bad " << nBad << " All " << nHit << " Reduced " << hitMap.size();
511 #endif
512  std::map<unsigned int, double>::const_iterator it = hitMap.begin();
513  for (; it != hitMap.end(); it++) {
514  double time = it->second;
516  DetId id(it->first);
517  if (indx != 2)
519  int idx = -1;
520  if (indx < 3) {
521  if (indx == 0) {
522  if ((id & 0x20000) != 0)
523  idx = indx + 1;
524  else if ((id & 0x40000) != 0)
525  idx = indx + 1;
526  else
527  idx = indx;
528  } else {
529  idx = indx + 2;
530  }
531  if (idx >= 0 && idx < 5)
532  timeAll_[idx]->Fill(time);
533  } else if (indx == 3) {
534  int idx(-1), subdet(0);
535  if (testNumber_) {
536  int ieta(0), phi(0), z(0), lay(0), depth(0);
537  HcalTestNumbering::unpackHcalIndex(id.rawId(), subdet, z, depth, ieta, phi, lay);
538  } else {
539  subdet = HcalDetId(id).subdet();
540  }
541  if (subdet == static_cast<int>(HcalBarrel)) {
542  idx = indx + 2;
543  } else if (subdet == static_cast<int>(HcalEndcap)) {
544  idx = indx + 3;
545  } else if (subdet == static_cast<int>(HcalOuter)) {
546  idx = indx + 4;
547  } else if (subdet == static_cast<int>(HcalForward)) {
548  idx = indx + 5;
549  }
550  if (idx > 0) {
551  timeAll_[idx]->Fill(time);
552  eta_[idx]->Fill(point.eta());
553  phi_[idx]->Fill(point.phi());
554  }
555  }
556  }
557 }
558 
560  int nHit = 0;
561  edm::PSimHitContainer::const_iterator ihit;
562  std::string label(" ");
563  if (indx >= 0 && indx < 6)
564  label = tkHighLab_[indx];
565  else if (indx >= 6 && indx < 12)
566  label = tkLowLab_[indx - 6];
567  else if (indx >= 12 && indx < nTrack_)
568  label = muonLab_[indx - 12];
569  for (ihit = hits->begin(); ihit != hits->end(); ihit++) {
570  edepTk_[indx]->Fill(ihit->energyLoss());
571  tofTk_[indx]->Fill(ihit->timeOfFlight());
572  nHit++;
573  }
574  hitTk_[indx]->Fill(float(nHit));
575 #ifdef EDM_ML_DEBUG
576  edm::LogVerbatim("HitStudy") << "CaloSimHitStudy::analyzeHits: for " << label << " Index " << indx << " # of Hits "
577  << nHit;
578 #endif
579 }
580 
581 void CaloSimHitStudy::analyzeHits(std::vector<PCaloHit>& ebHits,
582  std::vector<PCaloHit>& eeHits,
583  std::vector<PCaloHit>& hcHits) {
584  double edepEB = 0, edepEBT = 0;
585  for (unsigned int i = 0; i < ebHits.size(); i++) {
586  double edep = ebHits[i].energy();
587  double time = ebHits[i].time();
588  if (((ebHits[i].depth()) & PCaloHit::kEcalDepthIdMask) == 0) {
589  edepEB += edep;
590  if (time < tmax_)
591  edepEBT += edep;
592  }
593  }
594  double edepEE = 0, edepEET = 0;
595  for (unsigned int i = 0; i < eeHits.size(); i++) {
596  double edep = eeHits[i].energy();
597  double time = eeHits[i].time();
598  edepEE += edep;
599  if (time < tmax_)
600  edepEET += edep;
601  }
602  double edepH = 0, edepHT = 0, edepHO = 0, edepHOT = 0;
603  for (unsigned int i = 0; i < hcHits.size(); i++) {
604  double edep = hcHits[i].energy();
605  double time = hcHits[i].time();
606  int subdet(0);
607  if (testNumber_) {
608  int ieta(0), phi(0), z(0), lay(0), depth(0);
609  HcalTestNumbering::unpackHcalIndex(hcHits[i].id(), subdet, z, depth, ieta, phi, lay);
610  } else {
611  HcalDetId id = HcalDetId(hcHits[i].id());
612  subdet = id.subdet();
613  }
614  if (subdet == static_cast<int>(HcalBarrel) || subdet == static_cast<int>(HcalEndcap)) {
615  edepH += edep;
616  if (time < tmax_)
617  edepHT += edep;
618  } else if (subdet == static_cast<int>(HcalOuter)) {
619  edepHO += edep;
620  if (time < tmax_)
621  edepHOT += edep;
622  }
623  }
624  double edepE = edepEB + edepEE;
625  double edepET = edepEBT + edepEET;
626  double edepHC = edepH + edepHO;
627  double edepHCT = edepHT + edepHOT;
628 #ifdef EDM_ML_DEBUG
629  edm::LogVerbatim("HitStudy") << "CaloSimHitStudy::energy in EB " << edepEB << " (" << edepEBT << ") from "
630  << ebHits.size() << " hits; "
631  << " energy in EE " << edepEE << " (" << edepEET << ") from " << eeHits.size()
632  << " hits; energy in HC " << edepH << ", " << edepHO << " (" << edepHT << ", " << edepHOT
633  << ") from " << hcHits.size() << " hits";
634 #endif
635  edepC_[6]->Fill(edepE);
636  edepT_[6]->Fill(edepET);
637  edepC_[7]->Fill(edepH);
638  edepT_[7]->Fill(edepHT);
639  edepC_[8]->Fill(edepHC);
640  edepT_[8]->Fill(edepHCT);
641  if (edepE < eMIP_) {
642  edepC_[0]->Fill(edepE);
643  edepC_[1]->Fill(edepH);
644  edepC_[2]->Fill(edepHC);
645  } else {
646  edepC_[3]->Fill(edepE);
647  edepC_[4]->Fill(edepH);
648  edepC_[5]->Fill(edepHC);
649  }
650  if (edepET < eMIP_) {
651  edepT_[0]->Fill(edepET);
652  edepT_[1]->Fill(edepHT);
653  edepT_[2]->Fill(edepHCT);
654  } else {
655  edepT_[3]->Fill(edepET);
656  edepT_[4]->Fill(edepHT);
657  edepT_[5]->Fill(edepHCT);
658  }
659 }
660 
661 //define this as a plug-in
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
const edm::EDGetTokenT< edm::HepMCProduct > tok_evt_
std::vector< edm::EDGetTokenT< edm::PSimHitContainer > > toks_tkLow_
const std::vector< std::string > tkLowLab_
static constexpr int nCalo_
TH1F * edepT_[nCalo_]
void analyzeHits(std::vector< PCaloHit > &, int)
TH1F * timeAll_[nCalo_]
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
TH1F * hitTk_[nTrack_]
void beginRun(edm::Run const &, edm::EventSetup const &) override
TH1F * edep_[nCalo_]
constexpr HcalSubdetector subdet() const
get the subdetector
Definition: HcalDetId.h:138
char const * label
const double eMIP_
TH1F * etot_[nCalo_]
GlobalPoint getPosition(const DetId &id) const
Get the position of a given detector id.
Definition: CaloGeometry.cc:50
const std::vector< std::string > muonLab_
const std::string g4Label_
Definition: tfile.py:1
const std::vector< edm::EDGetTokenT< edm::PCaloHitContainer > > toks_calo_
TH1F * tofTk_[nTrack_]
const edm::ESGetToken< CaloGeometry, CaloGeometryRecord > tokGeom_
static void unpackHcalIndex(const uint32_t &idx, int &det, int &z, int &depth, int &eta, int &phi, int &lay)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
TH1F * edepHad_[nCalo_]
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
CaloSimHitStudy(const edm::ParameterSet &ps)
const double maxEnergy_
const HepMC::GenEvent * GetEvent() const
Definition: HepMCProduct.h:37
#define M_PI
Definition: DetId.h:17
TH1F * etotg_[nCalo_]
void beginJob() override
TH1F * phi_[nCalo_]
TH1F * time_[nCalo_]
const std::vector< std::string > tkHighLab_
void add(std::string const &label, ParameterSetDescription const &psetDescription)
const bool testNumber_
TH1F * edepTk_[nTrack_]
bool isValid() const
Definition: HandleBase.h:70
~CaloSimHitStudy() override
TH1F * hit_[nCalo_]
void endRun(edm::Run const &, edm::EventSetup const &) override
HLT enums.
std::vector< edm::EDGetTokenT< edm::PSimHitContainer > > toks_tkHigh_
const HcalGeometry * hcalGeom_
std::vector< edm::EDGetTokenT< edm::PSimHitContainer > > toks_track_
TH1F * edepC_[nCalo_]
const bool storeRL_
const double tmax_
void analyze(edm::Event const &, edm::EventSetup const &) override
const std::vector< std::string > hitLab_
TH1F * edepEM_[nCalo_]
TH1F * eta_[nCalo_]
const CaloSubdetectorGeometry * getSubdetectorGeometry(const DetId &id) const
access the subdetector geometry for the given subdetector directly
Definition: CaloGeometry.cc:34
const CaloGeometry * caloGeometry_
static constexpr int nTrack_
*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
Definition: Run.h:45
static const int kEcalDepthIdMask
Definition: PCaloHit.h:60