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