CMS 3D CMS Logo

ECalSD.cc
Go to the documentation of this file.
1 // File: ECalSD.cc
3 // Description: Sensitive Detector class for electromagnetic calorimeters
17 
21 
22 #include "G4LogicalVolumeStore.hh"
23 #include "G4LogicalVolume.hh"
24 #include "G4Step.hh"
25 #include "G4Track.hh"
26 #include "G4VProcess.hh"
27 #include "G4SystemOfUnits.hh"
28 
29 #include <algorithm>
30 
31 //#define EDM_ML_DEBUG
32 
33 using namespace geant_units::operators;
34 
35 template <class T>
36 bool any(const std::vector<T>& v, const T& what) {
37  return std::find(v.begin(), v.end(), what) != v.end();
38 }
39 
41  const EcalSimulationParameters* ecpar,
42  const SensitiveDetectorCatalog& clg,
43  edm::ParameterSet const& p,
44  const SimTrackManager* manager)
45  : CaloSD(name,
46  clg,
47  p,
48  manager,
49  (float)(p.getParameter<edm::ParameterSet>("ECalSD").getParameter<double>("TimeSliceUnit")),
50  p.getParameter<edm::ParameterSet>("ECalSD").getParameter<bool>("IgnoreTrackID")),
51  ecalSimParameters_(ecpar) {
52  numberingScheme_.reset(nullptr);
53  // static SimpleConfigurable<bool> on1(false, "ECalSD:UseBirkLaw");
54  // static SimpleConfigurable<double> bk1(0.00463,"ECalSD:BirkC1");
55  // static SimpleConfigurable<double> bk2(-0.03, "ECalSD:BirkC2");
56  // static SimpleConfigurable<double> bk3(1.0, "ECalSD:BirkC3");
57  // Values from NIM A484 (2002) 239-244: as implemented in Geant3
58  // useBirk = on1.value();
59  // birk1 = bk1.value()*(g/(MeV*cm2));
60  // birk2 = bk2.value()*(g/(MeV*cm2))*(g/(MeV*cm2));
61  edm::ParameterSet m_EC = p.getParameter<edm::ParameterSet>("ECalSD");
62  useBirk = m_EC.getParameter<bool>("UseBirkLaw");
63  useBirkL3 = m_EC.getParameter<bool>("BirkL3Parametrization");
64  double bunit = (CLHEP::g / (CLHEP::MeV * CLHEP::cm2));
65  birk1 = m_EC.getParameter<double>("BirkC1") * bunit;
66  birk2 = m_EC.getParameter<double>("BirkC2");
67  birk3 = m_EC.getParameter<double>("BirkC3");
68  birkSlope = m_EC.getParameter<double>("BirkSlope");
69  birkCut = m_EC.getParameter<double>("BirkCut");
70  slopeLY = m_EC.getParameter<double>("SlopeLightYield");
71  storeTrack = m_EC.getParameter<bool>("StoreSecondary");
72  crystalMat = m_EC.getUntrackedParameter<std::string>("XtalMat", "E_PbWO4");
73  bool isItTB = m_EC.getUntrackedParameter<bool>("TestBeam", false);
74  bool nullNS = m_EC.getUntrackedParameter<bool>("NullNumbering", false);
75  storeRL = m_EC.getUntrackedParameter<bool>("StoreRadLength", false);
76  scaleRL = m_EC.getUntrackedParameter<double>("ScaleRadLength", 1.0);
77  int dumpGeom = m_EC.getUntrackedParameter<int>("DumpGeometry", 0);
78 
79  //Changes for improved timing simulation
80  storeLayerTimeSim = m_EC.getUntrackedParameter<bool>("StoreLayerTimeSim", false);
81 
82  ageingWithSlopeLY = m_EC.getUntrackedParameter<bool>("AgeingWithSlopeLY", false);
84  ageing.setLumies(p.getParameter<edm::ParameterSet>("ECalSD").getParameter<double>("DelivLuminosity"),
85  p.getParameter<edm::ParameterSet>("ECalSD").getParameter<double>("InstLuminosity"));
86 
87  if (ecalSimParameters_ == nullptr) {
88  edm::LogError("EcalSim") << "ECalSD : Cannot find EcalSimulationParameters for " << name;
89  throw cms::Exception("Unknown", "ECalSD") << "Cannot find EcalSimulationParameters for " << name << "\n";
90  }
91 
92  // Use of Weight
94 #ifdef EDM_ML_DEBUG
95  edm::LogVerbatim("EcalSim") << "ECalSD:: useWeight " << useWeight;
96 #endif
99 #ifdef EDM_ML_DEBUG
100  edm::LogVerbatim("EcalSim") << "Names (Depth 1):" << depth1Name << " (Depth 2):" << depth2Name << std::endl;
101 #endif
102  int type(-1);
103  bool dump(false);
104  EcalNumberingScheme* scheme = nullptr;
105  if (nullNS) {
106  scheme = nullptr;
107  } else if (name == "EcalHitsEB") {
108  scheme = dynamic_cast<EcalNumberingScheme*>(new EcalBarrelNumberingScheme());
109  type = 0;
110  dump = ((dumpGeom % 10) > 0);
111  } else if (name == "EcalHitsEE") {
112  scheme = dynamic_cast<EcalNumberingScheme*>(new EcalEndcapNumberingScheme());
113  type = 1;
114  dump = (((dumpGeom / 10) % 10) > 0);
115  } else if (name == "EcalHitsES") {
116  if (isItTB)
117  scheme = dynamic_cast<EcalNumberingScheme*>(new ESTBNumberingScheme());
118  else
120  useWeight = false;
121  type = 2;
122  dump = (((dumpGeom / 100) % 10) > 0);
123  } else {
124  edm::LogWarning("EcalSim") << "ECalSD: ReadoutName not supported";
125  }
126  int type0 = dumpGeom / 1000;
127  type += (10 * type0);
128 
129  if (nullptr != scheme) {
131  }
132 #ifdef EDM_ML_DEBUG
133  edm::LogVerbatim("EcalSim") << "Constructing a ECalSD with name " << GetName();
134 #endif
135  if (useWeight) {
136  edm::LogVerbatim("EcalSim") << "ECalSD:: Use of Birks law is set to " << useBirk
137  << " with three constants kB = " << birk1 / bunit << ", C1 = " << birk2
138  << ", C2 = " << birk3 << "\n Use of L3 parametrization " << useBirkL3
139  << " with slope " << birkSlope << " and cut off " << birkCut << "\n"
140  << " Slope for Light yield is set to " << slopeLY;
141  } else {
142  edm::LogVerbatim("EcalSim") << "ECalSD:: energy deposit is not corrected "
143  << " by Birk or light yield curve";
144  }
145 
146  edm::LogVerbatim("EcalSim") << "ECalSD:: Suppression Flag " << suppressHeavy << "\tprotons below "
147  << kmaxProton / CLHEP::MeV << " MeV,\tneutrons below " << kmaxNeutron / CLHEP::MeV
148  << " MeV,\tions below " << kmaxIon / CLHEP::MeV << " MeV \n\tDepth1 Name = " << depth1Name
149  << "\tDepth2 Name = " << depth2Name << "\n\tstoreRL " << storeRL << ":" << scaleRL
150  << "\tstoreLayerTimeSim " << storeLayerTimeSim << "\n\ttime Granularity "
151  << p.getParameter<edm::ParameterSet>("ECalSD").getParameter<double>("TimeSliceUnit")
152  << " ns";
153  if (useWeight)
154  initMap();
155 #ifdef plotDebug
157  if (tfile.isAvailable()) {
158  TFileDirectory ecDir = tfile->mkdir("ProfileFromECalSD");
159  static const std::string ctype[4] = {"EB", "EBref", "EE", "EERef"};
160  for (int k = 0; k < 4; ++k) {
161  std::string name = "ECLL_" + ctype[k];
162  std::string title = "Local vs Global for " + ctype[k];
163  double xmin = (k > 1) ? 3000.0 : 1000.0;
164  g2L_[k] = ecDir.make<TH2F>(name.c_str(), title.c_str(), 100, xmin, xmin + 1000., 100, 0.0, 3000.);
165  }
166  } else {
167  for (int k = 0; k < 4; ++k)
168  g2L_[k] = 0;
169  }
170 #endif
171  if (dump) {
172  const auto& lvNames = clg.logicalNames(name);
174  geom.update();
175  }
176 }
177 
179 
180 double ECalSD::getEnergyDeposit(const G4Step* aStep) {
181  const G4StepPoint* preStepPoint = aStep->GetPreStepPoint();
182  const G4Track* theTrack = aStep->GetTrack();
183  double edep = aStep->GetTotalEnergyDeposit();
184 
185  // take into account light collection curve for crystals
186  double weight = 1.;
187  if (suppressHeavy) {
188  TrackInformation* trkInfo = (TrackInformation*)(theTrack->GetUserInformation());
189  if (trkInfo) {
190  int pdg = theTrack->GetDefinition()->GetPDGEncoding();
191  if (!(trkInfo->isPrimary())) { // Only secondary particles
192  double ke = theTrack->GetKineticEnergy();
193  if (((pdg / 1000000000 == 1 && ((pdg / 10000) % 100) > 0 && ((pdg / 10) % 100) > 0)) && (ke < kmaxIon))
194  weight = 0;
195  if ((pdg == 2212) && (ke < kmaxProton))
196  weight = 0;
197  if ((pdg == 2112) && (ke < kmaxNeutron))
198  weight = 0;
199  }
200  }
201  }
202  const G4LogicalVolume* lv = preStepPoint->GetTouchable()->GetVolume(0)->GetLogicalVolume();
203  double wt1 = 1.0;
204  if (useWeight && !any(noWeight, lv)) {
205  weight *= curve_LY(lv);
206  if (useBirk) {
207  if (useBirkL3)
208  weight *= getBirkL3(aStep);
209  else
210  weight *= getAttenuation(aStep, birk1, birk2, birk3);
211  }
212  wt1 = getResponseWt(theTrack);
213  }
214  edep *= weight * wt1;
215  // Russian Roulette
216  double wt2 = theTrack->GetWeight();
217  if (wt2 > 0.0) {
218  edep *= wt2;
219  }
220 #ifdef EDM_ML_DEBUG
221  edm::LogVerbatim("EcalSim") << lv->GetName() << " "
222  << DD4hep2DDDName::noNameSpace(static_cast<std::string>(lv->GetName()))
223  << " Light Collection Efficiency " << weight << ":" << wt1 << " wt2= " << wt2
224  << " Weighted Energy Deposit " << edep / CLHEP::MeV << " MeV at "
225  << preStepPoint->GetPosition();
226 #endif
227  return edep;
228 }
229 
230 double ECalSD::EnergyCorrected(const G4Step& step, const G4Track* track) {
231  const G4StepPoint* hitPoint = step.GetPreStepPoint();
232  const G4LogicalVolume* lv = hitPoint->GetTouchable()->GetVolume(0)->GetLogicalVolume();
233  if (lv->GetSensitiveDetector() != this)
234  return 0.0;
235 
236  double edep = step.GetTotalEnergyDeposit();
237  if (useWeight && !any(noWeight, lv)) {
238  currentLocalPoint = setToLocal(hitPoint->GetPosition(), hitPoint->GetTouchable());
239  auto ite = xtalLMap.find(lv);
240  crystalLength = (ite == xtalLMap.end()) ? 230._mm : std::abs(ite->second);
241  crystalDepth = (ite == xtalLMap.end()) ? 0.0 : (std::abs(0.5 * (ite->second) + currentLocalPoint.z()));
242  edep *= curve_LY(lv) * getResponseWt(track);
243  }
244  return edep;
245 }
246 
247 int ECalSD::getTrackID(const G4Track* aTrack) {
248  int primaryID(0);
249  if (storeTrack && depth > 0) {
250  forceSave = true;
251  primaryID = aTrack->GetTrackID();
252  } else {
253  primaryID = CaloSD::getTrackID(aTrack);
254  }
255  return primaryID;
256 }
257 
258 uint16_t ECalSD::getDepth(const G4Step* aStep) {
259  // this method should be called first at a step
260  const G4StepPoint* hitPoint = aStep->GetPreStepPoint();
261  currentLocalPoint = setToLocal(hitPoint->GetPosition(), hitPoint->GetTouchable());
262  const G4LogicalVolume* lv = hitPoint->GetTouchable()->GetVolume(0)->GetLogicalVolume();
263 
264  auto ite = xtalLMap.find(lv);
265  crystalLength = (ite == xtalLMap.end()) ? 230._mm : std::abs(ite->second);
266  crystalDepth = (ite == xtalLMap.end()) ? 0.0 : (std::abs(0.5 * (ite->second) + currentLocalPoint.z()));
267  depth = any(useDepth1, lv) ? 1 : (any(useDepth2, lv) ? 2 : 0);
268  uint16_t depth1(0), depth2(0);
269  if (storeRL) {
270  depth1 = (ite == xtalLMap.end()) ? 0 : (((ite->second) >= 0) ? 0 : PCaloHit::kEcalDepthRefz);
271  depth2 = getRadiationLength(hitPoint, lv);
272  depth |= (((depth2 & PCaloHit::kEcalDepthMask) << PCaloHit::kEcalDepthOffset) | depth1);
273  } else if (storeLayerTimeSim) {
274  depth2 = getLayerIDForTimeSim();
276  }
277 #ifdef EDM_ML_DEBUG
278  if (isXtal(lv))
279  edm::LogVerbatim("EcalSimX") << "ECalSD::Volume " << lv->GetName() << " DetId " << std::hex << setDetUnitId(aStep)
280  << std::dec << " Global " << (hitPoint->GetPosition()).rho() << ":"
281  << (hitPoint->GetPosition()).z() << " Local Z " << currentLocalPoint.z() << " Depth "
282  << crystalDepth;
283  edm::LogVerbatim("EcalSim") << "ECalSD::Depth " << std::hex << depth1 << ":" << depth2 << ":" << depth << std::dec
284  << " L " << (ite == xtalLMap.end()) << ":" << ite->second << " local "
285  << currentLocalPoint << " Crystal length " << crystalLength << ":" << crystalDepth;
286 #endif
287  return depth;
288 }
289 
290 uint16_t ECalSD::getRadiationLength(const G4StepPoint* hitPoint, const G4LogicalVolume* lv) {
291  uint16_t thisX0 = 0;
292  if (useWeight) {
293  double radl = hitPoint->GetMaterial()->GetRadlen();
294  thisX0 = (uint16_t)floor(scaleRL * crystalDepth / radl);
295 #ifdef plotDebug
296  const std::string lvname = DD4hep2DDDName::noNameSpace(static_cast<std::string>(lv->GetName()));
297  int k1 = (lvname.find("EFRY") != std::string::npos) ? 2 : 0;
298  int k2 = (lvname.find("refl") != std::string::npos) ? 1 : 0;
299  int kk = k1 + k2;
300  double rz = (k1 == 0) ? (hitPoint->GetPosition()).rho() : std::abs((hitPoint->GetPosition()).z());
301  edm::LogVerbatim("EcalSim") << lvname << " # " << k1 << ":" << k2 << ":" << kk << " rz " << rz << " D " << thisX0;
302  g2L_[kk]->Fill(rz, thisX0);
303 #endif
304 #ifdef EDM_ML_DEBUG
305  G4ThreeVector localPoint = setToLocal(hitPoint->GetPosition(), hitPoint->GetTouchable());
306  edm::LogVerbatim("EcalSim") << lv->GetName() << " "
307  << DD4hep2DDDName::noNameSpace(static_cast<std::string>(lv->GetName())) << " Global "
308  << hitPoint->GetPosition() << ":" << (hitPoint->GetPosition()).rho() << " Local "
309  << localPoint << " Crystal Length " << crystalLength << " Radl " << radl
310  << " crystalDepth " << crystalDepth << " Index " << thisX0 << " : "
312 #endif
313  }
314  return thisX0;
315 }
316 
318  const double invLayerSize = 0.1; //layer size in 1/mm
319  return (int)crystalDepth * invLayerSize;
320 }
321 
322 uint32_t ECalSD::setDetUnitId(const G4Step* aStep) {
323  if (numberingScheme_.get() == nullptr) {
324  return EBDetId(1, 1)();
325  } else {
326  getBaseNumber(aStep);
327  return numberingScheme_->getUnitID(theBaseNumber);
328  }
329 }
330 
332  if (scheme != nullptr) {
333  edm::LogVerbatim("EcalSim") << "EcalSD: updates numbering scheme for " << GetName();
334  numberingScheme_.reset(scheme);
335  }
336 }
337 
339  std::vector<const G4LogicalVolume*> lvused;
340  const G4LogicalVolumeStore* lvs = G4LogicalVolumeStore::GetInstance();
341  std::map<const std::string, const G4LogicalVolume*> nameMap;
342  for (auto lvi = lvs->begin(), lve = lvs->end(); lvi != lve; ++lvi)
343  nameMap.emplace(DD4hep2DDDName::noNameSpace(static_cast<std::string>((*lvi)->GetName())), *lvi);
344 
345  for (unsigned int it = 0; it < ecalSimParameters_->lvNames_.size(); ++it) {
346  const std::string& matname = ecalSimParameters_->matNames_[it];
347  const std::string& lvname = ecalSimParameters_->lvNames_[it];
348  const G4LogicalVolume* lv = nameMap[lvname];
349  int ibec = (lvname.find("EFRY") == std::string::npos) ? 0 : 1;
350  int iref = (lvname.find("refl") == std::string::npos) ? 0 : 1;
351  int type = (ibec + iref == 1) ? 1 : -1;
352  if (depth1Name != " ") {
353  if (strncmp(lvname.c_str(), depth1Name.c_str(), 4) == 0) {
354  if (!any(useDepth1, lv)) {
355  useDepth1.push_back(lv);
356 #ifdef EDM_ML_DEBUG
357  edm::LogVerbatim("EcalSim") << "ECalSD::initMap Logical Volume " << lvname << " in Depth 1 volume list";
358 #endif
359  }
360  const G4LogicalVolume* lvr = nameMap[lvname + "_refl"];
361  if (lvr != nullptr && !any(useDepth1, lvr)) {
362  useDepth1.push_back(lvr);
363 #ifdef EDM_ML_DEBUG
364  edm::LogVerbatim("EcalSim") << "ECalSD::initMap Logical Volume " << lvname << "_refl"
365  << " in Depth 1 volume list";
366 #endif
367  }
368  }
369  }
370  if (depth2Name != " ") {
371  if (strncmp(lvname.c_str(), depth2Name.c_str(), 4) == 0) {
372  if (!any(useDepth2, lv)) {
373  useDepth2.push_back(lv);
374 #ifdef EDM_ML_DEBUG
375  edm::LogVerbatim("EcalSim") << "ECalSD::initMap Logical Volume " << lvname << " in Depth 2 volume list";
376 #endif
377  }
378  const G4LogicalVolume* lvr = nameMap[lvname + "_refl"];
379  if (lvr != nullptr && !any(useDepth2, lvr)) {
380  useDepth2.push_back(lvr);
381 #ifdef EDM_ML_DEBUG
382  edm::LogVerbatim("EcalSim") << "ECalSD::initMap Logical Volume " << lvname << "_refl"
383  << " in Depth 2 volume list";
384 #endif
385  }
386  }
387  }
388  if (lv != nullptr) {
389  if (crystalMat.size() == matname.size() && !strcmp(crystalMat.c_str(), matname.c_str())) {
390  if (!any(lvused, lv)) {
391  lvused.push_back(lv);
392  double dz = ecalSimParameters_->dzs_[it];
393  xtalLMap.insert(std::pair<const G4LogicalVolume*, double>(lv, dz * type));
394  lv = nameMap[lvname + "_refl"];
395  if (lv != nullptr) {
396  xtalLMap.insert(std::pair<const G4LogicalVolume*, double>(lv, -dz * type));
397  }
398  }
399  } else {
400  if (!any(noWeight, lv)) {
401  noWeight.push_back(lv);
402 #ifdef EDM_ML_DEBUG
403  edm::LogVerbatim("EcalSim") << "ECalSD::initMap Logical Volume " << lvname << " Material " << matname
404  << " in noWeight list";
405 #endif
406  }
407  lv = nameMap[lvname];
408  if (lv != nullptr && !any(noWeight, lv)) {
409  noWeight.push_back(lv);
410 #ifdef EDM_ML_DEBUG
411  edm::LogVerbatim("EcalSim") << "ECalSD::initMap Logical Volume " << lvname << " Material " << matname
412  << " in noWeight list";
413 #endif
414  }
415  }
416  }
417  }
418 #ifdef EDM_ML_DEBUG
419  edm::LogVerbatim("EcalSim") << "ECalSD: Length Table:";
420  int i = 0;
421  for (auto ite : xtalLMap) {
422  std::string name("Unknown");
423  if (ite.first != nullptr)
424  name = DD4hep2DDDName::noNameSpace(static_cast<std::string>((ite.first)->GetName()));
425  edm::LogVerbatim("EcalSim") << " " << i << " " << ite.first << " " << name << " L = " << ite.second;
426  ++i;
427  }
428 #endif
429 }
430 
431 double ECalSD::curve_LY(const G4LogicalVolume* lv) {
432  double weight = 1.;
433  if (ageingWithSlopeLY) {
434  //position along the crystal in mm from 0 to 230 (in EB)
435  if (crystalDepth >= -0.1 || crystalDepth <= crystalLength + 0.1)
437  } else {
438  double dapd = crystalLength - crystalDepth;
439  if (dapd >= -0.1 || dapd <= crystalLength + 0.1) {
440  if (dapd <= 100.)
441  weight = 1.0 + slopeLY - dapd * 0.01 * slopeLY;
442  } else {
443  edm::LogWarning("EcalSim") << "ECalSD: light coll curve : wrong distance "
444  << "to APD " << dapd << " crlength = " << crystalLength << ":" << crystalDepth
445  << " crystal name = " << lv->GetName() << " "
446  << DD4hep2DDDName::noNameSpace(static_cast<std::string>(lv->GetName()))
447  << " z of localPoint = " << currentLocalPoint.z() << " take weight = " << weight;
448  }
449  }
450 #ifdef EDM_ML_DEBUG
451  edm::LogVerbatim("EcalSim") << "ECalSD: light coll curve : crlength = " << crystalLength << " Depth " << crystalDepth
452  << " crystal name = " << lv->GetName() << " "
453  << DD4hep2DDDName::noNameSpace(static_cast<std::string>(lv->GetName()))
454  << " z of localPoint = " << currentLocalPoint.z() << " take weight = " << weight;
455 #endif
456  return weight;
457 }
458 
459 void ECalSD::getBaseNumber(const G4Step* aStep) {
461  const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
462  int theSize = touch->GetHistoryDepth() + 1;
463  if (theBaseNumber.getCapacity() < theSize)
464  theBaseNumber.setSize(theSize);
465  //Get name and copy numbers
466  if (theSize > 1) {
467  for (int ii = 0; ii < theSize; ii++) {
468  std::string name = DD4hep2DDDName::noNameSpace(static_cast<std::string>(touch->GetVolume(ii)->GetName()));
469  theBaseNumber.addLevel(name, touch->GetReplicaNumber(ii));
470 #ifdef EDM_ML_DEBUG
471  edm::LogVerbatim("EcalSim") << "ECalSD::getBaseNumber(): Adding level " << ii << ": " << name << "["
472  << touch->GetReplicaNumber(ii) << "]";
473 #endif
474  }
475  }
476 }
477 
478 double ECalSD::getBirkL3(const G4Step* aStep) {
479  double weight = 1.;
480  const G4StepPoint* preStepPoint = aStep->GetPreStepPoint();
481  double charge = preStepPoint->GetCharge();
482 
483  if (charge != 0. && aStep->GetStepLength() > 0.) {
484  const G4Material* mat = preStepPoint->GetMaterial();
485  double density = mat->GetDensity();
486  double dedx = aStep->GetTotalEnergyDeposit() / aStep->GetStepLength();
487  double rkb = birk1 / density;
488  if (dedx > 0) {
489  weight = 1. - birkSlope * log(rkb * dedx);
490  if (weight < birkCut)
491  weight = birkCut;
492  else if (weight > 1.)
493  weight = 1.;
494  }
495 #ifdef EDM_ML_DEBUG
496  edm::LogVerbatim("EcalSim") << "ECalSD::getBirkL3 in "
497  << DD4hep2DDDName::noNameSpace(static_cast<std::string>(mat->GetName())) << " Charge "
498  << charge << " dE/dx " << dedx << " Birk Const " << rkb << " Weight = " << weight
499  << " dE " << aStep->GetTotalEnergyDeposit();
500 #endif
501  }
502  return weight;
503 }
504 
505 bool ECalSD::isXtal(const G4LogicalVolume* lv) { return (xtalLMap.find(lv) != xtalLMap.end()); }
Log< level::Info, true > LogVerbatim
double crystalLength
Definition: ECalSD.h:74
bool isPrimary() const
bool isXtal(const G4LogicalVolume *)
Definition: ECalSD.cc:505
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
bool useBirkL3
Definition: ECalSD.h:62
double kmaxNeutron
Definition: CaloSD.h:155
std::vector< std::string > lvNames_
double birkSlope
Definition: ECalSD.h:63
double slopeLY
Definition: ECalSD.h:64
void setLumies(double x, double y)
Definition: CaloSD.h:40
std::vector< const G4LogicalVolume * > useDepth2
Definition: ECalSD.h:67
void getBaseNumber(const G4Step *)
Definition: ECalSD.cc:459
static const int kEcalDepthRefz
Definition: PCaloHit.h:63
std::vector< const G4LogicalVolume * > noWeight
Definition: ECalSD.h:67
bool any(const std::vector< T > &v, const T &what)
Definition: ECalSD.cc:36
double birk1
Definition: ECalSD.h:63
double getEnergyDeposit(const G4Step *) override
Definition: ECalSD.cc:180
bool useWeight
Definition: ECalSD.h:61
G4ThreeVector setToLocal(const G4ThreeVector &, const G4VTouchable *) const
Definition: CaloSD.cc:454
Definition: weight.py:1
std::unique_ptr< EcalNumberingScheme > numberingScheme_
Definition: ECalSD.h:60
virtual int getTrackID(const G4Track *)
Definition: CaloSD.cc:873
double crystalDepth
Definition: ECalSD.h:75
Log< level::Error, false > LogError
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
double kmaxProton
Definition: CaloSD.h:155
const std::vector< std::string_view > logicalNames(const std::string &readoutName) const
bool storeLayerTimeSim
Definition: ECalSD.h:61
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
bool storeRL
Definition: ECalSD.h:61
~ECalSD() override
Definition: ECalSD.cc:178
double getBirkL3(const G4Step *)
Definition: ECalSD.cc:478
void initMap()
Definition: ECalSD.cc:338
T getUntrackedParameter(std::string const &, T const &) const
int getTrackID(const G4Track *) override
Definition: ECalSD.cc:247
std::vector< std::string > matNames_
bool forceSave
Definition: CaloSD.h:157
T * make(const Args &...args) const
make new ROOT object
std::string depth1Name
Definition: ECalSD.h:65
void addLevel(const std::string &name, const int &copyNumber)
std::map< const G4LogicalVolume *, double > xtalLMap
Definition: ECalSD.h:66
double getResponseWt(const G4Track *, int k=0)
Definition: CaloSD.cc:928
EcalBaseNumber theBaseNumber
Definition: ECalSD.h:68
bool storeTrack
Definition: ECalSD.h:61
double kmaxIon
Definition: CaloSD.h:155
bool suppressHeavy
Definition: CaloSD.h:154
Definition: tfile.py:1
static const int kEcalDepthMask
Definition: PCaloHit.h:61
double scaleRL
Definition: ECalSD.h:64
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
uint16_t getLayerIDForTimeSim()
Definition: ECalSD.cc:317
static const int kEcalDepthOffset
Definition: PCaloHit.h:62
const EcalSimulationParameters * ecalSimParameters_
Definition: ECalSD.h:59
ii
Definition: cuy.py:589
std::string depth2Name
Definition: ECalSD.h:65
G4ThreeVector currentLocalPoint
Definition: ECalSD.h:73
uint16_t getRadiationLength(const G4StepPoint *hitPoint, const G4LogicalVolume *lv)
Definition: ECalSD.cc:290
std::vector< const G4LogicalVolume * > useDepth1
Definition: ECalSD.h:67
bool ageingWithSlopeLY
Definition: ECalSD.h:70
bool useBirk
Definition: ECalSD.h:62
double getAttenuation(const G4Step *aStep, double birk1, double birk2, double birk3) const
Definition: CaloSD.cc:727
double calcLightCollectionEfficiencyWeighted(DetId id, double z)
HLT enums.
ECalSD(const std::string &, const EcalSimulationParameters *, const SensitiveDetectorCatalog &, edm::ParameterSet const &p, const SimTrackManager *)
Definition: ECalSD.cc:40
double birk2
Definition: ECalSD.h:63
uint16_t depth
Definition: ECalSD.h:76
CaloHitID currentID[2]
Definition: CaloSD.h:146
std::string noNameSpace(const std::string &name)
double curve_LY(const G4LogicalVolume *)
Definition: ECalSD.cc:431
step
Definition: StallMonitor.cc:83
EnergyResolutionVsLumi ageing
Definition: ECalSD.h:69
Log< level::Warning, false > LogWarning
void setNumberingScheme(EcalNumberingScheme *)
Definition: ECalSD.cc:331
uint32_t setDetUnitId(const G4Step *) override
Definition: ECalSD.cc:322
void setSize(const int &size)
long double T
double EnergyCorrected(const G4Step &, const G4Track *) override
Definition: ECalSD.cc:230
double birkCut
Definition: ECalSD.h:63
std::string crystalMat
Definition: ECalSD.h:65
double birk3
Definition: ECalSD.h:63
uint16_t getDepth(const G4Step *) override
Definition: ECalSD.cc:258