CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
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");
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 
144  HcalNumberingScheme* scheme;
145  if (testNumber || forTBH2) {
146  scheme = dynamic_cast<HcalNumberingScheme*>(new HcalTestNumberingScheme(forTBH2));
147  } else {
148  scheme = new HcalNumberingScheme();
149  }
150  setNumberingScheme(scheme);
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) {
231  readWeightFromFile(file);
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:
560  threshold = eminHitHB;
561  break;
562  case HcalEndcap:
563  threshold = eminHitHE;
564  break;
565  case HcalOuter:
566  threshold = eminHitHO;
567  break;
568  case HcalForward:
569  threshold = eminHitHF;
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 (numberingFromDDD.get()) {
580  //get the ID's as eta, phi, depth, ... indices
582  numberingFromDDD->unitID(det, math::XYZVectorD(pos.x(), pos.y(), pos.z()), depth, lay);
583  id = setDetUnitId(tmp);
584  }
585  return id;
586 }
587 
589  modifyDepth(tmp);
590  uint32_t id = (numberingScheme.get()) ? numberingScheme->getUnitID(tmp) : 0;
591  if ((!testNumber) && m_HcalTestNS.get()) {
592  bool ok = m_HcalTestNS->compare(tmp, id);
593  if (!ok)
594  edm::LogWarning("HcalSim") << "Det ID from HCalSD " << HcalDetId(id) << " " << std::hex << id << std::dec
595  << " does not match one from relabller for " << tmp.subdet << ":" << tmp.etaR << ":"
596  << tmp.phi << ":" << tmp.phis << ":" << tmp.depth << ":" << tmp.lay << std::endl;
597  }
598  return id;
599 }
600 
601 bool HCalSD::isItHF(const G4Step* aStep) {
602  const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
603  int levels = (touch->GetHistoryDepth()) + 1;
604  for (unsigned int it = 0; it < hfNames.size(); ++it) {
605  if (levels >= hfLevels[it]) {
606  const G4LogicalVolume* lv = touch->GetVolume(levels - hfLevels[it])->GetLogicalVolume();
607  if (lv == hfLV[it])
608  return true;
609  }
610  }
611  return false;
612 }
613 
614 bool HCalSD::isItHF(const G4String& name) {
615  for (const auto& nam : hfNames)
616  if (name == static_cast<G4String>(nam)) {
617  return true;
618  }
619  return false;
620 }
621 
622 bool HCalSD::isItFibre(const G4LogicalVolume* lv) {
623  for (auto lvol : fibreLV)
624  if (lv == lvol) {
625  return true;
626  }
627  return false;
628 }
629 
630 bool HCalSD::isItFibre(const G4String& name) {
631  for (const auto& nam : fibreNames)
632  if (name == static_cast<G4String>(nam)) {
633  return true;
634  }
635  return false;
636 }
637 
638 bool HCalSD::isItPMT(const G4LogicalVolume* lv) {
639  for (auto lvol : pmtLV)
640  if (lv == lvol) {
641  return true;
642  }
643  return false;
644 }
645 
646 bool HCalSD::isItStraightBundle(const G4LogicalVolume* lv) {
647  for (auto lvol : fibre1LV)
648  if (lv == lvol) {
649  return true;
650  }
651  return false;
652 }
653 
654 bool HCalSD::isItConicalBundle(const G4LogicalVolume* lv) {
655  for (auto lvol : fibre2LV)
656  if (lv == lvol) {
657  return true;
658  }
659  return false;
660 }
661 
662 bool HCalSD::isItScintillator(const G4Material* mat) {
663  for (auto amat : materials)
664  if (amat == mat) {
665  return true;
666  }
667  return false;
668 }
669 
670 bool HCalSD::isItinFidVolume(const G4ThreeVector& hitPoint) {
671  bool flag = true;
672  if (applyFidCut) {
673  int npmt = HFFibreFiducial::PMTNumber(hitPoint);
674 #ifdef EDM_ML_DEBUG
675  edm::LogVerbatim("HcalSim") << "HCalSD::isItinFidVolume:#PMT= " << npmt << " for hit point " << hitPoint;
676 #endif
677  if (npmt <= 0)
678  flag = false;
679  }
680 #ifdef EDM_ML_DEBUG
681  edm::LogVerbatim("HcalSim") << "HCalSD::isItinFidVolume: point " << hitPoint << " return flag " << flag;
682 #endif
683  return flag;
684 }
685 
686 void HCalSD::getFromHFLibrary(const G4Step* aStep, bool& isKilled) {
687  std::vector<HFShowerLibrary::Hit> hits = showerLibrary->getHits(aStep, isKilled, weight_, false);
688  if (!isKilled || hits.empty()) {
689  return;
690  }
691 
692  int primaryID = setTrackID(aStep);
693 
694  // Reset entry point for new primary
695  resetForNewPrimary(aStep);
696 
697  auto const theTrack = aStep->GetTrack();
698  int det = 5;
699 
701  edepositEM = 1. * GeV;
702  edepositHAD = 0.;
703  } else {
704  edepositEM = 0.;
705  edepositHAD = 1. * GeV;
706  }
707 #ifdef EDM_ML_DEBUG
708  edm::LogVerbatim("HcalSim") << "HCalSD::getFromLibrary " << hits.size() << " hits for " << GetName() << " of "
709  << primaryID << " with " << theTrack->GetDefinition()->GetParticleName() << " of "
710  << aStep->GetPreStepPoint()->GetKineticEnergy() / GeV << " GeV";
711 #endif
712  for (unsigned int i = 0; i < hits.size(); ++i) {
713  G4ThreeVector hitPoint = hits[i].position;
714  if (isItinFidVolume(hitPoint)) {
715  int depth = hits[i].depth;
716  double time = hits[i].time;
717  unsigned int unitID = setDetUnitId(det, hitPoint, depth);
718  currentID.setID(unitID, time, primaryID, 0);
719 #ifdef plotDebug
720  plotProfile(aStep, hitPoint, 1.0 * GeV, time, depth);
721  bool emType = G4TrackToParticleID::isGammaElectronPositron(theTrack->GetDefinition()->GetPDGEncoding());
722  plotHF(hitPoint, emType);
723 #endif
724  processHit(aStep);
725  }
726  }
727 }
728 
729 void HCalSD::hitForFibre(const G4Step* aStep) { // if not ParamShower
730 
731  std::vector<HFShower::Hit> hits = hfshower->getHits(aStep, weight_);
732  if (hits.empty()) {
733  return;
734  }
735 
736  auto const theTrack = aStep->GetTrack();
737  int primaryID = setTrackID(aStep);
738  int det = 5;
739 
741  edepositEM = 1. * GeV;
742  edepositHAD = 0.;
743  } else {
744  edepositEM = 0.;
745  edepositHAD = 1. * GeV;
746  }
747 
748 #ifdef EDM_ML_DEBUG
749  edm::LogVerbatim("HcalSim") << "HCalSD::hitForFibre " << hits.size() << " hits for " << GetName() << " of "
750  << primaryID << " with " << theTrack->GetDefinition()->GetParticleName() << " of "
751  << aStep->GetPreStepPoint()->GetKineticEnergy() / GeV << " GeV in detector type " << det;
752 #endif
753 
754  for (unsigned int i = 0; i < hits.size(); ++i) {
755  G4ThreeVector hitPoint = hits[i].position;
756  if (isItinFidVolume(hitPoint)) {
757  int depth = hits[i].depth;
758  double time = hits[i].time;
759  unsigned int unitID = setDetUnitId(det, hitPoint, depth);
760  currentID.setID(unitID, time, primaryID, 0);
761 #ifdef plotDebug
762  plotProfile(aStep, hitPoint, edepositEM, time, depth);
763  bool emType = (edepositEM > 0.) ? true : false;
764  plotHF(hitPoint, emType);
765 #endif
766  processHit(aStep);
767  }
768  }
769 }
770 
771 void HCalSD::getFromParam(const G4Step* aStep, bool& isKilled) {
772  std::vector<HFShowerParam::Hit> hits = showerParam->getHits(aStep, weight_, isKilled);
773  if (!isKilled || hits.empty()) {
774  return;
775  }
776 
777  int primaryID = setTrackID(aStep);
778  int det = 5;
779 
780 #ifdef EDM_ML_DEBUG
781  edm::LogVerbatim("HcalSim") << "HCalSD::getFromParam " << hits.size() << " hits for " << GetName() << " of "
782  << primaryID << " with " << aStep->GetTrack()->GetDefinition()->GetParticleName()
783  << " of " << aStep->GetPreStepPoint()->GetKineticEnergy() / GeV
784  << " GeV in detector type " << det;
785 #endif
786  for (unsigned int i = 0; i < hits.size(); ++i) {
787  G4ThreeVector hitPoint = hits[i].position;
788  int depth = hits[i].depth;
789  double time = hits[i].time;
790  unsigned int unitID = setDetUnitId(det, hitPoint, depth);
791  currentID.setID(unitID, time, primaryID, 0);
792  edepositEM = hits[i].edep * GeV;
793  edepositHAD = 0.;
794 #ifdef plotDebug
795  plotProfile(aStep, hitPoint, edepositEM, time, depth);
796 #endif
797  processHit(aStep);
798  }
799 }
800 
801 void HCalSD::getHitPMT(const G4Step* aStep) {
802  auto const preStepPoint = aStep->GetPreStepPoint();
803  auto const theTrack = aStep->GetTrack();
804  double edep = showerPMT->getHits(aStep);
805 
806  if (edep >= 0.) {
807  edep *= GeV;
808  double etrack = preStepPoint->GetKineticEnergy();
809  int primaryID = 0;
810  if (etrack >= energyCut) {
811  primaryID = theTrack->GetTrackID();
812  } else {
813  primaryID = theTrack->GetParentID();
814  if (primaryID == 0)
815  primaryID = theTrack->GetTrackID();
816  }
817  // Reset entry point for new primary
818  resetForNewPrimary(aStep);
819  //
820  int det = static_cast<int>(HcalForward);
821  const G4ThreeVector& hitPoint = preStepPoint->GetPosition();
822  double rr = (hitPoint.x() * hitPoint.x() + hitPoint.y() * hitPoint.y());
823  double phi = (rr == 0. ? 0. : atan2(hitPoint.y(), hitPoint.x()));
824  double etaR = showerPMT->getRadius();
825  int depth = 1;
826  if (etaR < 0) {
827  depth = 2;
828  etaR = -etaR;
829  }
830  if (hitPoint.z() < 0)
831  etaR = -etaR;
832 #ifdef EDM_ML_DEBUG
833  edm::LogVerbatim("HcalSim") << "HCalSD::Hit for Detector " << det << " etaR " << etaR << " phi " << phi / deg
834  << " depth " << depth;
835 #endif
836  double time = (aStep->GetPostStepPoint()->GetGlobalTime());
837  uint32_t unitID = 0;
838  if (numberingFromDDD) {
839  HcalNumberingFromDDD::HcalID tmp = numberingFromDDD->unitID(det, etaR, phi, depth, 1);
840  unitID = setDetUnitId(tmp);
841  }
842  currentID.setID(unitID, time, primaryID, 1);
843 
844  edepositHAD = aStep->GetTotalEnergyDeposit();
845  edepositEM = -edepositHAD + edep;
846 #ifdef plotDebug
847  plotProfile(aStep, hitPoint, edep, time, depth);
848 #endif
849 #ifdef EDM_ML_DEBUG
850  double beta = preStepPoint->GetBeta();
851  edm::LogVerbatim("HcalSim") << "HCalSD::getHitPMT 1 hit for " << GetName() << " of " << primaryID << " with "
852  << theTrack->GetDefinition()->GetParticleName() << " of "
853  << preStepPoint->GetKineticEnergy() / GeV << " GeV with velocity " << beta << " UnitID "
854  << std::hex << unitID << std::dec;
855 #endif
856  processHit(aStep);
857  }
858 }
859 
860 void HCalSD::getHitFibreBundle(const G4Step* aStep, bool type) {
861  auto const preStepPoint = aStep->GetPreStepPoint();
862  auto const theTrack = aStep->GetTrack();
863  double edep = showerBundle->getHits(aStep, type);
864 
865  if (edep >= 0.0) {
866  edep *= GeV;
867  double etrack = preStepPoint->GetKineticEnergy();
868  int primaryID = 0;
869  if (etrack >= energyCut) {
870  primaryID = theTrack->GetTrackID();
871  } else {
872  primaryID = theTrack->GetParentID();
873  if (primaryID == 0)
874  primaryID = theTrack->GetTrackID();
875  }
876  // Reset entry point for new primary
877  resetForNewPrimary(aStep);
878  //
879  int det = static_cast<int>(HcalForward);
880  const G4ThreeVector& hitPoint = preStepPoint->GetPosition();
881  double rr = hitPoint.x() * hitPoint.x() + hitPoint.y() * hitPoint.y();
882  double phi = rr == 0. ? 0. : atan2(hitPoint.y(), hitPoint.x());
883  double etaR = showerBundle->getRadius();
884  int depth = 1;
885  if (etaR < 0.) {
886  depth = 2;
887  etaR = -etaR;
888  }
889  if (hitPoint.z() < 0.)
890  etaR = -etaR;
891 #ifdef EDM_ML_DEBUG
892  edm::LogVerbatim("HcalSim") << "HCalSD::Hit for Detector " << det << " etaR " << etaR << " phi " << phi / deg
893  << " depth " << depth;
894 #endif
895  double time = (aStep->GetPostStepPoint()->GetGlobalTime());
896  uint32_t unitID = 0;
897  if (numberingFromDDD) {
898  HcalNumberingFromDDD::HcalID tmp = numberingFromDDD->unitID(det, etaR, phi, depth, 1);
899  unitID = setDetUnitId(tmp);
900  }
901  if (type)
902  currentID.setID(unitID, time, primaryID, 3);
903  else
904  currentID.setID(unitID, time, primaryID, 2);
905 
906  edepositHAD = aStep->GetTotalEnergyDeposit();
907  edepositEM = -edepositHAD + edep;
908 #ifdef plotDebug
909  plotProfile(aStep, hitPoint, edep, time, depth);
910 #endif
911 #ifdef EDM_ML_DEBUG
912  double beta = preStepPoint->GetBeta();
913  edm::LogVerbatim("HcalSim") << "HCalSD::getHitFibreBundle 1 hit for " << GetName() << " of " << primaryID
914  << " with " << theTrack->GetDefinition()->GetParticleName() << " of "
915  << preStepPoint->GetKineticEnergy() / GeV << " GeV with velocity " << beta << " UnitID "
916  << std::hex << unitID << std::dec;
917 #endif
918  processHit(aStep);
919  } // non-zero energy deposit
920 }
921 
923  std::ifstream infile;
924  int entry = 0;
925  infile.open(fName.c_str(), std::ios::in);
926  if (infile) {
927  int det, zside, etaR, phi, lay;
928  double wt;
929  while (infile >> det >> zside >> etaR >> phi >> lay >> wt) {
930  uint32_t id = HcalTestNumbering::packHcalIndex(det, zside, 1, etaR, phi, lay);
931  layerWeights.insert(std::pair<uint32_t, double>(id, wt));
932  ++entry;
933 #ifdef EDM_ML_DEBUG
934  edm::LogVerbatim("HcalSim") << "HCalSD::readWeightFromFile:Entry " << entry << " ID " << std::hex << id
935  << std::dec << " (" << det << "/" << zside << "/1/" << etaR << "/" << phi << "/"
936  << lay << ") Weight " << wt;
937 #endif
938  }
939  infile.close();
940  }
941  edm::LogVerbatim("HcalSim") << "HCalSD::readWeightFromFile: reads " << entry << " weights from " << fName;
942  if (entry <= 0)
943  useLayerWt = false;
944 }
945 
946 double HCalSD::layerWeight(int det, const G4ThreeVector& pos, int depth, int lay) {
947  double wt = 1.;
948  if (numberingFromDDD) {
949  //get the ID's as eta, phi, depth, ... indices
951  numberingFromDDD->unitID(det, math::XYZVectorD(pos.x(), pos.y(), pos.z()), depth, lay);
952  modifyDepth(tmp);
953  uint32_t id = HcalTestNumbering::packHcalIndex(tmp.subdet, tmp.zside, 1, tmp.etaR, tmp.phis, tmp.lay);
954  std::map<uint32_t, double>::const_iterator ite = layerWeights.find(id);
955  if (ite != layerWeights.end())
956  wt = ite->second;
957 #ifdef EDM_ML_DEBUG
958  edm::LogVerbatim("HcalSim") << "HCalSD::layerWeight: ID " << std::hex << id << std::dec << " (" << tmp.subdet << "/"
959  << tmp.zside << "/1/" << tmp.etaR << "/" << tmp.phis << "/" << tmp.lay << ") Weight "
960  << wt;
961 #endif
962  }
963  return wt;
964 }
965 
966 void HCalSD::plotProfile(const G4Step* aStep, const G4ThreeVector& global, double edep, double time, int id) {
967  const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
968  static const unsigned int names = 10;
969  static const G4String modName[names] = {
970  "HEModule", "HVQF", "HBModule", "MBAT", "MBBT", "MBBTC", "MBBT_R1P", "MBBT_R1M", "MBBT_R1PX", "MBBT_R1MX"};
971  G4ThreeVector local;
972  bool found = false;
973  double depth = -2000;
974  int idx = 4;
975  for (int n = 0; n < touch->GetHistoryDepth(); ++n) {
976  G4String name(static_cast<std::string>(dd4hep::dd::noNamespace(touch->GetVolume(n)->GetName())));
977 #ifdef EDM_ML_DEBUG
978  edm::LogVerbatim("HcalSim") << "plotProfile Depth " << n << " Name " << name;
979 #endif
980  for (unsigned int ii = 0; ii < names; ++ii) {
981  if (name == modName[ii]) {
982  found = true;
983  int dn = touch->GetHistoryDepth() - n;
984  local = touch->GetHistory()->GetTransform(dn).TransformPoint(global);
985  if (ii == 0) {
986  depth = local.z() - 4006.5;
987  idx = 1;
988  } else if (ii == 1) {
989  depth = local.z() + 825.0;
990  idx = 3;
991  } else if (ii == 2) {
992  depth = local.x() - 1775.;
993  idx = 0;
994  } else {
995  depth = local.y() + 15.;
996  idx = 2;
997  }
998  break;
999  }
1000  }
1001  if (found)
1002  break;
1003  }
1004  if (!found)
1005  depth = std::abs(global.z()) - 11500;
1006 #ifdef EDM_ML_DEBUG
1007  edm::LogVerbatim("HcalSim") << "plotProfile Found " << found << " Global " << global << " Local " << local
1008  << " depth " << depth << " ID " << id << " EDEP " << edep << " Time " << time;
1009 #endif
1010  if (hit_[idx] != nullptr)
1011  hit_[idx]->Fill(edep);
1012  if (time_[idx] != nullptr)
1013  time_[idx]->Fill(time, edep);
1014  if (dist_[idx] != nullptr)
1015  dist_[idx]->Fill(depth, edep);
1016  int jd = 2 * idx + id - 7;
1017  if (jd >= 0 && jd < 4) {
1018  jd += 5;
1019  if (hit_[jd] != nullptr)
1020  hit_[jd]->Fill(edep);
1021  if (time_[jd] != nullptr)
1022  time_[jd]->Fill(time, edep);
1023  if (dist_[jd] != nullptr)
1024  dist_[jd]->Fill(depth, edep);
1025  }
1026 }
1027 
1028 void HCalSD::plotHF(const G4ThreeVector& hitPoint, bool emType) {
1029  double zv = std::abs(hitPoint.z()) - gpar[4];
1030  if (emType) {
1031  if (hzvem != nullptr)
1032  hzvem->Fill(zv);
1033  } else {
1034  if (hzvhad != nullptr)
1035  hzvhad->Fill(zv);
1036  }
1037 }
1038 
1040  if (id.subdet == 4) {
1041  int ieta = (id.zside == 0) ? -id.etaR : id.etaR;
1042  if (hcalConstants_->maxHFDepth(ieta, id.phis) > 2) {
1043  if (id.depth <= 2) {
1044  if (G4UniformRand() > 0.5)
1045  id.depth += 2;
1046  }
1047  }
1048  } else if ((id.subdet == 1 || id.subdet == 2) && testNumber) {
1049  id.depth = (depth_ == 0) ? 1 : 2;
1050  }
1051 }
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:114
void readWeightFromFile(const std::string &)
Definition: HCalSD.cc:922
std::vector< const G4LogicalVolume * > hfLV
Definition: HCalSD.h:116
bool useParam
Definition: HCalSD.h:105
double eminHitHE
Definition: HCalSD.h:106
TH1F * time_[9]
Definition: HCalSD.h:118
const double GeV
Definition: MathUtil.h:16
double kmaxNeutron
Definition: CaloSD.h:149
void hitForFibre(const G4Step *step)
Definition: HCalSD.cc:729
void fillLogVolumeVector(const std::string &, const std::vector< std::string > &, std::vector< const G4LogicalVolume * > &)
Definition: HCalSD.cc:309
constexpr int ietaAbs() const
get the absolute value of the cell ieta
Definition: HcalDetId.h:148
static bool isMuon(int pdgCode)
uint16_t *__restrict__ id
constexpr int zside() const
get the z-side of the cell (1/-1)
Definition: HcalDetId.h:141
bool useLayerWt
Definition: HCalSD.h:102
Definition: CaloSD.h:40
std::vector< std::string > fibreNames
Definition: HCalSD.h:113
std::unique_ptr< HFShowerParam > showerParam
Definition: HCalSD.h:89
double weight_
Definition: HCalSD.h:108
std::vector< double > gpar
Definition: HCalSD.h:110
bool useFibreBundle
Definition: HCalSD.h:102
const HcalDDDSimConstants * hcalConstants_
Definition: HCalSD.h:93
double betaThr
Definition: HCalSD.h:104
std::vector< std::string > hfNames_
double deliveredLumi
Definition: HCalSD.h:107
float *__restrict__ zv
std::vector< std::string > hfFibreNames_
double eminHitHB
Definition: HCalSD.h:106
bool useShowerLibrary
Definition: HCalSD.h:105
std::vector< std::string > hfNames
Definition: HCalSD.h:112
void plotProfile(const G4Step *step, const G4ThreeVector &pos, double edep, double time, int id)
Definition: HCalSD.cc:966
double birk2
Definition: HCalSD.h:104
bool usePMTHit
Definition: HCalSD.h:102
std::vector< const G4LogicalVolume * > fibre2LV
Definition: HCalSD.h:116
int zside(DetId const &)
static const unsigned int numberOfZLayers
Definition: HFDarkening.h:24
void modifyDepth(HcalNumberingFromDDD::HcalID &id)
Definition: HCalSD.cc:1039
void update(const BeginOfJob *) override
This routine will be called when the appropriate signal arrives.
Definition: HCalSD.cc:551
const std::vector< std::string_view > logicalNames(const std::string &readoutName) const
int ii
Definition: cuy.py:589
void setNumberingScheme(HcalNumberingScheme *)
Definition: HCalSD.cc:544
void processHit(const G4Step *step)
Definition: CaloSD.h:113
bool agingFlagHE
Definition: HCalSD.h:101
double birk1
Definition: HCalSD.h:104
double kmaxProton
Definition: CaloSD.h:149
std::unique_ptr< HFShowerPMT > showerPMT
Definition: HCalSD.h:90
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:118
const HcalSimulationParameters * hcalsimpar() const
bool isItConicalBundle(const G4LogicalVolume *)
Definition: HCalSD.cc:654
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:106
const double MeV
static uint32_t packHcalIndex(int det, int z, int depth, int eta, int phi, int lay)
tuple dd4hep
Definition: dd4hep_cff.py:3
std::unique_ptr< HcalNumberingFromDDD > numberingFromDDD
Definition: HCalSD.h:85
TH1F * hit_[9]
Definition: HCalSD.h:118
ROOT::Math::DisplacementVector3D< ROOT::Math::Cartesian3D< double > > XYZVectorD
spatial vector with cartesian internal representation
Definition: Vector3D.h:8
int getMaxDepth(const int &type) const
bool agingFlagHB
Definition: HCalSD.h:101
bool testNumber
Definition: HCalSD.h:103
constexpr int iphi() const
get the cell iphi
Definition: HcalDetId.h:157
double kmaxIon
Definition: CaloSD.h:149
bool suppressHeavy
Definition: CaloSD.h:148
std::vector< std::string > hfFibreConicalNames_
bool useBirk
Definition: HCalSD.h:102
float edepositHAD
Definition: CaloSD.h:140
std::unique_ptr< HFDarkening > m_HFDarkening
Definition: HCalSD.h:97
double birk3
Definition: HCalSD.h:104
void resetForNewPrimary(const G4Step *)
Definition: CaloSD.cc:644
const HBHEDarkening * m_HEDarkening
Definition: HCalSD.h:96
bool isAvailable() const
Definition: Service.h:40
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:48
bool applyFidCut
Definition: HCalSD.h:105
TH1F * dist_[9]
Definition: HCalSD.h:118
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
void plotHF(const G4ThreeVector &pos, bool emType)
Definition: HCalSD.cc:1028
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)
T * make(const Args &...args) const
make new ROOT object
int getNumberOfHits()
Definition: CaloSD.cc:452
std::vector< int > hfLevels
Definition: HCalSD.h:111
void getFromParam(const G4Step *step, bool &isKilled)
Definition: HCalSD.cc:771
void getFromHFLibrary(const G4Step *step, bool &isKilled)
Definition: HCalSD.cc:686
uint32_t setDetUnitId(const G4Step *step) override
Definition: HCalSD.cc:532
double getLayer0Wt(const int &det, const int &phi, const int &zside) const
bool isItStraightBundle(const G4LogicalVolume *)
Definition: HCalSD.cc:646
bool testNS_
Definition: HCalSD.h:103
bool getFromLibrary(const G4Step *) override
Definition: HCalSD.cc:334
double eminHitHO
Definition: HCalSD.h:106
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:118
const std::vector< double > & getGparHF() const
std::vector< const G4LogicalVolume * > pmtLV
Definition: HCalSD.h:116
int PMTNumber(const G4ThreeVector &pe_effect)
Definition: DetId.h:17
static const unsigned int lowZLimit
Definition: HFDarkening.h:27
TFileDirectory mkdir(const std::string &dir, const std::string &descr="")
create a new subdirectory
Definition: TFileService.h:69
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:103
bool isHF
Definition: HCalSD.h:100
void getHitFibreBundle(const G4Step *step, bool type)
Definition: HCalSD.cc:860
std::vector< std::string > hcalMaterialNames_
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
double getEnergyDeposit(const G4Step *) override
Definition: HCalSD.cc:398
int ke
bool isPrimary() const
std::unique_ptr< HFShowerFibreBundle > showerBundle
Definition: HCalSD.h:91
std::vector< std::string > hfPMTNames_
bool isItPMT(const G4LogicalVolume *)
Definition: HCalSD.cc:638
bool isItinFidVolume(const G4ThreeVector &)
Definition: HCalSD.cc:670
bool isHF(int etabin, int depth)
void getHitPMT(const G4Step *step)
Definition: HCalSD.cc:801
double getAttenuation(const G4Step *aStep, double birk1, double birk2, double birk3) const
Definition: CaloSD.cc:657
std::unique_ptr< HFShowerLibrary > showerLibrary
Definition: HCalSD.h:87
double layerWeight(int, const G4ThreeVector &, int, int)
Definition: HCalSD.cc:946
const HcalSimulationConstants * hcalSimConstants_
Definition: HCalSD.h:94
int depth_
Definition: HCalSD.h:109
virtual int setTrackID(const G4Step *)
Definition: CaloSD.cc:814
std::map< uint32_t, double > layerWeights
Definition: HCalSD.h:117
bool isItHF(const G4Step *)
Definition: HCalSD.cc:601
static bool isGammaElectronPositron(int pdgCode)
bool isItScintillator(const G4Material *)
Definition: HCalSD.cc:662
std::vector< std::string > hfFibreStraightNames_
list entry
Definition: mps_splice.py:68
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:116
tuple tfile
Definition: compare.py:324
uint32_t getUnitID() const
Definition: CaloG4Hit.h:66
Log< level::Warning, false > LogWarning
double getResponseWt(const G4Track *)
Definition: CaloSD.cc:849
std::vector< const G4LogicalVolume * > fibre1LV
Definition: HCalSD.h:116
tmp
align.sh
Definition: createJobs.py:716
int maxHFDepth(const int &ieta, const int &iphi) const
std::vector< const G4Material * > materials
Definition: HCalSD.h:115
void initRun() override
Definition: HCalSD.cc:553
const HBHEDarkening * m_HBDarkening
Definition: HCalSD.h:95
float degradation(float intlumi, int ieta, int lay) const
#define EDM_ML_DEBUG
double getEnergyDeposit() const
Definition: CaloG4Hit.h:79
bool useHF
Definition: HCalSD.h:105
bool isItFibre(const G4LogicalVolume *)
Definition: HCalSD.cc:622
void setParameterized(bool val)
Definition: CaloSD.h:110