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