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 
13 
18 
19 #include "G4LogicalVolumeStore.hh"
20 #include "G4LogicalVolume.hh"
21 #include "G4Step.hh"
22 #include "G4Track.hh"
23 
24 #include "G4SystemOfUnits.hh"
25 #include "G4PhysicalConstants.hh"
26 #include "Randomize.hh"
27 
28 #include "DD4hep/Filter.h"
29 
30 #include <iostream>
31 #include <fstream>
32 #include <iomanip>
33 #include <sstream>
34 
35 //#define EDM_ML_DEBUG
36 //#define plotDebug
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  bool 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) {
179  G4String namv(static_cast<std::string>(dd4hep::dd::noNamespace(hfNames[i])));
180  lv = nullptr;
181  for (auto lvol : *lvs) {
182  if (dd4hep::dd::noNamespace(lvol->GetName()) == 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 (static_cast<std::string>(dd4hep::dd::noNamespace((*matite)->GetName())) == 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) {
317  G4String namv(static_cast<std::string>(dd4hep::dd::noNamespace(lvnames[i])));
318  lv = nullptr;
319  for (auto lvol : *lvs) {
320  if (dd4hep::dd::noNamespace(lvol->GetName()) == 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: isHF " << isHF << " darken " << (m_HFDarkening != nullptr)
342  << " useParam " << useParam << " useShowerLibrary " << useShowerLibrary << " Muon? "
343  << G4TrackToParticleID::isMuon(track) << " electron? "
344  << G4TrackToParticleID::isGammaElectronPositron(track) << " Stable Hadron? "
346 #endif
347  if (isHF) {
348  if (m_HFDarkening) {
349  G4ThreeVector hitPoint = aStep->GetPreStepPoint()->GetPosition();
350  const double invcm = 1. / CLHEP::cm;
351  double r = hitPoint.perp() * invcm;
352  double z = std::abs(hitPoint.z()) * invcm;
353  double dose_acquired = 0.;
355  unsigned int hfZLayer = (unsigned int)((z - HFDarkening::lowZLimit) / 5);
356  if (hfZLayer >= HFDarkening::upperZLimit)
357  hfZLayer = (HFDarkening::upperZLimit - 1);
358  float normalized_lumi = m_HFDarkening->int_lumi(deliveredLumi);
359  for (int i = hfZLayer; i != HFDarkening::numberOfZLayers; ++i) {
360  dose_acquired = m_HFDarkening->dose(i, r);
361  weight_ *= m_HFDarkening->degradation(normalized_lumi * dose_acquired);
362  }
363  }
364 #ifdef EDM_ML_DEBUG
365  edm::LogVerbatim("HcalSim") << "HCalSD::getFromLibrary: HFLumiDarkening at "
366  << "r= " << r << ", z= " << z << " Dose= " << dose_acquired << " weight= " << weight_;
367 #endif
368  }
369 
370  if (useParam) {
371  getFromParam(aStep, kill);
372 #ifdef EDM_ML_DEBUG
373  G4String nameVolume = aStep->GetPreStepPoint()->GetPhysicalVolume()->GetLogicalVolume()->GetName();
374  edm::LogVerbatim("HcalSim") << "HCalSD: " << getNumberOfHits() << " hits from parametrization in " << nameVolume
375  << " for Track " << track->GetTrackID() << " ("
376  << track->GetDefinition()->GetParticleName() << ")";
377 #endif
380 #ifdef EDM_ML_DEBUG
381  auto nameVolume = aStep->GetPreStepPoint()->GetPhysicalVolume()->GetLogicalVolume()->GetName();
382  edm::LogVerbatim("HcalSim") << "HCalSD: Starts shower library from " << nameVolume << " for Track "
383  << track->GetTrackID() << " (" << track->GetDefinition()->GetParticleName() << ")";
384 
385 #endif
386  getFromHFLibrary(aStep, kill);
387  }
388  }
389  }
390 #ifdef EDM_ML_DEBUG
391  edm::LogVerbatim("HcalSim") << "HCalSD::getFromLibrary ID= " << track->GetTrackID() << " ("
392  << track->GetDefinition()->GetParticleName() << ") kill= " << kill
393  << " weight= " << weight_ << " depth= " << depth_ << " isHF: " << isHF;
394 #endif
395  return kill;
396 }
397 
398 double HCalSD::getEnergyDeposit(const G4Step* aStep) {
399  double destep(0.0);
400  auto const lv = aStep->GetPreStepPoint()->GetPhysicalVolume()->GetLogicalVolume();
401  auto const theTrack = aStep->GetTrack();
402 
403  if (isHF) {
404  if (useShowerLibrary && G4TrackToParticleID::isMuon(theTrack) && isItFibre(lv)) {
405 #ifdef EDM_ML_DEBUG
406  edm::LogVerbatim("HcalSim") << "HCalSD: Hit at Fibre in LV " << lv->GetName() << " for track "
407  << aStep->GetTrack()->GetTrackID() << " ("
408  << aStep->GetTrack()->GetDefinition()->GetParticleName() << ")";
409 #endif
410  hitForFibre(aStep);
411  }
412  return destep;
413  }
414 
415  if (isItPMT(lv)) {
416  if (usePMTHit && showerPMT) {
417  getHitPMT(aStep);
418  }
419 #ifdef EDM_ML_DEBUG
420  edm::LogVerbatim("HcalSim") << "HCalSD: Hit from PMT parametrization in LV " << lv->GetName() << " for Track "
421  << aStep->GetTrack()->GetTrackID() << " ("
422  << aStep->GetTrack()->GetDefinition()->GetParticleName() << ")";
423 #endif
424  return destep;
425 
426  } else if (isItStraightBundle(lv)) {
427  if (useFibreBundle && showerBundle) {
428  getHitFibreBundle(aStep, false);
429  }
430 #ifdef EDM_ML_DEBUG
431  edm::LogVerbatim("HcalSim") << "HCalSD: Hit from straight FibreBundle in LV: " << lv->GetName() << " for track "
432  << aStep->GetTrack()->GetTrackID() << " ("
433  << aStep->GetTrack()->GetDefinition()->GetParticleName() << ")";
434 #endif
435  return destep;
436 
437  } else if (isItConicalBundle(lv)) {
438  if (useFibreBundle && showerBundle) {
439  getHitFibreBundle(aStep, true);
440  }
441 #ifdef EDM_ML_DEBUG
442  edm::LogVerbatim("HcalSim") << "HCalSD: Hit from conical FibreBundle PV: " << lv->GetName() << " for track "
443  << aStep->GetTrack()->GetTrackID() << " ("
444  << aStep->GetTrack()->GetDefinition()->GetParticleName() << ")";
445 #endif
446  return destep;
447  }
448 
449  // normal hit
450  destep = aStep->GetTotalEnergyDeposit();
451 
452  const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
453  uint32_t detid = setDetUnitId(aStep);
454  int det(0), ieta(0), phi(0), z(0), lay, depth(-1);
455  if (testNumber) {
456  HcalTestNumbering::unpackHcalIndex(detid, det, z, depth, ieta, phi, lay);
457  if (z == 0) {
458  z = -1;
459  }
460  } else {
461  HcalDetId hcid(detid);
462  det = hcid.subdetId();
463  ieta = hcid.ietaAbs();
464  phi = hcid.iphi();
465  z = hcid.zside();
466  }
467  lay = (touch->GetReplicaNumber(0) / 10) % 100 + 1;
468 #ifdef EDM_ML_DEBUG
469  edm::LogVerbatim("HcalSim") << "HCalSD: det: " << det << " ieta: " << ieta << " iphi: " << phi << " zside " << z
470  << " lay: " << lay - 2;
471 #endif
472  if (depth_ == 0 && (det == 1 || det == 2) && ((!testNumber) || neutralDensity))
474  if (useLayerWt) {
475  G4ThreeVector hitPoint = aStep->GetPreStepPoint()->GetPosition();
476  weight_ = layerWeight(det + 2, hitPoint, depth_, lay);
477  }
478 
479  if (agingFlagHB && m_HBDarkening && det == 1) {
480  double dweight = m_HBDarkening->degradation(deliveredLumi, ieta, lay);
481  weight_ *= dweight;
482 #ifdef EDM_ML_DEBUG
483  edm::LogVerbatim("HcalSim") << "HCalSD: HB Lumi: " << deliveredLumi << " coefficient = " << dweight
484  << " Weight= " << weight_;
485 #endif
486  }
487 
488  if (agingFlagHE && m_HEDarkening && det == 2) {
489  double dweight = m_HEDarkening->degradation(deliveredLumi, ieta, lay);
490  weight_ *= dweight;
491 #ifdef EDM_ML_DEBUG
492  edm::LogVerbatim("HcalSim") << "HCalSD: HB Lumi: " << deliveredLumi << " coefficient = " << dweight
493  << " Weight= " << weight_;
494 #endif
495  }
496 
497  if (suppressHeavy) {
498  TrackInformation* trkInfo = (TrackInformation*)(theTrack->GetUserInformation());
499  if (trkInfo) {
500  int pdg = theTrack->GetDefinition()->GetPDGEncoding();
501  if (!(trkInfo->isPrimary())) { // Only secondary particles
502  double ke = theTrack->GetKineticEnergy();
503  if (pdg / 1000000000 == 1 && (pdg / 10000) % 100 > 0 && (pdg / 10) % 100 > 0 && ke < kmaxIon)
504  weight_ = 0;
505  if ((pdg == 2212) && (ke < kmaxProton))
506  weight_ = 0;
507  if ((pdg == 2112) && (ke < kmaxNeutron))
508  weight_ = 0;
509  }
510  }
511  }
512  double wt0(1.0);
513  if (useBirk) {
514  const G4Material* mat = aStep->GetPreStepPoint()->GetMaterial();
515  if (isItScintillator(mat))
516  wt0 = getAttenuation(aStep, birk1, birk2, birk3);
517  }
518  weight_ *= wt0;
519  double wt1 = getResponseWt(theTrack);
520  double wt2 = theTrack->GetWeight();
521  double edep = weight_ * wt1 * destep;
522  if (wt2 > 0.0) {
523  edep *= wt2;
524  }
525 #ifdef EDM_ML_DEBUG
526  edm::LogVerbatim("HcalSim") << "HCalSD: edep= " << edep << " Det: " << det + 2 << " depth= " << depth_
527  << " weight= " << weight_ << " wt0= " << wt0 << " wt1= " << wt1 << " wt2= " << wt2;
528 #endif
529  return edep;
530 }
531 
532 uint32_t HCalSD::setDetUnitId(const G4Step* aStep) {
533  auto const prePoint = aStep->GetPreStepPoint();
534  auto const touch = prePoint->GetTouchable();
535  const G4ThreeVector& hitPoint = prePoint->GetPosition();
536 
537  int depth = (touch->GetReplicaNumber(0)) % 10 + 1;
538  int lay = (touch->GetReplicaNumber(0) / 10) % 100 + 1;
539  int det = (touch->GetReplicaNumber(1)) / 1000;
540 
541  return setDetUnitId(det, hitPoint, depth, lay);
542 }
543 
545  if (scheme != nullptr) {
546  edm::LogVerbatim("HcalSim") << "HCalSD: updates numbering scheme for " << GetName();
547  numberingScheme.reset(scheme);
548  }
549 }
550 
551 void HCalSD::update(const BeginOfJob* job) {}
552 
554 
555 bool HCalSD::filterHit(CaloG4Hit* aHit, double time) {
556  double threshold = 0;
557  DetId theId(aHit->getUnitID());
558  switch (theId.subdetId()) {
559  case HcalBarrel:
561  break;
562  case HcalEndcap:
564  break;
565  case HcalOuter:
567  break;
568  case HcalForward:
570  break;
571  default:
572  break;
573  }
574  return ((time <= tmaxHit) && (aHit->getEnergyDeposit() > threshold));
575 }
576 
577 uint32_t HCalSD::setDetUnitId(int det, const G4ThreeVector& pos, int depth, int lay = 1) {
578  uint32_t id = 0;
579  if (det == 0) {
580  if (std::abs(pos.z()) > maxZ_) {
581  det = 5;
582  } else if (!(hcalConstants_->isHE())) {
583  det = 3;
584  } else {
585  double minR = minRoff_ + slopeHE_ * std::abs(pos.z());
586  double maxR = maxRoff_ + slopeHE_ * std::abs(pos.z());
587  det = ((pos.perp() > minR) && (pos.perp() < maxR)) ? 4 : 3;
588  }
589 #ifdef EDM_ML_DEBUG
590  edm::LogVerbatim("HcalSim") << "Position " << pos.perp() << ":" << std::abs(pos.z()) << " Limits "
591  << !(hcalConstants_->isHE()) << ":" << maxZ_ << " det " << det;
592 #endif
593  }
594  if (numberingFromDDD.get()) {
595  //get the ID's as eta, phi, depth, ... indices
597  numberingFromDDD->unitID(det, math::XYZVectorD(pos.x(), pos.y(), pos.z()), depth, lay);
598  id = setDetUnitId(tmp);
599  }
600  return id;
601 }
602 
604  modifyDepth(tmp);
605  uint32_t id = (numberingScheme.get()) ? numberingScheme->getUnitID(tmp) : 0;
606  if ((!testNumber) && m_HcalTestNS.get()) {
607  bool ok = m_HcalTestNS->compare(tmp, id);
608  if (!ok)
609  edm::LogWarning("HcalSim") << "Det ID from HCalSD " << HcalDetId(id) << " " << std::hex << id << std::dec
610  << " does not match one from relabller for " << tmp.subdet << ":" << tmp.etaR << ":"
611  << tmp.phi << ":" << tmp.phis << ":" << tmp.depth << ":" << tmp.lay << std::endl;
612  }
613  return id;
614 }
615 
616 bool HCalSD::isItHF(const G4Step* aStep) {
617  const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
618  int levels = (touch->GetHistoryDepth()) + 1;
619  for (unsigned int it = 0; it < hfNames.size(); ++it) {
620  if (levels >= hfLevels[it]) {
621  const G4LogicalVolume* lv = touch->GetVolume(levels - hfLevels[it])->GetLogicalVolume();
622  if (lv == hfLV[it])
623  return true;
624  }
625  }
626  return false;
627 }
628 
629 bool HCalSD::isItHF(const G4String& name) {
630  for (const auto& nam : hfNames)
631  if (name == static_cast<G4String>(nam)) {
632  return true;
633  }
634  return false;
635 }
636 
637 bool HCalSD::isItFibre(const G4LogicalVolume* lv) {
638  for (auto lvol : fibreLV)
639  if (lv == lvol) {
640  return true;
641  }
642  return false;
643 }
644 
645 bool HCalSD::isItFibre(const G4String& name) {
646  for (const auto& nam : fibreNames)
647  if (name == static_cast<G4String>(nam)) {
648  return true;
649  }
650  return false;
651 }
652 
653 bool HCalSD::isItPMT(const G4LogicalVolume* lv) {
654  for (auto lvol : pmtLV)
655  if (lv == lvol) {
656  return true;
657  }
658  return false;
659 }
660 
661 bool HCalSD::isItStraightBundle(const G4LogicalVolume* lv) {
662  for (auto lvol : fibre1LV)
663  if (lv == lvol) {
664  return true;
665  }
666  return false;
667 }
668 
669 bool HCalSD::isItConicalBundle(const G4LogicalVolume* lv) {
670  for (auto lvol : fibre2LV)
671  if (lv == lvol) {
672  return true;
673  }
674  return false;
675 }
676 
677 bool HCalSD::isItScintillator(const G4Material* mat) {
678  for (auto amat : materials)
679  if (amat == mat) {
680  return true;
681  }
682  return false;
683 }
684 
685 bool HCalSD::isItinFidVolume(const G4ThreeVector& hitPoint) {
686  bool flag = true;
687  if (applyFidCut) {
688  int npmt = HFFibreFiducial::PMTNumber(hitPoint);
689 #ifdef EDM_ML_DEBUG
690  edm::LogVerbatim("HcalSim") << "HCalSD::isItinFidVolume:#PMT= " << npmt << " for hit point " << hitPoint;
691 #endif
692  if (npmt <= 0)
693  flag = false;
694  }
695 #ifdef EDM_ML_DEBUG
696  edm::LogVerbatim("HcalSim") << "HCalSD::isItinFidVolume: point " << hitPoint << " return flag " << flag;
697 #endif
698  return flag;
699 }
700 
701 void HCalSD::getFromHFLibrary(const G4Step* aStep, bool& isKilled) {
702  std::vector<HFShowerLibrary::Hit> hits = showerLibrary->getHits(aStep, isKilled, weight_, false);
703  if (!isKilled || hits.empty()) {
704  return;
705  }
706 
707  int primaryID = setTrackID(aStep);
708 
709  // Reset entry point for new primary
710  resetForNewPrimary(aStep);
711 
712  auto const theTrack = aStep->GetTrack();
713  int det = 5;
714 
716  edepositEM = 1. * GeV;
717  edepositHAD = 0.;
718  } else {
719  edepositEM = 0.;
720  edepositHAD = 1. * GeV;
721  }
722 #ifdef EDM_ML_DEBUG
723  edm::LogVerbatim("HcalSim") << "HCalSD::getFromLibrary " << hits.size() << " hits for " << GetName() << " of "
724  << primaryID << " with " << theTrack->GetDefinition()->GetParticleName() << " of "
725  << aStep->GetPreStepPoint()->GetKineticEnergy() / GeV << " GeV";
726 #endif
727  for (unsigned int i = 0; i < hits.size(); ++i) {
728  G4ThreeVector hitPoint = hits[i].position;
729  if (isItinFidVolume(hitPoint)) {
730  int depth = hits[i].depth;
731  double time = hits[i].time;
732  unsigned int unitID = setDetUnitId(det, hitPoint, depth);
733  currentID.setID(unitID, time, primaryID, 0);
734 #ifdef plotDebug
735  plotProfile(aStep, hitPoint, 1.0 * GeV, time, depth);
736  bool emType = G4TrackToParticleID::isGammaElectronPositron(theTrack->GetDefinition()->GetPDGEncoding());
737  plotHF(hitPoint, emType);
738 #endif
739  processHit(aStep);
740  }
741  }
742 }
743 
744 void HCalSD::hitForFibre(const G4Step* aStep) { // if not ParamShower
745 
746  std::vector<HFShower::Hit> hits = hfshower->getHits(aStep, weight_);
747  if (hits.empty()) {
748  return;
749  }
750 
751  auto const theTrack = aStep->GetTrack();
752  int primaryID = setTrackID(aStep);
753  int det = 5;
754 
756  edepositEM = 1. * GeV;
757  edepositHAD = 0.;
758  } else {
759  edepositEM = 0.;
760  edepositHAD = 1. * GeV;
761  }
762 
763 #ifdef EDM_ML_DEBUG
764  edm::LogVerbatim("HcalSim") << "HCalSD::hitForFibre " << hits.size() << " hits for " << GetName() << " of "
765  << primaryID << " with " << theTrack->GetDefinition()->GetParticleName() << " of "
766  << aStep->GetPreStepPoint()->GetKineticEnergy() / GeV << " GeV in detector type " << det;
767 #endif
768 
769  for (unsigned int i = 0; i < hits.size(); ++i) {
770  G4ThreeVector hitPoint = hits[i].position;
771  if (isItinFidVolume(hitPoint)) {
772  int depth = hits[i].depth;
773  double time = hits[i].time;
774  unsigned int unitID = setDetUnitId(det, hitPoint, depth);
775  currentID.setID(unitID, time, primaryID, 0);
776 #ifdef plotDebug
777  plotProfile(aStep, hitPoint, edepositEM, time, depth);
778  bool emType = (edepositEM > 0.) ? true : false;
779  plotHF(hitPoint, emType);
780 #endif
781  processHit(aStep);
782  }
783  }
784 }
785 
786 void HCalSD::getFromParam(const G4Step* aStep, bool& isKilled) {
787  std::vector<HFShowerParam::Hit> hits = showerParam->getHits(aStep, weight_, isKilled);
788  if (!isKilled || hits.empty()) {
789  return;
790  }
791 
792  int primaryID = setTrackID(aStep);
793  int det = 5;
794 
795 #ifdef EDM_ML_DEBUG
796  edm::LogVerbatim("HcalSim") << "HCalSD::getFromParam " << hits.size() << " hits for " << GetName() << " of "
797  << primaryID << " with " << aStep->GetTrack()->GetDefinition()->GetParticleName()
798  << " of " << aStep->GetPreStepPoint()->GetKineticEnergy() / GeV
799  << " GeV in detector type " << det;
800 #endif
801  for (unsigned int i = 0; i < hits.size(); ++i) {
802  G4ThreeVector hitPoint = hits[i].position;
803  int depth = hits[i].depth;
804  double time = hits[i].time;
805  unsigned int unitID = setDetUnitId(det, hitPoint, depth);
806  currentID.setID(unitID, time, primaryID, 0);
807  edepositEM = hits[i].edep * GeV;
808  edepositHAD = 0.;
809 #ifdef plotDebug
810  plotProfile(aStep, hitPoint, edepositEM, time, depth);
811 #endif
812  processHit(aStep);
813  }
814 }
815 
816 void HCalSD::getHitPMT(const G4Step* aStep) {
817  auto const preStepPoint = aStep->GetPreStepPoint();
818  auto const theTrack = aStep->GetTrack();
819  double edep = showerPMT->getHits(aStep);
820 
821  if (edep >= 0.) {
822  edep *= GeV;
823  double etrack = preStepPoint->GetKineticEnergy();
824  int primaryID = 0;
825  if (etrack >= energyCut) {
826  primaryID = theTrack->GetTrackID();
827  } else {
828  primaryID = theTrack->GetParentID();
829  if (primaryID == 0)
830  primaryID = theTrack->GetTrackID();
831  }
832  // Reset entry point for new primary
833  resetForNewPrimary(aStep);
834  //
835  int det = static_cast<int>(HcalForward);
836  const G4ThreeVector& hitPoint = preStepPoint->GetPosition();
837  double rr = (hitPoint.x() * hitPoint.x() + hitPoint.y() * hitPoint.y());
838  double phi = (rr == 0. ? 0. : atan2(hitPoint.y(), hitPoint.x()));
839  double etaR = showerPMT->getRadius();
840  int depth = 1;
841  if (etaR < 0) {
842  depth = 2;
843  etaR = -etaR;
844  }
845  if (hitPoint.z() < 0)
846  etaR = -etaR;
847 #ifdef EDM_ML_DEBUG
848  edm::LogVerbatim("HcalSim") << "HCalSD::Hit for Detector " << det << " etaR " << etaR << " phi " << phi / deg
849  << " depth " << depth;
850 #endif
851  double time = (aStep->GetPostStepPoint()->GetGlobalTime());
852  uint32_t unitID = 0;
853  if (numberingFromDDD) {
854  HcalNumberingFromDDD::HcalID tmp = numberingFromDDD->unitID(det, etaR, phi, depth, 1);
855  unitID = setDetUnitId(tmp);
856  }
857  currentID.setID(unitID, time, primaryID, 1);
858 
859  edepositHAD = aStep->GetTotalEnergyDeposit();
860  edepositEM = -edepositHAD + edep;
861 #ifdef plotDebug
862  plotProfile(aStep, hitPoint, edep, time, depth);
863 #endif
864 #ifdef EDM_ML_DEBUG
865  double beta = preStepPoint->GetBeta();
866  edm::LogVerbatim("HcalSim") << "HCalSD::getHitPMT 1 hit for " << GetName() << " of " << primaryID << " with "
867  << theTrack->GetDefinition()->GetParticleName() << " of "
868  << preStepPoint->GetKineticEnergy() / GeV << " GeV with velocity " << beta << " UnitID "
869  << std::hex << unitID << std::dec;
870 #endif
871  processHit(aStep);
872  }
873 }
874 
875 void HCalSD::getHitFibreBundle(const G4Step* aStep, bool type) {
876  auto const preStepPoint = aStep->GetPreStepPoint();
877  auto const theTrack = aStep->GetTrack();
878  double edep = showerBundle->getHits(aStep, type);
879 
880  if (edep >= 0.0) {
881  edep *= GeV;
882  double etrack = preStepPoint->GetKineticEnergy();
883  int primaryID = 0;
884  if (etrack >= energyCut) {
885  primaryID = theTrack->GetTrackID();
886  } else {
887  primaryID = theTrack->GetParentID();
888  if (primaryID == 0)
889  primaryID = theTrack->GetTrackID();
890  }
891  // Reset entry point for new primary
892  resetForNewPrimary(aStep);
893  //
894  int det = static_cast<int>(HcalForward);
895  const G4ThreeVector& hitPoint = preStepPoint->GetPosition();
896  double rr = hitPoint.x() * hitPoint.x() + hitPoint.y() * hitPoint.y();
897  double phi = rr == 0. ? 0. : atan2(hitPoint.y(), hitPoint.x());
898  double etaR = showerBundle->getRadius();
899  int depth = 1;
900  if (etaR < 0.) {
901  depth = 2;
902  etaR = -etaR;
903  }
904  if (hitPoint.z() < 0.)
905  etaR = -etaR;
906 #ifdef EDM_ML_DEBUG
907  edm::LogVerbatim("HcalSim") << "HCalSD::Hit for Detector " << det << " etaR " << etaR << " phi " << phi / deg
908  << " depth " << depth;
909 #endif
910  double time = (aStep->GetPostStepPoint()->GetGlobalTime());
911  uint32_t unitID = 0;
912  if (numberingFromDDD) {
913  HcalNumberingFromDDD::HcalID tmp = numberingFromDDD->unitID(det, etaR, phi, depth, 1);
914  unitID = setDetUnitId(tmp);
915  }
916  if (type)
917  currentID.setID(unitID, time, primaryID, 3);
918  else
919  currentID.setID(unitID, time, primaryID, 2);
920 
921  edepositHAD = aStep->GetTotalEnergyDeposit();
922  edepositEM = -edepositHAD + edep;
923 #ifdef plotDebug
924  plotProfile(aStep, hitPoint, edep, time, depth);
925 #endif
926 #ifdef EDM_ML_DEBUG
927  double beta = preStepPoint->GetBeta();
928  edm::LogVerbatim("HcalSim") << "HCalSD::getHitFibreBundle 1 hit for " << GetName() << " of " << primaryID
929  << " with " << theTrack->GetDefinition()->GetParticleName() << " of "
930  << preStepPoint->GetKineticEnergy() / GeV << " GeV with velocity " << beta << " UnitID "
931  << std::hex << unitID << std::dec;
932 #endif
933  processHit(aStep);
934  } // non-zero energy deposit
935 }
936 
938  std::ifstream infile;
939  int entry = 0;
940  infile.open(fName.c_str(), std::ios::in);
941  if (infile) {
942  int det, zside, etaR, phi, lay;
943  double wt;
944  while (infile >> det >> zside >> etaR >> phi >> lay >> wt) {
945  uint32_t id = HcalTestNumbering::packHcalIndex(det, zside, 1, etaR, phi, lay);
946  layerWeights.insert(std::pair<uint32_t, double>(id, wt));
947  ++entry;
948 #ifdef EDM_ML_DEBUG
949  edm::LogVerbatim("HcalSim") << "HCalSD::readWeightFromFile:Entry " << entry << " ID " << std::hex << id
950  << std::dec << " (" << det << "/" << zside << "/1/" << etaR << "/" << phi << "/"
951  << lay << ") Weight " << wt;
952 #endif
953  }
954  infile.close();
955  }
956  edm::LogVerbatim("HcalSim") << "HCalSD::readWeightFromFile: reads " << entry << " weights from " << fName;
957  if (entry <= 0)
958  useLayerWt = false;
959 }
960 
961 double HCalSD::layerWeight(int det, const G4ThreeVector& pos, int depth, int lay) {
962  double wt = 1.;
963  if (numberingFromDDD) {
964  //get the ID's as eta, phi, depth, ... indices
966  numberingFromDDD->unitID(det, math::XYZVectorD(pos.x(), pos.y(), pos.z()), depth, lay);
967  modifyDepth(tmp);
968  uint32_t id = HcalTestNumbering::packHcalIndex(tmp.subdet, tmp.zside, 1, tmp.etaR, tmp.phis, tmp.lay);
969  std::map<uint32_t, double>::const_iterator ite = layerWeights.find(id);
970  if (ite != layerWeights.end())
971  wt = ite->second;
972 #ifdef EDM_ML_DEBUG
973  edm::LogVerbatim("HcalSim") << "HCalSD::layerWeight: ID " << std::hex << id << std::dec << " (" << tmp.subdet << "/"
974  << tmp.zside << "/1/" << tmp.etaR << "/" << tmp.phis << "/" << tmp.lay << ") Weight "
975  << wt;
976 #endif
977  }
978  return wt;
979 }
980 
981 void HCalSD::plotProfile(const G4Step* aStep, const G4ThreeVector& global, double edep, double time, int id) {
982  const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
983  static const unsigned int names = 10;
984  static const G4String modName[names] = {
985  "HEModule", "HVQF", "HBModule", "MBAT", "MBBT", "MBBTC", "MBBT_R1P", "MBBT_R1M", "MBBT_R1PX", "MBBT_R1MX"};
986  G4ThreeVector local;
987  bool found = false;
988  double depth = -2000;
989  int idx = 4;
990  for (int n = 0; n < touch->GetHistoryDepth(); ++n) {
991  G4String name(static_cast<std::string>(dd4hep::dd::noNamespace(touch->GetVolume(n)->GetName())));
992 #ifdef EDM_ML_DEBUG
993  edm::LogVerbatim("HcalSim") << "plotProfile Depth " << n << " Name " << name;
994 #endif
995  for (unsigned int ii = 0; ii < names; ++ii) {
996  if (name == modName[ii]) {
997  found = true;
998  int dn = touch->GetHistoryDepth() - n;
999  local = touch->GetHistory()->GetTransform(dn).TransformPoint(global);
1000  if (ii == 0) {
1001  depth = local.z() - 4006.5;
1002  idx = 1;
1003  } else if (ii == 1) {
1004  depth = local.z() + 825.0;
1005  idx = 3;
1006  } else if (ii == 2) {
1007  depth = local.x() - 1775.;
1008  idx = 0;
1009  } else {
1010  depth = local.y() + 15.;
1011  idx = 2;
1012  }
1013  break;
1014  }
1015  }
1016  if (found)
1017  break;
1018  }
1019  if (!found)
1020  depth = std::abs(global.z()) - 11500;
1021 #ifdef EDM_ML_DEBUG
1022  edm::LogVerbatim("HcalSim") << "plotProfile Found " << found << " Global " << global << " Local " << local
1023  << " depth " << depth << " ID " << id << " EDEP " << edep << " Time " << time;
1024 #endif
1025  if (hit_[idx] != nullptr)
1026  hit_[idx]->Fill(edep);
1027  if (time_[idx] != nullptr)
1028  time_[idx]->Fill(time, edep);
1029  if (dist_[idx] != nullptr)
1030  dist_[idx]->Fill(depth, edep);
1031  int jd = 2 * idx + id - 7;
1032  if (jd >= 0 && jd < 4) {
1033  jd += 5;
1034  if (hit_[jd] != nullptr)
1035  hit_[jd]->Fill(edep);
1036  if (time_[jd] != nullptr)
1037  time_[jd]->Fill(time, edep);
1038  if (dist_[jd] != nullptr)
1039  dist_[jd]->Fill(depth, edep);
1040  }
1041 }
1042 
1043 void HCalSD::plotHF(const G4ThreeVector& hitPoint, bool emType) {
1044  double zv = std::abs(hitPoint.z()) - gpar[4];
1045  if (emType) {
1046  if (hzvem != nullptr)
1047  hzvem->Fill(zv);
1048  } else {
1049  if (hzvhad != nullptr)
1050  hzvhad->Fill(zv);
1051  }
1052 }
1053 
1055  if (id.subdet == 4) {
1056  int ieta = (id.zside == 0) ? -id.etaR : id.etaR;
1057  if (hcalConstants_->maxHFDepth(ieta, id.phis) > 2) {
1058  if (id.depth <= 2) {
1059  if (G4UniformRand() > 0.5)
1060  id.depth += 2;
1061  }
1062  }
1063  } else if ((id.subdet == 1 || id.subdet == 2) && testNumber) {
1064  id.depth = (depth_ == 0) ? 1 : 2;
1065  }
1066 }
float edepositEM
Definition: CaloSD.h:140
double energyCut
Definition: CaloSD.h:144
Log< level::Info, true > LogVerbatim
std::vector< std::string > matNames
Definition: HCalSD.h:118
void readWeightFromFile(const std::string &)
Definition: HCalSD.cc:937
bool isPrimary() const
std::vector< const G4LogicalVolume * > hfLV
Definition: HCalSD.h:120
bool useParam
Definition: HCalSD.h:109
double eminHitHE
Definition: HCalSD.h:110
TH1F * time_[9]
Definition: HCalSD.h:122
double kmaxNeutron
Definition: CaloSD.h:149
void hitForFibre(const G4Step *step)
Definition: HCalSD.cc:744
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:106
Definition: CaloSD.h:40
std::vector< std::string > fibreNames
Definition: HCalSD.h:117
static constexpr double slopeHE_
Definition: HCalSD.h:103
std::unique_ptr< HFShowerParam > showerParam
Definition: HCalSD.h:89
double weight_
Definition: HCalSD.h:112
std::vector< double > gpar
Definition: HCalSD.h:114
bool useFibreBundle
Definition: HCalSD.h:106
constexpr int ietaAbs() const
get the absolute value of the cell ieta
Definition: HcalDetId.h:148
const HcalDDDSimConstants * hcalConstants_
Definition: HCalSD.h:93
double betaThr
Definition: HCalSD.h:108
std::vector< std::string > hfNames_
double deliveredLumi
Definition: HCalSD.h:111
int maxHFDepth(const int &ieta, const int &iphi) const
float degradation(float intlumi, int ieta, int lay) const
float *__restrict__ zv
std::vector< std::string > hfFibreNames_
double eminHitHB
Definition: HCalSD.h:110
bool useShowerLibrary
Definition: HCalSD.h:109
std::vector< std::string > hfNames
Definition: HCalSD.h:116
void plotProfile(const G4Step *step, const G4ThreeVector &pos, double edep, double time, int id)
Definition: HCalSD.cc:981
double birk2
Definition: HCalSD.h:108
bool usePMTHit
Definition: HCalSD.h:106
double getLayer0Wt(const int &det, const int &phi, const int &zside) const
std::vector< const G4LogicalVolume * > fibre2LV
Definition: HCalSD.h:120
int zside(DetId const &)
static const unsigned int numberOfZLayers
Definition: HFDarkening.h:24
void modifyDepth(HcalNumberingFromDDD::HcalID &id)
Definition: HCalSD.cc:1054
void update(const BeginOfJob *) override
This routine will be called when the appropriate signal arrives.
Definition: HCalSD.cc:551
void setNumberingScheme(HcalNumberingScheme *)
Definition: HCalSD.cc:544
void processHit(const G4Step *step)
Definition: CaloSD.h:113
bool agingFlagHE
Definition: HCalSD.h:105
double birk1
Definition: HCalSD.h:108
double kmaxProton
Definition: CaloSD.h:149
const std::vector< std::string_view > logicalNames(const std::string &readoutName) const
std::unique_ptr< HFShowerPMT > showerPMT
Definition: HCalSD.h:90
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:122
bool isItConicalBundle(const G4LogicalVolume *)
Definition: HCalSD.cc:669
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:110
const HcalSimulationParameters * hcalsimpar() const
T * make(const Args &...args) const
make new ROOT object
static constexpr double minRoff_
Definition: HCalSD.h:101
static uint32_t packHcalIndex(int det, int z, int depth, int eta, int phi, int lay)
std::unique_ptr< HcalNumberingFromDDD > numberingFromDDD
Definition: HCalSD.h:85
TH1F * hit_[9]
Definition: HCalSD.h:122
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:105
bool testNumber
Definition: HCalSD.h:107
double kmaxIon
Definition: CaloSD.h:149
bool suppressHeavy
Definition: CaloSD.h:148
std::vector< std::string > hfFibreConicalNames_
bool useBirk
Definition: HCalSD.h:106
float edepositHAD
Definition: CaloSD.h:140
std::unique_ptr< HFDarkening > m_HFDarkening
Definition: HCalSD.h:97
double birk3
Definition: HCalSD.h:108
void resetForNewPrimary(const G4Step *)
Definition: CaloSD.cc:653
Definition: tfile.py:1
const HBHEDarkening * m_HEDarkening
Definition: HCalSD.h:96
bool applyFidCut
Definition: HCalSD.h:109
TH1F * dist_[9]
Definition: HCalSD.h:122
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
void plotHF(const G4ThreeVector &pos, bool emType)
Definition: HCalSD.cc:1043
std::unique_ptr< HcalTestNS > m_HcalTestNS
Definition: HCalSD.h:98
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
int getNumberOfHits()
Definition: CaloSD.cc:461
Definition: value.py:1
std::vector< int > hfLevels
Definition: HCalSD.h:115
void getFromParam(const G4Step *step, bool &isKilled)
Definition: HCalSD.cc:786
void getFromHFLibrary(const G4Step *step, bool &isKilled)
Definition: HCalSD.cc:701
uint32_t setDetUnitId(const G4Step *step) override
Definition: HCalSD.cc:532
bool isItStraightBundle(const G4LogicalVolume *)
Definition: HCalSD.cc:661
bool testNS_
Definition: HCalSD.h:107
ii
Definition: cuy.py:589
bool getFromLibrary(const G4Step *) override
Definition: HCalSD.cc:334
double eminHitHO
Definition: HCalSD.h:110
double tmaxHit
Definition: CaloSD.h:144
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:122
std::vector< const G4LogicalVolume * > pmtLV
Definition: HCalSD.h:120
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:555
CaloHitID currentID
Definition: CaloSD.h:142
static const unsigned int upperZLimit
Definition: HFDarkening.h:28
bool neutralDensity
Definition: HCalSD.h:107
bool isHF
Definition: HCalSD.h:104
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:100
void getHitFibreBundle(const G4Step *step, bool type)
Definition: HCalSD.cc:875
std::vector< std::string > hcalMaterialNames_
double getEnergyDeposit(const G4Step *) override
Definition: HCalSD.cc:398
double getAttenuation(const G4Step *aStep, double birk1, double birk2, double birk3) const
Definition: CaloSD.cc:666
double getEnergyDeposit() const
Definition: CaloG4Hit.h:79
std::unique_ptr< HFShowerFibreBundle > showerBundle
Definition: HCalSD.h:91
std::vector< std::string > hfPMTNames_
bool isItPMT(const G4LogicalVolume *)
Definition: HCalSD.cc:653
bool isItinFidVolume(const G4ThreeVector &)
Definition: HCalSD.cc:685
uint32_t getUnitID() const
Definition: CaloG4Hit.h:66
bool isHF(int etabin, int depth)
void getHitPMT(const G4Step *step)
Definition: HCalSD.cc:816
std::unique_ptr< HFShowerLibrary > showerLibrary
Definition: HCalSD.h:87
double layerWeight(int, const G4ThreeVector &, int, int)
Definition: HCalSD.cc:961
const HcalSimulationConstants * hcalSimConstants_
Definition: HCalSD.h:94
int depth_
Definition: HCalSD.h:113
HLT enums.
virtual int setTrackID(const G4Step *)
Definition: CaloSD.cc:823
std::map< uint32_t, double > layerWeights
Definition: HCalSD.h:121
bool isItHF(const G4Step *)
Definition: HCalSD.cc:616
static bool isGammaElectronPositron(int pdgCode)
bool isItScintillator(const G4Material *)
Definition: HCalSD.cc:677
std::vector< std::string > hfFibreStraightNames_
std::unique_ptr< HcalNumberingScheme > numberingScheme
Definition: HCalSD.h:86
std::unique_ptr< HFShower > hfshower
Definition: HCalSD.h:88
std::vector< const G4LogicalVolume * > fibreLV
Definition: HCalSD.h:120
Log< level::Warning, false > LogWarning
double getResponseWt(const G4Track *)
Definition: CaloSD.cc:858
std::vector< const G4LogicalVolume * > fibre1LV
Definition: HCalSD.h:120
tmp
align.sh
Definition: createJobs.py:716
std::vector< const G4Material * > materials
Definition: HCalSD.h:119
int getMaxDepth(const int &type) const
void initRun() override
Definition: HCalSD.cc:553
const HBHEDarkening * m_HBDarkening
Definition: HCalSD.h:95
constexpr int iphi() const
get the cell iphi
Definition: HcalDetId.h:157
static constexpr double maxRoff_
Definition: HCalSD.h:102
bool useHF
Definition: HCalSD.h:109
bool isItFibre(const G4LogicalVolume *)
Definition: HCalSD.cc:637
void setParameterized(bool val)
Definition: CaloSD.h:110