CMS 3D CMS Logo

HCalSD.cc
Go to the documentation of this file.
1 // File: HCalSD.cc
3 // Description: Sensitive Detector class for calorimeters
5 
14 
19 
20 #include "G4LogicalVolumeStore.hh"
21 #include "G4LogicalVolume.hh"
22 #include "G4Step.hh"
23 #include "G4Track.hh"
24 
25 #include "G4SystemOfUnits.hh"
26 #include "G4PhysicalConstants.hh"
27 #include "Randomize.hh"
28 
29 #include <iostream>
30 #include <fstream>
31 #include <iomanip>
32 #include <sstream>
33 
34 //#define EDM_ML_DEBUG
35 //#define plotDebug
36 //#define printDebug
37 
38 #ifdef plotDebug
39 #include <TH1F.h>
40 #endif
41 
43  const HcalDDDSimConstants* hcns,
44  const HcalDDDRecConstants* hcnr,
45  const HcalSimulationConstants* hscs,
46  const HBHEDarkening* hbd,
47  const HBHEDarkening* hed,
48  const SensitiveDetectorCatalog& clg,
49  edm::ParameterSet const& p,
50  const SimTrackManager* manager)
51  : CaloSD(name,
52  clg,
53  p,
54  manager,
55  (float)(p.getParameter<edm::ParameterSet>("HCalSD").getParameter<double>("TimeSliceUnit")),
56  p.getParameter<edm::ParameterSet>("HCalSD").getParameter<bool>("IgnoreTrackID")),
57  hcalConstants_(hcns),
58  hcalSimConstants_(hscs),
59  m_HBDarkening(hbd),
60  m_HEDarkening(hed),
61  isHF(false),
62  weight_(1.0),
63  depth_(1) {
64  numberingFromDDD.reset(nullptr);
65  numberingScheme.reset(nullptr);
66  showerLibrary.reset(nullptr);
67  hfshower.reset(nullptr);
68  showerParam.reset(nullptr);
69  showerPMT.reset(nullptr);
70  showerBundle.reset(nullptr);
71  m_HFDarkening.reset(nullptr);
72  m_HcalTestNS.reset(nullptr);
73 
74  //static SimpleConfigurable<double> bk1(0.013, "HCalSD:BirkC1");
75  //static SimpleConfigurable<double> bk2(0.0568,"HCalSD:BirkC2");
76  //static SimpleConfigurable<double> bk3(1.75, "HCalSD:BirkC3");
77  // Values from NIM 80 (1970) 239-244: as implemented in Geant3
78 
79  dd4hep_ = p.getParameter<bool>("g4GeometryDD4hepSource");
80  edm::ParameterSet m_HC = p.getParameter<edm::ParameterSet>("HCalSD");
81  useBirk = m_HC.getParameter<bool>("UseBirkLaw");
82  double bunit = (CLHEP::g / (CLHEP::MeV * CLHEP::cm2));
83  birk1 = m_HC.getParameter<double>("BirkC1") * bunit;
84  birk2 = m_HC.getParameter<double>("BirkC2");
85  birk3 = m_HC.getParameter<double>("BirkC3");
86  useShowerLibrary = m_HC.getParameter<bool>("UseShowerLibrary");
87  useParam = m_HC.getParameter<bool>("UseParametrize");
88  testNumber = m_HC.getParameter<bool>("TestNumberingScheme");
89  neutralDensity = m_HC.getParameter<bool>("doNeutralDensityFilter");
90  usePMTHit = m_HC.getParameter<bool>("UsePMTHits");
91  betaThr = m_HC.getParameter<double>("BetaThreshold");
92  eminHitHB = m_HC.getParameter<double>("EminHitHB") * MeV;
93  eminHitHE = m_HC.getParameter<double>("EminHitHE") * MeV;
94  eminHitHO = m_HC.getParameter<double>("EminHitHO") * MeV;
95  eminHitHF = m_HC.getParameter<double>("EminHitHF") * MeV;
96  useFibreBundle = m_HC.getParameter<bool>("UseFibreBundleHits");
97  deliveredLumi = m_HC.getParameter<double>("DelivLuminosity");
98  agingFlagHB = m_HC.getParameter<bool>("HBDarkening");
99  agingFlagHE = m_HC.getParameter<bool>("HEDarkening");
100  bool agingFlagHF = m_HC.getParameter<bool>("HFDarkening");
101  useHF = m_HC.getUntrackedParameter<bool>("UseHF", true);
102  bool forTBHC = m_HC.getUntrackedParameter<bool>("ForTBHCAL", false);
103  bool forTBH2 = m_HC.getUntrackedParameter<bool>("ForTBH2", false);
104  useLayerWt = m_HC.getUntrackedParameter<bool>("UseLayerWt", false);
105  std::string file = m_HC.getUntrackedParameter<std::string>("WtFile", "None");
106  testNS_ = m_HC.getUntrackedParameter<bool>("TestNS", false);
107  edm::ParameterSet m_HF = p.getParameter<edm::ParameterSet>("HFShower");
108  applyFidCut = m_HF.getParameter<bool>("ApplyFiducialCut");
109  bool dumpGeom = m_HC.getUntrackedParameter<bool>("DumpGeometry", false);
110 
111 #ifdef EDM_ML_DEBUG
112  edm::LogVerbatim("HcalSim") << "***************************************************"
113  << "\n"
114  << "* Constructing a HCalSD with name " << name << "\n"
115  << "\n"
116  << "***************************************************";
117 #endif
118  edm::LogVerbatim("HcalSim") << "HCalSD:: Use of HF code is set to " << useHF
119  << "\n Use of shower parametrization set to " << useParam
120  << "\n Use of shower library is set to " << useShowerLibrary
121  << "\n Use PMT Hit is set to " << usePMTHit << " with beta Threshold " << betaThr
122  << "\n USe of FibreBundle Hit set to " << useFibreBundle
123  << "\n Use of Birks law is set to " << useBirk
124  << " with three constants kB = " << birk1 / bunit << ", C1 = " << birk2
125  << ", C2 = " << birk3;
126  edm::LogVerbatim("HcalSim") << "HCalSD:: Suppression Flag " << suppressHeavy << " protons below " << kmaxProton
127  << " MeV,"
128  << " neutrons below " << kmaxNeutron << " MeV and"
129  << " ions below " << kmaxIon << " MeV\n"
130  << " Threshold for storing hits in HB: " << eminHitHB << " HE: " << eminHitHE
131  << " HO: " << eminHitHO << " HF: " << eminHitHF << "\n"
132  << "Delivered luminosity for Darkening " << deliveredLumi << " Flag (HE) " << agingFlagHE
133  << " Flag (HB) " << agingFlagHB << " Flag (HF) " << agingFlagHF << "\n"
134  << "Application of Fiducial Cut " << applyFidCut
135  << "Flag for test number|neutral density filter " << testNumber << " " << neutralDensity;
136 
137  if (forTBHC) {
138  useHF = false;
139  matNames.emplace_back("Scintillator");
140  } else {
142  }
143 
145  if (testNumber || forTBH2) {
146  scheme = dynamic_cast<HcalNumberingScheme*>(new HcalTestNumberingScheme(forTBH2));
147  } else {
148  scheme = new HcalNumberingScheme();
149  }
151 
152  // always call getFromLibrary() method to identify HF region
153  setParameterized(true);
154 
155  const G4LogicalVolumeStore* lvs = G4LogicalVolumeStore::GetInstance();
156  // std::vector<G4LogicalVolume *>::const_iterator lvcite;
157  const G4LogicalVolume* lv;
158  std::string attribute, value;
159 
160  if (useHF) {
161  if (useParam) {
162  showerParam = std::make_unique<HFShowerParam>(name, hcalConstants_, hcalSimConstants_->hcalsimpar(), p);
163  } else {
164  if (useShowerLibrary) {
165  showerLibrary = std::make_unique<HFShowerLibrary>(name, hcalConstants_, hcalSimConstants_->hcalsimpar(), p);
166  }
167  hfshower = std::make_unique<HFShower>(name, hcalConstants_, hcalSimConstants_->hcalsimpar(), p, 0);
168  }
169 
170  // HF volume names
172  const std::vector<int>& temp = hcalSimConstants_->hcalsimpar()->hfLevels_;
173 #ifdef EDM_ML_DEBUG
174  std::stringstream ss0;
175  ss0 << "HCalSD: Names to be tested for Volume = HF has " << hfNames.size() << " elements";
176 #endif
177  int addlevel = dd4hep_ ? 1 : 0;
178  for (unsigned int i = 0; i < hfNames.size(); ++i) {
180  lv = nullptr;
181  for (auto lvol : *lvs) {
182  if (DD4hep2DDDName::nameMatterLV(lvol->GetName(), dd4hep_) == namv) {
183  lv = lvol;
184  break;
185  }
186  }
187  hfLV.emplace_back(lv);
188  hfLevels.emplace_back(temp[i] + addlevel);
189 #ifdef EDM_ML_DEBUG
190  ss0 << "\n HF[" << i << "] = " << namv << " LV " << lv << " at level " << (temp[i] + addlevel);
191 #endif
192  }
193 #ifdef EDM_ML_DEBUG
194  edm::LogVerbatim("HcalSim") << ss0.str();
195 #endif
196  // HF Fibre volume names
199  const std::vector<std::string>& pmtNames = hcalSimConstants_->hcalsimpar()->hfPMTNames_;
200  fillLogVolumeVector("HFPMT", pmtNames, pmtLV);
201  const std::vector<std::string>& straightNames = hcalSimConstants_->hcalsimpar()->hfFibreStraightNames_;
202  fillLogVolumeVector("HFFibreBundleStraight", straightNames, fibre1LV);
203  const std::vector<std::string>& conicalNames = hcalSimConstants_->hcalsimpar()->hfFibreConicalNames_;
204  fillLogVolumeVector("HFFibreBundleConical", conicalNames, fibre2LV);
205  }
206 
207  //Material list for HB/HE/HO sensitive detectors
208  const G4MaterialTable* matTab = G4Material::GetMaterialTable();
209  std::vector<G4Material*>::const_iterator matite;
210  for (auto const& namx : matNames) {
211  const G4Material* mat = nullptr;
212  for (matite = matTab->begin(); matite != matTab->end(); ++matite) {
213  if (DD4hep2DDDName::nameMatterLV((*matite)->GetName(), dd4hep_) == namx) {
214  mat = (*matite);
215  break;
216  }
217  }
218  materials.emplace_back(mat);
219  }
220 #ifdef EDM_ML_DEBUG
221  std::stringstream ss1;
222  for (unsigned int i = 0; i < matNames.size(); ++i) {
223  if (i / 10 * 10 == i) {
224  ss1 << "\n";
225  }
226  ss1 << " " << matNames[i];
227  }
228  edm::LogVerbatim("HcalSim") << "HCalSD: Material names for HCAL: " << ss1.str();
229 #endif
230  if (useLayerWt) {
232  }
233  numberingFromDDD = std::make_unique<HcalNumberingFromDDD>(hcalConstants_);
234 
235  //Special Geometry parameters
237 #ifdef EDM_ML_DEBUG
238  std::stringstream ss2;
239  for (unsigned int ig = 0; ig < gpar.size(); ig++) {
240  ss2 << "\n gpar[" << ig << "] = " << gpar[ig] / cm << " cm";
241  }
242  edm::LogVerbatim("HcalSim") << "Maximum depth for HF " << hcalConstants_->getMaxDepth(2) << gpar.size()
243  << " gpar (cm)" << ss2.str();
244 #endif
245 
246  //Test Hcal Numbering Scheme
247  if (testNS_)
248  m_HcalTestNS = std::make_unique<HcalTestNS>(hcnr);
249 
250  for (int i = 0; i < 9; ++i) {
251  hit_[i] = time_[i] = dist_[i] = nullptr;
252  }
253  hzvem = hzvhad = nullptr;
254 
255  if (agingFlagHF) {
256  m_HFDarkening = std::make_unique<HFDarkening>(m_HC.getParameter<edm::ParameterSet>("HFDarkeningParameterBlock"));
257  }
258 #ifdef plotDebug
260 
261  if (tfile.isAvailable()) {
262  static const char* const labels[] = {"HB",
263  "HE",
264  "HO",
265  "HF Absorber",
266  "HF PMT",
267  "HF Absorber Long",
268  "HF Absorber Short",
269  "HF PMT Long",
270  "HF PMT Short"};
271  TFileDirectory hcDir = tfile->mkdir("ProfileFromHCalSD");
272  char name[20], title[60];
273  for (int i = 0; i < 9; ++i) {
274  sprintf(title, "Hit energy in %s", labels[i]);
275  sprintf(name, "HCalSDHit%d", i);
276  hit_[i] = hcDir.make<TH1F>(name, title, 2000, 0., 2000.);
277  sprintf(title, "Energy (MeV)");
278  hit_[i]->GetXaxis()->SetTitle(title);
279  hit_[i]->GetYaxis()->SetTitle("Hits");
280  sprintf(title, "Time of the hit in %s", labels[i]);
281  sprintf(name, "HCalSDTime%d", i);
282  time_[i] = hcDir.make<TH1F>(name, title, 2000, 0., 2000.);
283  sprintf(title, "Time (ns)");
284  time_[i]->GetXaxis()->SetTitle(title);
285  time_[i]->GetYaxis()->SetTitle("Hits");
286  sprintf(title, "Longitudinal profile in %s", labels[i]);
287  sprintf(name, "HCalSDDist%d", i);
288  dist_[i] = hcDir.make<TH1F>(name, title, 2000, 0., 2000.);
289  sprintf(title, "Distance (mm)");
290  dist_[i]->GetXaxis()->SetTitle(title);
291  dist_[i]->GetYaxis()->SetTitle("Hits");
292  }
293  if (useHF && (!useParam)) {
294  hzvem = hcDir.make<TH1F>("hzvem", "Longitudinal Profile (EM Part)", 330, 0.0, 1650.0);
295  hzvem->GetXaxis()->SetTitle("Longitudinal Profile (EM Part)");
296  hzvhad = hcDir.make<TH1F>("hzvhad", "Longitudinal Profile (Had Part)", 330, 0.0, 1650.0);
297  hzvhad->GetXaxis()->SetTitle("Longitudinal Profile (Hadronic Part)");
298  }
299  }
300 #endif
301  if (dumpGeom) {
302  std::unique_ptr<HcalNumberingFromDDD> hcn = std::make_unique<HcalNumberingFromDDD>(hcalConstants_);
303  const auto& lvNames = clg.logicalNames(name);
304  HcalDumpGeometry geom(lvNames, hcn.get(), testNumber, false);
305  geom.update();
306  }
307 }
308 
310  const std::vector<std::string>& lvnames,
311  std::vector<const G4LogicalVolume*>& lvvec) {
312  const G4LogicalVolumeStore* lvs = G4LogicalVolumeStore::GetInstance();
313  const G4LogicalVolume* lv;
314  std::stringstream ss3;
315  ss3 << "HCalSD: " << lvnames.size() << " names to be tested for Volume <" << value << ">:";
316  for (unsigned int i = 0; i < lvnames.size(); ++i) {
318  lv = nullptr;
319  for (auto lvol : *lvs) {
320  if (DD4hep2DDDName::nameMatterLV(lvol->GetName(), dd4hep_) == namv) {
321  lv = lvol;
322  break;
323  }
324  }
325  lvvec.emplace_back(lv);
326  if (i / 10 * 10 == i) {
327  ss3 << "\n";
328  }
329  ss3 << " " << namv;
330  }
331  edm::LogVerbatim("HcalSim") << ss3.str();
332 }
333 
334 bool HCalSD::getFromLibrary(const G4Step* aStep) {
335  auto const track = aStep->GetTrack();
336  depth_ = (aStep->GetPreStepPoint()->GetTouchable()->GetReplicaNumber(0)) % 10;
337  weight_ = 1.0;
338  bool kill(false);
339  isHF = isItHF(aStep);
340 #ifdef EDM_ML_DEBUG
341  edm::LogVerbatim("HcalSim") << "GetFromLibrary: "
342  << (aStep->GetPreStepPoint()->GetTouchable()->GetVolume(0)->GetName())
343 
344  << " isHF " << isHF << " darken " << (m_HFDarkening != nullptr) << " useParam "
345  << useParam << " useShowerLibrary " << useShowerLibrary << " Muon? "
346  << G4TrackToParticleID::isMuon(track) << " electron? "
347  << G4TrackToParticleID::isGammaElectronPositron(track) << " Stable Hadron? "
349 #endif
350  if (isHF) {
351  if (m_HFDarkening) {
352  G4ThreeVector hitPoint = aStep->GetPreStepPoint()->GetPosition();
353  const double invcm = 1. / CLHEP::cm;
354  double r = hitPoint.perp() * invcm;
355  double z = std::abs(hitPoint.z()) * invcm;
356  double dose_acquired = 0.;
358  unsigned int hfZLayer = (unsigned int)((z - HFDarkening::lowZLimit) / 5);
359  if (hfZLayer >= HFDarkening::upperZLimit)
360  hfZLayer = (HFDarkening::upperZLimit - 1);
361  float normalized_lumi = m_HFDarkening->int_lumi(deliveredLumi);
362  for (int i = hfZLayer; i != HFDarkening::numberOfZLayers; ++i) {
363  dose_acquired = m_HFDarkening->dose(i, r);
364  weight_ *= m_HFDarkening->degradation(normalized_lumi * dose_acquired);
365  }
366  }
367 #ifdef EDM_ML_DEBUG
368  edm::LogVerbatim("HcalSim") << "HCalSD::getFromLibrary: HFLumiDarkening at "
369  << "r= " << r << ", z= " << z << " Dose= " << dose_acquired << " weight= " << weight_;
370 #endif
371  }
372 
373  if (useParam) {
374  getFromParam(aStep, kill);
375 #ifdef EDM_ML_DEBUG
376  G4String nameVolume = aStep->GetPreStepPoint()->GetPhysicalVolume()->GetLogicalVolume()->GetName();
377  edm::LogVerbatim("HcalSim") << "HCalSD: " << getNumberOfHits() << " hits from parametrization in " << nameVolume
378  << " for Track " << track->GetTrackID() << " ("
379  << track->GetDefinition()->GetParticleName() << ")";
380 #endif
383 #ifdef EDM_ML_DEBUG
384  auto nameVolume = aStep->GetPreStepPoint()->GetPhysicalVolume()->GetLogicalVolume()->GetName();
385  edm::LogVerbatim("HcalSim") << "HCalSD: Starts shower library from " << nameVolume << " for Track "
386  << track->GetTrackID() << " (" << track->GetDefinition()->GetParticleName() << ")";
387 
388 #endif
389  getFromHFLibrary(aStep, kill);
390  }
391  }
392  }
393 #ifdef EDM_ML_DEBUG
394  edm::LogVerbatim("HcalSim") << "HCalSD::getFromLibrary ID= " << track->GetTrackID() << " ("
395  << track->GetDefinition()->GetParticleName() << ") kill= " << kill
396  << " weight= " << weight_ << " depth= " << depth_ << " isHF: " << isHF;
397 #endif
398  return kill;
399 }
400 
401 double HCalSD::getEnergyDeposit(const G4Step* aStep) {
402  double destep(0.0);
403  auto const lv = aStep->GetPreStepPoint()->GetPhysicalVolume()->GetLogicalVolume();
404  auto const theTrack = aStep->GetTrack();
405 
406  if (isHF) {
407  if (useShowerLibrary && G4TrackToParticleID::isMuon(theTrack) && isItFibre(lv)) {
408 #ifdef EDM_ML_DEBUG
409  edm::LogVerbatim("HcalSim") << "HCalSD: Hit at Fibre in LV " << lv->GetName() << " for track "
410  << aStep->GetTrack()->GetTrackID() << " ("
411  << aStep->GetTrack()->GetDefinition()->GetParticleName() << ")";
412 #endif
413  hitForFibre(aStep);
414  }
415  return destep;
416  }
417 
418  if (isItPMT(lv)) {
419  if (usePMTHit && showerPMT) {
420  getHitPMT(aStep);
421  }
422 #ifdef EDM_ML_DEBUG
423  edm::LogVerbatim("HcalSim") << "HCalSD: Hit from PMT parametrization in LV " << lv->GetName() << " for Track "
424  << aStep->GetTrack()->GetTrackID() << " ("
425  << aStep->GetTrack()->GetDefinition()->GetParticleName() << ")";
426 #endif
427  return destep;
428 
429  } else if (isItStraightBundle(lv)) {
430  if (useFibreBundle && showerBundle) {
431  getHitFibreBundle(aStep, false);
432  }
433 #ifdef EDM_ML_DEBUG
434  edm::LogVerbatim("HcalSim") << "HCalSD: Hit from straight FibreBundle in LV: " << lv->GetName() << " for track "
435  << aStep->GetTrack()->GetTrackID() << " ("
436  << aStep->GetTrack()->GetDefinition()->GetParticleName() << ")";
437 #endif
438  return destep;
439 
440  } else if (isItConicalBundle(lv)) {
441  if (useFibreBundle && showerBundle) {
442  getHitFibreBundle(aStep, true);
443  }
444 #ifdef EDM_ML_DEBUG
445  edm::LogVerbatim("HcalSim") << "HCalSD: Hit from conical FibreBundle PV: " << lv->GetName() << " for track "
446  << aStep->GetTrack()->GetTrackID() << " ("
447  << aStep->GetTrack()->GetDefinition()->GetParticleName() << ")";
448 #endif
449  return destep;
450  }
451 
452  // normal hit
453  destep = aStep->GetTotalEnergyDeposit();
454 
455  const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
456  uint32_t detid = setDetUnitId(aStep);
457  int det(0), ieta(0), phi(0), z(0), lay, depth(-1);
458  if (testNumber) {
459  HcalTestNumbering::unpackHcalIndex(detid, det, z, depth, ieta, phi, lay);
460  if (z == 0) {
461  z = -1;
462  }
463  } else {
464  HcalDetId hcid(detid);
465  det = hcid.subdetId();
466  ieta = hcid.ietaAbs();
467  phi = hcid.iphi();
468  z = hcid.zside();
469  }
470  lay = (touch->GetReplicaNumber(0) / 10) % 100 + 1;
471 #ifdef EDM_ML_DEBUG
472  edm::LogVerbatim("HcalSim") << "HCalSD: det: " << det << " ieta: " << ieta << " iphi: " << phi << " zside " << z
473  << " lay: " << lay - 2;
474 #endif
475  if (depth_ == 0 && (det == 1 || det == 2) && ((!testNumber) || neutralDensity))
477  if (useLayerWt) {
478  G4ThreeVector hitPoint = aStep->GetPreStepPoint()->GetPosition();
479  weight_ = layerWeight(det + 2, hitPoint, depth_, lay);
480  }
481 
482  if (agingFlagHB && m_HBDarkening && det == 1) {
483  double dweight = m_HBDarkening->degradation(deliveredLumi, ieta, lay);
484  weight_ *= dweight;
485 #ifdef EDM_ML_DEBUG
486  edm::LogVerbatim("HcalSim") << "HCalSD: HB Lumi: " << deliveredLumi << " coefficient = " << dweight
487  << " Weight= " << weight_;
488 #endif
489  }
490 
491  if (agingFlagHE && m_HEDarkening && det == 2) {
492  double dweight = m_HEDarkening->degradation(deliveredLumi, ieta, lay);
493  weight_ *= dweight;
494 #ifdef EDM_ML_DEBUG
495  edm::LogVerbatim("HcalSim") << "HCalSD: HB Lumi: " << deliveredLumi << " coefficient = " << dweight
496  << " Weight= " << weight_;
497 #endif
498  }
499 
500  if (suppressHeavy) {
501  TrackInformation* trkInfo = (TrackInformation*)(theTrack->GetUserInformation());
502  if (trkInfo) {
503  int pdg = theTrack->GetDefinition()->GetPDGEncoding();
504  if (!(trkInfo->isPrimary())) { // Only secondary particles
505  double ke = theTrack->GetKineticEnergy();
506  if (pdg / 1000000000 == 1 && (pdg / 10000) % 100 > 0 && (pdg / 10) % 100 > 0 && ke < kmaxIon)
507  weight_ = 0;
508  if ((pdg == 2212) && (ke < kmaxProton))
509  weight_ = 0;
510  if ((pdg == 2112) && (ke < kmaxNeutron))
511  weight_ = 0;
512  }
513  }
514  }
515  double wt0(1.0);
516  if (useBirk) {
517  const G4Material* mat = aStep->GetPreStepPoint()->GetMaterial();
518  if (isItScintillator(mat))
519  wt0 = getAttenuation(aStep, birk1, birk2, birk3);
520  }
521  weight_ *= wt0;
522  double wt1 = getResponseWt(theTrack);
523  double wt2 = theTrack->GetWeight();
524  double edep = weight_ * wt1 * destep;
525  if (wt2 > 0.0) {
526  edep *= wt2;
527  }
528 #ifdef EDM_ML_DEBUG
529  edm::LogVerbatim("HcalSim") << "HCalSD: edep= " << edep << " Det: " << det + 2 << " depth= " << depth_
530  << " weight= " << weight_ << " wt0= " << wt0 << " wt1= " << wt1 << " wt2= " << wt2;
531 #endif
532  return edep;
533 }
534 
535 uint32_t HCalSD::setDetUnitId(const G4Step* aStep) {
536  auto const prePoint = aStep->GetPreStepPoint();
537  auto const touch = prePoint->GetTouchable();
538  const G4ThreeVector& hitPoint = prePoint->GetPosition();
539 
540  int depth = (touch->GetReplicaNumber(0)) % 10 + 1;
541  int lay = (touch->GetReplicaNumber(0) / 10) % 100 + 1;
542  int det = (touch->GetReplicaNumber(1)) / 1000;
543 
544  uint32_t idx = setDetUnitId(det, hitPoint, depth, lay);
545 #ifdef EDM_ML_DEBUG
546  if (depth == 1) {
547  edm::LogVerbatim("HcalSim") << "HCalSD: Check for " << det << ":" << depth << ":" << lay << " ID " << std::hex
548  << idx << std::dec;
549  int det0, z0, depth0, eta0, phi0, lay0(-1);
550  if (testNumber) {
551  HcalTestNumbering::unpackHcalIndex(idx, det0, z0, depth0, eta0, phi0, lay0);
552  } else {
553  HcalDetId hcid0(idx);
554  det0 = hcid0.subdetId();
555  eta0 = hcid0.ietaAbs();
556  phi0 = hcid0.iphi();
557  z0 = hcid0.zside();
558  depth0 = hcid0.depth();
559  }
560  edm::LogVerbatim("HcalSim") << "HCalSD: det|z|depth|eta|phi|lay " << det0 << ":" << z0 << ":" << depth0 << ":"
561  << eta0 << ":" << phi0 << ":" << lay0;
562  printVolume(touch);
563  }
564 #endif
565  return idx;
566 }
567 
569  if (scheme != nullptr) {
570  edm::LogVerbatim("HcalSim") << "HCalSD: updates numbering scheme for " << GetName();
571  numberingScheme.reset(scheme);
572  }
573 }
574 
575 void HCalSD::update(const BeginOfJob* job) {}
576 
578 
579 bool HCalSD::filterHit(CaloG4Hit* aHit, double time) {
580  double threshold = 0;
581  DetId theId(aHit->getUnitID());
582  switch (theId.subdetId()) {
583  case HcalBarrel:
585  break;
586  case HcalEndcap:
588  break;
589  case HcalOuter:
591  break;
592  case HcalForward:
594  break;
595  default:
596  break;
597  }
598  return ((time <= tmaxHit) && (aHit->getEnergyDeposit() > threshold));
599 }
600 
601 uint32_t HCalSD::setDetUnitId(int det, const G4ThreeVector& pos, int depth, int lay = 1) {
602  uint32_t id = 0;
603  if (det == 0) {
604 #ifdef printDebug
605  double eta = std::abs(pos.eta());
606 #endif
607  if (std::abs(pos.z()) > maxZ_) {
608  det = 5;
609 #ifdef printDebug
610  if (eta < 2.868)
611  ++detNull_[2];
612 #endif
613  } else if (!(hcalConstants_->isHE())) {
614  det = 3;
615 #ifdef printDebug
616  ++detNull_[0];
617 #endif
618  } else {
619  double minR = minRoff_ + slopeHE_ * std::abs(pos.z());
620  double maxR = maxRoff_ + slopeHE_ * std::abs(pos.z());
621  det = ((pos.perp() > minR) && (pos.perp() < maxR)) ? 4 : 3;
622 #ifdef printDebug
623  ++detNull_[det - 3];
624 #endif
625  }
626 #ifdef printDEBUG
627  edm::LogVerbatim("HcalSim") << "Position " << pos.perp() << ":" << std::abs(pos.z()) << " Limits "
628  << !(hcalConstants_->isHE()) << ":" << maxZ_ << " det " << det;
629  } else {
630  ++detNull_[3];
631 #endif
632  }
633 
634  if (numberingFromDDD.get()) {
635  //get the ID's as eta, phi, depth, ... indices
637  numberingFromDDD->unitID(det, math::XYZVectorD(pos.x(), pos.y(), pos.z()), depth, lay);
638  id = setDetUnitId(tmp);
639  }
640  return id;
641 }
642 
644  modifyDepth(tmp);
645  uint32_t id = (numberingScheme.get()) ? numberingScheme->getUnitID(tmp) : 0;
646  if ((!testNumber) && m_HcalTestNS.get()) {
647  bool ok = m_HcalTestNS->compare(tmp, id);
648  if (!ok)
649  edm::LogWarning("HcalSim") << "Det ID from HCalSD " << HcalDetId(id) << " " << std::hex << id << std::dec
650  << " does not match one from relabller for " << tmp.subdet << ":" << tmp.etaR << ":"
651  << tmp.phi << ":" << tmp.phis << ":" << tmp.depth << ":" << tmp.lay << std::endl;
652  }
653  return id;
654 }
655 
656 bool HCalSD::isItHF(const G4Step* aStep) {
657  const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
658  int levels = (touch->GetHistoryDepth()) + 1;
659  for (unsigned int it = 0; it < hfNames.size(); ++it) {
660  if (levels >= hfLevels[it]) {
661  const G4LogicalVolume* lv = touch->GetVolume(levels - hfLevels[it])->GetLogicalVolume();
662  if (lv == hfLV[it])
663  return true;
664  }
665  }
666  return false;
667 }
668 
669 bool HCalSD::isItHF(const G4String& name) {
670  for (const auto& nam : hfNames)
671  if (name == static_cast<G4String>(nam)) {
672  return true;
673  }
674  return false;
675 }
676 
677 bool HCalSD::isItFibre(const G4LogicalVolume* lv) {
678  for (auto lvol : fibreLV)
679  if (lv == lvol) {
680  return true;
681  }
682  return false;
683 }
684 
685 bool HCalSD::isItFibre(const G4String& name) {
686  for (const auto& nam : fibreNames)
687  if (name == static_cast<G4String>(nam)) {
688  return true;
689  }
690  return false;
691 }
692 
693 bool HCalSD::isItPMT(const G4LogicalVolume* lv) {
694  for (auto lvol : pmtLV)
695  if (lv == lvol) {
696  return true;
697  }
698  return false;
699 }
700 
701 bool HCalSD::isItStraightBundle(const G4LogicalVolume* lv) {
702  for (auto lvol : fibre1LV)
703  if (lv == lvol) {
704  return true;
705  }
706  return false;
707 }
708 
709 bool HCalSD::isItConicalBundle(const G4LogicalVolume* lv) {
710  for (auto lvol : fibre2LV)
711  if (lv == lvol) {
712  return true;
713  }
714  return false;
715 }
716 
717 bool HCalSD::isItScintillator(const G4Material* mat) {
718  for (auto amat : materials)
719  if (amat == mat) {
720  return true;
721  }
722  return false;
723 }
724 
725 bool HCalSD::isItinFidVolume(const G4ThreeVector& hitPoint) {
726  bool flag = true;
727  if (applyFidCut) {
728  int npmt = HFFibreFiducial::PMTNumber(hitPoint);
729 #ifdef EDM_ML_DEBUG
730  edm::LogVerbatim("HcalSim") << "HCalSD::isItinFidVolume:#PMT= " << npmt << " for hit point " << hitPoint;
731 #endif
732  if (npmt <= 0)
733  flag = false;
734  }
735 #ifdef EDM_ML_DEBUG
736  edm::LogVerbatim("HcalSim") << "HCalSD::isItinFidVolume: point " << hitPoint << " return flag " << flag;
737 #endif
738  return flag;
739 }
740 
741 void HCalSD::getFromHFLibrary(const G4Step* aStep, bool& isKilled) {
742  std::vector<HFShowerLibrary::Hit> hits = showerLibrary->getHits(aStep, isKilled, weight_, false);
743  if (!isKilled || hits.empty()) {
744  return;
745  }
746 
747  int primaryID = setTrackID(aStep);
748 
749  // Reset entry point for new primary
750  resetForNewPrimary(aStep);
751 
752  auto const theTrack = aStep->GetTrack();
753  int det = 5;
754 
756  edepositEM = 1. * GeV;
757  edepositHAD = 0.;
758  } else {
759  edepositEM = 0.;
760  edepositHAD = 1. * GeV;
761  }
762 #ifdef EDM_ML_DEBUG
763  edm::LogVerbatim("HcalSim") << "HCalSD::getFromLibrary " << hits.size() << " hits for " << GetName() << " of "
764  << primaryID << " with " << theTrack->GetDefinition()->GetParticleName() << " of "
765  << aStep->GetPreStepPoint()->GetKineticEnergy() / GeV << " GeV";
766 #endif
767  for (unsigned int i = 0; i < hits.size(); ++i) {
768  G4ThreeVector hitPoint = hits[i].position;
769  if (isItinFidVolume(hitPoint)) {
770  int depth = hits[i].depth;
771  double time = hits[i].time;
772  unsigned int unitID = setDetUnitId(det, hitPoint, depth);
773  currentID[0].setID(unitID, time, primaryID, 0);
774 #ifdef plotDebug
775  plotProfile(aStep, hitPoint, 1.0 * GeV, time, depth);
776  bool emType = G4TrackToParticleID::isGammaElectronPositron(theTrack->GetDefinition()->GetPDGEncoding());
777  plotHF(hitPoint, emType);
778 #endif
779  processHit(aStep);
780  }
781  }
782 }
783 
784 void HCalSD::hitForFibre(const G4Step* aStep) { // if not ParamShower
785 
786  std::vector<HFShower::Hit> hits = hfshower->getHits(aStep, weight_);
787  if (hits.empty()) {
788  return;
789  }
790 
791  auto const theTrack = aStep->GetTrack();
792  int primaryID = setTrackID(aStep);
793  int det = 5;
794 
796  edepositEM = 1. * GeV;
797  edepositHAD = 0.;
798  } else {
799  edepositEM = 0.;
800  edepositHAD = 1. * GeV;
801  }
802 
803 #ifdef EDM_ML_DEBUG
804  edm::LogVerbatim("HcalSim") << "HCalSD::hitForFibre " << hits.size() << " hits for " << GetName() << " of "
805  << primaryID << " with " << theTrack->GetDefinition()->GetParticleName() << " of "
806  << aStep->GetPreStepPoint()->GetKineticEnergy() / GeV << " GeV in detector type " << det;
807 #endif
808 
809  for (unsigned int i = 0; i < hits.size(); ++i) {
810  G4ThreeVector hitPoint = hits[i].position;
811  if (isItinFidVolume(hitPoint)) {
812  int depth = hits[i].depth;
813  double time = hits[i].time;
814  unsigned int unitID = setDetUnitId(det, hitPoint, depth);
815  currentID[0].setID(unitID, time, primaryID, 0);
816 #ifdef plotDebug
817  plotProfile(aStep, hitPoint, edepositEM, time, depth);
818  bool emType = (edepositEM > 0.) ? true : false;
819  plotHF(hitPoint, emType);
820 #endif
821  processHit(aStep);
822  }
823  }
824 }
825 
826 void HCalSD::getFromParam(const G4Step* aStep, bool& isKilled) {
827  std::vector<HFShowerParam::Hit> hits = showerParam->getHits(aStep, weight_, isKilled);
828  if (!isKilled || hits.empty()) {
829  return;
830  }
831 
832  int primaryID = setTrackID(aStep);
833  int det = 5;
834 
835 #ifdef EDM_ML_DEBUG
836  edm::LogVerbatim("HcalSim") << "HCalSD::getFromParam " << hits.size() << " hits for " << GetName() << " of "
837  << primaryID << " with " << aStep->GetTrack()->GetDefinition()->GetParticleName()
838  << " of " << aStep->GetPreStepPoint()->GetKineticEnergy() / GeV
839  << " GeV in detector type " << det;
840 #endif
841  for (unsigned int i = 0; i < hits.size(); ++i) {
842  G4ThreeVector hitPoint = hits[i].position;
843  int depth = hits[i].depth;
844  double time = hits[i].time;
845  unsigned int unitID = setDetUnitId(det, hitPoint, depth);
846  currentID[0].setID(unitID, time, primaryID, 0);
847  edepositEM = hits[i].edep * GeV;
848  edepositHAD = 0.;
849 #ifdef plotDebug
850  plotProfile(aStep, hitPoint, edepositEM, time, depth);
851 #endif
852  processHit(aStep);
853  }
854 }
855 
856 void HCalSD::getHitPMT(const G4Step* aStep) {
857  auto const preStepPoint = aStep->GetPreStepPoint();
858  auto const theTrack = aStep->GetTrack();
859  double edep = showerPMT->getHits(aStep);
860 
861  if (edep >= 0.) {
862  edep *= GeV;
863  double etrack = preStepPoint->GetKineticEnergy();
864  int primaryID = 0;
865  if (etrack >= energyCut) {
866  primaryID = theTrack->GetTrackID();
867  } else {
868  primaryID = theTrack->GetParentID();
869  if (primaryID == 0)
870  primaryID = theTrack->GetTrackID();
871  }
872  // Reset entry point for new primary
873  resetForNewPrimary(aStep);
874  //
875  int det = static_cast<int>(HcalForward);
876  const G4ThreeVector& hitPoint = preStepPoint->GetPosition();
877  double rr = (hitPoint.x() * hitPoint.x() + hitPoint.y() * hitPoint.y());
878  double phi = (rr == 0. ? 0. : atan2(hitPoint.y(), hitPoint.x()));
879  double etaR = showerPMT->getRadius();
880  int depth = 1;
881  if (etaR < 0) {
882  depth = 2;
883  etaR = -etaR;
884  }
885  if (hitPoint.z() < 0)
886  etaR = -etaR;
887 #ifdef EDM_ML_DEBUG
888  edm::LogVerbatim("HcalSim") << "HCalSD::Hit for Detector " << det << " etaR " << etaR << " phi " << phi / deg
889  << " depth " << depth;
890 #endif
891  double time = (aStep->GetPostStepPoint()->GetGlobalTime());
892  uint32_t unitID = 0;
893  if (numberingFromDDD) {
894  HcalNumberingFromDDD::HcalID tmp = numberingFromDDD->unitID(det, etaR, phi, depth, 1);
895  unitID = setDetUnitId(tmp);
896  }
897  currentID[0].setID(unitID, time, primaryID, 1);
898 
899  edepositHAD = aStep->GetTotalEnergyDeposit();
900  edepositEM = -edepositHAD + edep;
901 #ifdef plotDebug
902  plotProfile(aStep, hitPoint, edep, time, depth);
903 #endif
904 #ifdef EDM_ML_DEBUG
905  double beta = preStepPoint->GetBeta();
906  edm::LogVerbatim("HcalSim") << "HCalSD::getHitPMT 1 hit for " << GetName() << " of " << primaryID << " with "
907  << theTrack->GetDefinition()->GetParticleName() << " of "
908  << preStepPoint->GetKineticEnergy() / GeV << " GeV with velocity " << beta << " UnitID "
909  << std::hex << unitID << std::dec;
910 #endif
911  processHit(aStep);
912  }
913 }
914 
915 void HCalSD::getHitFibreBundle(const G4Step* aStep, bool type) {
916  auto const preStepPoint = aStep->GetPreStepPoint();
917  auto const theTrack = aStep->GetTrack();
918  double edep = showerBundle->getHits(aStep, type);
919 
920  if (edep >= 0.0) {
921  edep *= GeV;
922  double etrack = preStepPoint->GetKineticEnergy();
923  int primaryID = 0;
924  if (etrack >= energyCut) {
925  primaryID = theTrack->GetTrackID();
926  } else {
927  primaryID = theTrack->GetParentID();
928  if (primaryID == 0)
929  primaryID = theTrack->GetTrackID();
930  }
931  // Reset entry point for new primary
932  resetForNewPrimary(aStep);
933  //
934  int det = static_cast<int>(HcalForward);
935  const G4ThreeVector& hitPoint = preStepPoint->GetPosition();
936  double rr = hitPoint.x() * hitPoint.x() + hitPoint.y() * hitPoint.y();
937  double phi = rr == 0. ? 0. : atan2(hitPoint.y(), hitPoint.x());
938  double etaR = showerBundle->getRadius();
939  int depth = 1;
940  if (etaR < 0.) {
941  depth = 2;
942  etaR = -etaR;
943  }
944  if (hitPoint.z() < 0.)
945  etaR = -etaR;
946 #ifdef EDM_ML_DEBUG
947  edm::LogVerbatim("HcalSim") << "HCalSD::Hit for Detector " << det << " etaR " << etaR << " phi " << phi / deg
948  << " depth " << depth;
949 #endif
950  double time = (aStep->GetPostStepPoint()->GetGlobalTime());
951  uint32_t unitID = 0;
952  if (numberingFromDDD) {
953  HcalNumberingFromDDD::HcalID tmp = numberingFromDDD->unitID(det, etaR, phi, depth, 1);
954  unitID = setDetUnitId(tmp);
955  }
956  if (type)
957  currentID[0].setID(unitID, time, primaryID, 3);
958  else
959  currentID[0].setID(unitID, time, primaryID, 2);
960 
961  edepositHAD = aStep->GetTotalEnergyDeposit();
962  edepositEM = -edepositHAD + edep;
963 #ifdef plotDebug
964  plotProfile(aStep, hitPoint, edep, time, depth);
965 #endif
966 #ifdef EDM_ML_DEBUG
967  double beta = preStepPoint->GetBeta();
968  edm::LogVerbatim("HcalSim") << "HCalSD::getHitFibreBundle 1 hit for " << GetName() << " of " << primaryID
969  << " with " << theTrack->GetDefinition()->GetParticleName() << " of "
970  << preStepPoint->GetKineticEnergy() / GeV << " GeV with velocity " << beta << " UnitID "
971  << std::hex << unitID << std::dec;
972 #endif
973  processHit(aStep);
974  } // non-zero energy deposit
975 }
976 
978  std::ifstream infile;
979  int entry = 0;
980  infile.open(fName.c_str(), std::ios::in);
981  if (infile) {
982  int det, zside, etaR, phi, lay;
983  double wt;
984  while (infile >> det >> zside >> etaR >> phi >> lay >> wt) {
985  uint32_t id = HcalTestNumbering::packHcalIndex(det, zside, 1, etaR, phi, lay);
986  layerWeights.insert(std::pair<uint32_t, double>(id, wt));
987  ++entry;
988 #ifdef EDM_ML_DEBUG
989  edm::LogVerbatim("HcalSim") << "HCalSD::readWeightFromFile:Entry " << entry << " ID " << std::hex << id
990  << std::dec << " (" << det << "/" << zside << "/1/" << etaR << "/" << phi << "/"
991  << lay << ") Weight " << wt;
992 #endif
993  }
994  infile.close();
995  }
996  edm::LogVerbatim("HcalSim") << "HCalSD::readWeightFromFile: reads " << entry << " weights from " << fName;
997  if (entry <= 0)
998  useLayerWt = false;
999 }
1000 
1001 double HCalSD::layerWeight(int det, const G4ThreeVector& pos, int depth, int lay) {
1002  double wt = 1.;
1003  if (numberingFromDDD) {
1004  //get the ID's as eta, phi, depth, ... indices
1006  numberingFromDDD->unitID(det, math::XYZVectorD(pos.x(), pos.y(), pos.z()), depth, lay);
1007  modifyDepth(tmp);
1008  uint32_t id = HcalTestNumbering::packHcalIndex(tmp.subdet, tmp.zside, 1, tmp.etaR, tmp.phis, tmp.lay);
1009  std::map<uint32_t, double>::const_iterator ite = layerWeights.find(id);
1010  if (ite != layerWeights.end())
1011  wt = ite->second;
1012 #ifdef EDM_ML_DEBUG
1013  edm::LogVerbatim("HcalSim") << "HCalSD::layerWeight: ID " << std::hex << id << std::dec << " (" << tmp.subdet << "/"
1014  << tmp.zside << "/1/" << tmp.etaR << "/" << tmp.phis << "/" << tmp.lay << ") Weight "
1015  << wt;
1016 #endif
1017  }
1018  return wt;
1019 }
1020 
1021 void HCalSD::plotProfile(const G4Step* aStep, const G4ThreeVector& global, double edep, double time, int id) {
1022  const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
1023  static const unsigned int names = 10;
1024  static const G4String modName[names] = {
1025  "HEModule", "HVQF", "HBModule", "MBAT", "MBBT", "MBBTC", "MBBT_R1P", "MBBT_R1M", "MBBT_R1PX", "MBBT_R1MX"};
1026  G4ThreeVector local;
1027  bool found = false;
1028  double depth = -2000;
1029  int idx = 4;
1030  for (int n = 0; n < touch->GetHistoryDepth(); ++n) {
1031  G4String name(static_cast<G4String>(DD4hep2DDDName::nameMatterLV(touch->GetVolume(n)->GetName(), dd4hep_)));
1032 #ifdef EDM_ML_DEBUG
1033  edm::LogVerbatim("HcalSim") << "plotProfile Depth " << n << " Name " << name;
1034 #endif
1035  for (unsigned int ii = 0; ii < names; ++ii) {
1036  if (name == modName[ii]) {
1037  found = true;
1038  int dn = touch->GetHistoryDepth() - n;
1039  local = touch->GetHistory()->GetTransform(dn).TransformPoint(global);
1040  if (ii == 0) {
1041  depth = local.z() - 4006.5;
1042  idx = 1;
1043  } else if (ii == 1) {
1044  depth = local.z() + 825.0;
1045  idx = 3;
1046  } else if (ii == 2) {
1047  depth = local.x() - 1775.;
1048  idx = 0;
1049  } else {
1050  depth = local.y() + 15.;
1051  idx = 2;
1052  }
1053  break;
1054  }
1055  }
1056  if (found)
1057  break;
1058  }
1059  if (!found)
1060  depth = std::abs(global.z()) - 11500;
1061 #ifdef EDM_ML_DEBUG
1062  edm::LogVerbatim("HcalSim") << "plotProfile Found " << found << " Global " << global << " Local " << local
1063  << " depth " << depth << " ID " << id << " EDEP " << edep << " Time " << time;
1064 #endif
1065  if (hit_[idx] != nullptr)
1066  hit_[idx]->Fill(edep);
1067  if (time_[idx] != nullptr)
1068  time_[idx]->Fill(time, edep);
1069  if (dist_[idx] != nullptr)
1070  dist_[idx]->Fill(depth, edep);
1071  int jd = 2 * idx + id - 7;
1072  if (jd >= 0 && jd < 4) {
1073  jd += 5;
1074  if (hit_[jd] != nullptr)
1075  hit_[jd]->Fill(edep);
1076  if (time_[jd] != nullptr)
1077  time_[jd]->Fill(time, edep);
1078  if (dist_[jd] != nullptr)
1079  dist_[jd]->Fill(depth, edep);
1080  }
1081 }
1082 
1083 void HCalSD::plotHF(const G4ThreeVector& hitPoint, bool emType) {
1084  double zv = std::abs(hitPoint.z()) - gpar[4];
1085  if (emType) {
1086  if (hzvem != nullptr)
1087  hzvem->Fill(zv);
1088  } else {
1089  if (hzvhad != nullptr)
1090  hzvhad->Fill(zv);
1091  }
1092 }
1093 
1095  if (id.subdet == 4) {
1096  int ieta = (id.zside == 0) ? -id.etaR : id.etaR;
1097  if (hcalConstants_->maxHFDepth(ieta, id.phis) > 2) {
1098  if (id.depth <= 2) {
1099  if (G4UniformRand() > 0.5)
1100  id.depth += 2;
1101  }
1102  }
1103  } else if ((id.subdet == 1 || id.subdet == 2) && testNumber) {
1104  id.depth = (depth_ == 0) ? 1 : 2;
1105  }
1106 }
1107 
1109 #ifdef printDebug
1110  detNull_ = {0, 0, 0, 0};
1111 #endif
1112 }
1113 
1115 #ifdef printDebug
1116  int sum = detNull_[0] + detNull_[1] + detNull_[2];
1117  if (sum > 0)
1118  edm::LogVerbatim("HcalSim") << "NullDets " << detNull_[0] << " " << detNull_[1] << " " << detNull_[2] << " "
1119  << detNull_[3] << " " << (static_cast<float>(sum) / (sum + detNull_[3]));
1120 #endif
1121 }
1122 
1123 void HCalSD::printVolume(const G4VTouchable* touch) const {
1124  if (touch) {
1125 #ifdef EDM_ML_DEBUG
1126  int level = ((touch->GetHistoryDepth()) + 1);
1127  edm::LogVerbatim("CaloSimX") << "HCalSD::printVolume with " << level << " levels";
1128  static const std::string unknown("Unknown");
1129  //Get name and copy numbers
1130  for (int ii = 0; ii < level; ii++) {
1131  int i = level - ii - 1;
1132  G4VPhysicalVolume* pv = touch->GetVolume(i);
1133  G4String name = (pv != nullptr) ? pv->GetName() : unknown;
1134  G4int copyno = touch->GetReplicaNumber(i);
1135  edm::LogVerbatim("HcalSim") << "[" << ii << "] " << name << ":" << copyno;
1136  }
1137 #endif
1138  }
1139 }
float edepositEM
Definition: CaloSD.h:144
int getNumberOfHits(int k=0)
Definition: CaloSD.cc:521
double energyCut
Definition: CaloSD.h:148
Log< level::Info, true > LogVerbatim
std::vector< std::string > matNames
Definition: HCalSD.h:122
void readWeightFromFile(const std::string &)
Definition: HCalSD.cc:977
bool isPrimary() const
std::vector< const G4LogicalVolume * > hfLV
Definition: HCalSD.h:124
bool useParam
Definition: HCalSD.h:112
double eminHitHE
Definition: HCalSD.h:113
TH1F * time_[9]
Definition: HCalSD.h:126
double kmaxNeutron
Definition: CaloSD.h:155
void hitForFibre(const G4Step *step)
Definition: HCalSD.cc:784
void initEvent(const BeginOfEvent *) override
Definition: HCalSD.cc:1108
static bool isMuon(int pdgCode)
constexpr int zside() const
get the z-side of the cell (1/-1)
Definition: HcalDetId.h:141
bool useLayerWt
Definition: HCalSD.h:109
Definition: CaloSD.h:40
std::vector< std::string > fibreNames
Definition: HCalSD.h:121
static constexpr double slopeHE_
Definition: HCalSD.h:106
std::unique_ptr< HFShowerParam > showerParam
Definition: HCalSD.h:92
double weight_
Definition: HCalSD.h:115
std::vector< double > gpar
Definition: HCalSD.h:118
bool useFibreBundle
Definition: HCalSD.h:109
constexpr int ietaAbs() const
get the absolute value of the cell ieta
Definition: HcalDetId.h:148
const HcalDDDSimConstants * hcalConstants_
Definition: HCalSD.h:96
double betaThr
Definition: HCalSD.h:111
std::vector< std::string > hfNames_
double deliveredLumi
Definition: HCalSD.h:114
int maxHFDepth(const int &ieta, const int &iphi) const
float degradation(float intlumi, int ieta, int lay) const
std::vector< std::string > hfFibreNames_
void endEvent() override
Definition: HCalSD.cc:1114
double eminHitHB
Definition: HCalSD.h:113
bool useShowerLibrary
Definition: HCalSD.h:112
std::vector< std::string > hfNames
Definition: HCalSD.h:120
void plotProfile(const G4Step *step, const G4ThreeVector &pos, double edep, double time, int id)
Definition: HCalSD.cc:1021
double birk2
Definition: HCalSD.h:111
bool usePMTHit
Definition: HCalSD.h:109
double getLayer0Wt(const int &det, const int &phi, const int &zside) const
std::vector< const G4LogicalVolume * > fibre2LV
Definition: HCalSD.h:124
int zside(DetId const &)
static const unsigned int numberOfZLayers
Definition: HFDarkening.h:24
void modifyDepth(HcalNumberingFromDDD::HcalID &id)
Definition: HCalSD.cc:1094
void update(const BeginOfJob *) override
This routine will be called when the appropriate signal arrives.
Definition: HCalSD.cc:575
void setNumberingScheme(HcalNumberingScheme *)
Definition: HCalSD.cc:568
void processHit(const G4Step *step)
Definition: CaloSD.h:117
bool agingFlagHE
Definition: HCalSD.h:108
double birk1
Definition: HCalSD.h:111
double kmaxProton
Definition: CaloSD.h:155
const std::vector< std::string_view > logicalNames(const std::string &readoutName) const
std::unique_ptr< HFShowerPMT > showerPMT
Definition: HCalSD.h:93
const std::vector< double > & getGparHF() const
const std::string names[nVars_]
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e g
Definition: Activities.doc:4
TH1F * hzvhad
Definition: HCalSD.h:126
bool isItConicalBundle(const G4LogicalVolume *)
Definition: HCalSD.cc:709
HCalSD(const std::string &, const HcalDDDSimConstants *, const HcalDDDRecConstants *, const HcalSimulationConstants *, const HBHEDarkening *, const HBHEDarkening *, const SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *)
Definition: HCalSD.cc:42
double eminHitHF
Definition: HCalSD.h:113
const HcalSimulationParameters * hcalsimpar() const
T * make(const Args &...args) const
make new ROOT object
static constexpr double minRoff_
Definition: HCalSD.h:104
static uint32_t packHcalIndex(int det, int z, int depth, int eta, int phi, int lay)
std::unique_ptr< HcalNumberingFromDDD > numberingFromDDD
Definition: HCalSD.h:88
TH1F * hit_[9]
Definition: HCalSD.h:126
double getResponseWt(const G4Track *, int k=0)
Definition: CaloSD.cc:928
ROOT::Math::DisplacementVector3D< ROOT::Math::Cartesian3D< double > > XYZVectorD
spatial vector with cartesian internal representation
Definition: Vector3D.h:8
#define EDM_ML_DEBUG
Definition: MPUnpacker.cc:1
bool agingFlagHB
Definition: HCalSD.h:108
bool testNumber
Definition: HCalSD.h:110
double kmaxIon
Definition: CaloSD.h:155
bool suppressHeavy
Definition: CaloSD.h:154
bool dd4hep_
Definition: HCalSD.h:117
std::vector< std::string > hfFibreConicalNames_
bool useBirk
Definition: HCalSD.h:109
float edepositHAD
Definition: CaloSD.h:144
std::unique_ptr< HFDarkening > m_HFDarkening
Definition: HCalSD.h:100
double birk3
Definition: HCalSD.h:111
void resetForNewPrimary(const G4Step *)
Definition: CaloSD.cc:714
Definition: tfile.py:1
const HBHEDarkening * m_HEDarkening
Definition: HCalSD.h:99
bool applyFidCut
Definition: HCalSD.h:112
TH1F * dist_[9]
Definition: HCalSD.h:126
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
void plotHF(const G4ThreeVector &pos, bool emType)
Definition: HCalSD.cc:1083
std::unique_ptr< HcalTestNS > m_HcalTestNS
Definition: HCalSD.h:101
static void unpackHcalIndex(const uint32_t &idx, int &det, int &z, int &depth, int &eta, int &phi, int &lay)
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:48
Definition: value.py:1
std::vector< int > hfLevels
Definition: HCalSD.h:119
void getFromParam(const G4Step *step, bool &isKilled)
Definition: HCalSD.cc:826
void getFromHFLibrary(const G4Step *step, bool &isKilled)
Definition: HCalSD.cc:741
uint32_t setDetUnitId(const G4Step *step) override
Definition: HCalSD.cc:535
bool isItStraightBundle(const G4LogicalVolume *)
Definition: HCalSD.cc:701
bool testNS_
Definition: HCalSD.h:110
ii
Definition: cuy.py:589
bool getFromLibrary(const G4Step *) override
Definition: HCalSD.cc:334
double eminHitHO
Definition: HCalSD.h:113
double tmaxHit
Definition: CaloSD.h:148
static bool isStableHadronIon(const G4Track *)
void setID(uint32_t unitID, double timeSlice, int trackID, uint16_t depth=0)
Definition: CaloHitID.cc:41
TH1F * hzvem
Definition: HCalSD.h:126
std::vector< const G4LogicalVolume * > pmtLV
Definition: HCalSD.h:124
int PMTNumber(const G4ThreeVector &pe_effect)
Definition: DetId.h:17
static const unsigned int lowZLimit
Definition: HFDarkening.h:27
bool filterHit(CaloG4Hit *, double) override
Definition: HCalSD.cc:579
static const unsigned int upperZLimit
Definition: HFDarkening.h:28
bool neutralDensity
Definition: HCalSD.h:110
bool isHF
Definition: HCalSD.h:107
void fillLogVolumeVector(const std::string &, const std::vector< std::string > &, std::vector< const G4LogicalVolume *> &)
Definition: HCalSD.cc:309
static constexpr double maxZ_
Definition: HCalSD.h:103
void getHitFibreBundle(const G4Step *step, bool type)
Definition: HCalSD.cc:915
std::vector< std::string > hcalMaterialNames_
double getEnergyDeposit(const G4Step *) override
Definition: HCalSD.cc:401
double getAttenuation(const G4Step *aStep, double birk1, double birk2, double birk3) const
Definition: CaloSD.cc:727
double getEnergyDeposit() const
Definition: CaloG4Hit.h:79
std::unique_ptr< HFShowerFibreBundle > showerBundle
Definition: HCalSD.h:94
std::vector< std::string > hfPMTNames_
bool isItPMT(const G4LogicalVolume *)
Definition: HCalSD.cc:693
bool isItinFidVolume(const G4ThreeVector &)
Definition: HCalSD.cc:725
uint32_t getUnitID() const
Definition: CaloG4Hit.h:66
bool isHF(int etabin, int depth)
void printVolume(const G4VTouchable *touch) const
Definition: HCalSD.cc:1123
void getHitPMT(const G4Step *step)
Definition: HCalSD.cc:856
std::vector< int > detNull_
Definition: HCalSD.h:127
std::unique_ptr< HFShowerLibrary > showerLibrary
Definition: HCalSD.h:90
double layerWeight(int, const G4ThreeVector &, int, int)
Definition: HCalSD.cc:1001
const HcalSimulationConstants * hcalSimConstants_
Definition: HCalSD.h:97
int depth_
Definition: HCalSD.h:116
HLT enums.
virtual int setTrackID(const G4Step *)
Definition: CaloSD.cc:891
std::map< uint32_t, double > layerWeights
Definition: HCalSD.h:125
bool isItHF(const G4Step *)
Definition: HCalSD.cc:656
static bool isGammaElectronPositron(int pdgCode)
bool isItScintillator(const G4Material *)
Definition: HCalSD.cc:717
std::vector< std::string > hfFibreStraightNames_
std::unique_ptr< HcalNumberingScheme > numberingScheme
Definition: HCalSD.h:89
CaloHitID currentID[2]
Definition: CaloSD.h:146
std::unique_ptr< HFShower > hfshower
Definition: HCalSD.h:91
std::vector< const G4LogicalVolume * > fibreLV
Definition: HCalSD.h:124
Log< level::Warning, false > LogWarning
std::string nameMatterLV(const std::string &name, bool dd4hep)
std::vector< const G4LogicalVolume * > fibre1LV
Definition: HCalSD.h:124
tmp
align.sh
Definition: createJobs.py:716
std::vector< const G4Material * > materials
Definition: HCalSD.h:123
int getMaxDepth(const int &type) const
void initRun() override
Definition: HCalSD.cc:577
const HBHEDarkening * m_HBDarkening
Definition: HCalSD.h:98
constexpr int iphi() const
get the cell iphi
Definition: HcalDetId.h:157
static constexpr double maxRoff_
Definition: HCalSD.h:105
bool useHF
Definition: HCalSD.h:112
bool isItFibre(const G4LogicalVolume *)
Definition: HCalSD.cc:677
void setParameterized(bool val)
Definition: CaloSD.h:114
constexpr int depth() const
get the tower depth
Definition: HcalDetId.h:164