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