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