CMS 3D CMS Logo

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