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 "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  ageingWithSlopeLY = m_EC.getUntrackedParameter<bool>("AgeingWithSlopeLY", false);
71  if(ageingWithSlopeLY) ageing.setLumies(p.getParameter<edm::ParameterSet>("ECalSD").getParameter<double>("DelivLuminosity"),
72  p.getParameter<edm::ParameterSet>("ECalSD").getParameter<double>("InstLuminosity"));
73 
74  //Material list for HB/HE/HO sensitive detectors
75  std::string attribute = "ReadOutName";
77  DDValue ddv(attribute,name,0);
78  filter.setCriteria(ddv,DDCompOp::equals);
79  DDFilteredView fv(cpv);
80  fv.addFilter(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") scheme = dynamic_cast<EcalNumberingScheme*>(new EcalBarrelNumberingScheme());
97  else if (name == "EcalHitsEE") scheme = dynamic_cast<EcalNumberingScheme*>(new EcalEndcapNumberingScheme());
98  else if (name == "EcalHitsES") {
99  if (isItTB) scheme = dynamic_cast<EcalNumberingScheme*>(new ESTBNumberingScheme());
100  else scheme = dynamic_cast<EcalNumberingScheme*>(new EcalPreshowerNumberingScheme());
101  useWeight = false;
102  } else {edm::LogWarning("EcalSim") << "ECalSD: ReadoutName not supported\n";}
103 
104  if (scheme) setNumberingScheme(scheme);
105 #ifdef DebugLog
106  LogDebug("EcalSim") << "Constructing a ECalSD with name " << GetName();
107 #endif
108  if (useWeight) {
109  edm::LogInfo("EcalSim") << "ECalSD:: Use of Birks law is set to "
110  << useBirk << " with three constants kB = "
111  << birk1 << ", C1 = " << birk2 << ", C2 = "
112  << birk3 <<"\n Use of L3 parametrization "
113  << useBirkL3 << " with slope " << birkSlope
114  << " and cut off " << birkCut << "\n"
115  << " Slope for Light yield is set to "
116  << slopeLY;
117  } else {
118  edm::LogInfo("EcalSim") << "ECalSD:: energy deposit is not corrected "
119  << " by Birk or light yield curve";
120  }
121 
122  edm::LogInfo("EcalSim") << "ECalSD:: Suppression Flag " << suppressHeavy
123  << " protons below " << kmaxProton << " MeV,"
124  << " neutrons below " << kmaxNeutron << " MeV and"
125  << " ions below " << kmaxIon << " MeV\n"
126  << " Depth1 Name = " << depth1Name
127  << " and Depth2 Name = " << depth2Name;
128 
129  if (useWeight) initMap(name,cpv);
130 
131 }
132 
134  if (numberingScheme) delete numberingScheme;
135 }
136 
137 double ECalSD::getEnergyDeposit(G4Step * aStep) {
138 
139  if (aStep == NULL) {
140  return 0;
141  } else {
142  preStepPoint = aStep->GetPreStepPoint();
143  G4Track* theTrack = aStep->GetTrack();
144  double wt2 = theTrack->GetWeight();
145  G4String nameVolume = preStepPoint->GetPhysicalVolume()->GetName();
146 
147  // take into account light collection curve for crystals
148  double weight = 1.;
149  if (suppressHeavy) {
150  TrackInformation * trkInfo = (TrackInformation *)(theTrack->GetUserInformation());
151  if (trkInfo) {
152  int pdg = theTrack->GetDefinition()->GetPDGEncoding();
153  if (!(trkInfo->isPrimary())) { // Only secondary particles
154  double ke = theTrack->GetKineticEnergy()/MeV;
155  if (((pdg/1000000000 == 1 && ((pdg/10000)%100) > 0 &&
156  ((pdg/10)%100) > 0)) && (ke<kmaxIon)) weight = 0;
157  if ((pdg == 2212) && (ke < kmaxProton)) weight = 0;
158  if ((pdg == 2112) && (ke < kmaxNeutron)) weight = 0;
159 #ifdef DebugLog
160  if (weight == 0)
161  LogDebug("EcalSim") << "Ignore Track " << theTrack->GetTrackID()
162  << " Type " << theTrack->GetDefinition()->GetParticleName()
163  << " Kinetic Energy " << ke << " MeV";
164 #endif
165  }
166  }
167  }
168  G4LogicalVolume* lv = aStep->GetPreStepPoint()->GetTouchable()->GetVolume(0)->GetLogicalVolume();
169  if (useWeight && !any(noWeight,lv)) {
170  weight *= curve_LY(aStep);
171  if (useBirk) {
172  if (useBirkL3) weight *= getBirkL3(aStep);
173  else weight *= getAttenuation(aStep, birk1, birk2, birk3);
174  }
175  }
176  double wt1 = getResponseWt(theTrack);
177  double edep = aStep->GetTotalEnergyDeposit()*weight*wt1;
178  /*
179  if(wt2 != 1.0) {
180  std::cout << "ECalSD:: " << nameVolume
181  <<" LightWeight= " <<weight << " wt1= " <<wt1
182  << " wt2= " << wt2 << " "
183  << " Weighted Energy Deposit " << edep/MeV << " MeV"
184  << std::endl;
185  std::cout << theTrack->GetDefinition()->GetParticleName()
186  << " " << theTrack->GetKineticEnergy()
187  << " Id=" << theTrack->GetTrackID()
188  << " IdP=" << theTrack->GetParentID();
189  const G4VProcess* pr = theTrack->GetCreatorProcess();
190  if(pr) std::cout << " from " << pr->GetProcessName();
191  std::cout << std::endl;
192  }
193  */
194  if(wt2 > 0.0) { edep *= wt2; }
195 #ifdef DebugLog
196  LogDebug("EcalSim") << "ECalSD:: " << nameVolume
197  <<" Light Collection Efficiency " <<weight << ":" <<wt1
198  << " Weighted Energy Deposit " << edep/MeV << " MeV";
199 #endif
200  return edep;
201  }
202 }
203 
204 int ECalSD::getTrackID(G4Track* aTrack) {
205 
206  int primaryID(0);
207  bool flag(false);
208  if (storeTrack) {
209  G4LogicalVolume* lv = preStepPoint->GetTouchable()->GetVolume(0)->GetLogicalVolume();
210  if (any(useDepth1,lv)) {
211  flag = true;
212  } else if (any(useDepth2,lv)) {
213  flag = true;
214  }
215  }
216  if (flag) {
217  forceSave = true;
218  primaryID = aTrack->GetTrackID();
219  } else {
220  primaryID = CaloSD::getTrackID(aTrack);
221  }
222  return primaryID;
223 }
224 
225 uint16_t ECalSD::getDepth(G4Step * aStep) {
226  G4LogicalVolume* lv = aStep->GetPreStepPoint()->GetTouchable()->GetVolume(0)->GetLogicalVolume();
227  uint16_t ret = 0;
228  if (any(useDepth1,lv)) ret = 1;
229  else if (any(useDepth2,lv)) ret = 2;
230  else if (storeRL) ret = getRadiationLength(aStep);
231 #ifdef DebugLog
232  LogDebug("EcalSim") << "Volume " << lv->GetName() << " Depth " << ret;
233 #endif
234  return ret;
235 }
236 
237 uint16_t ECalSD::getRadiationLength(G4Step * aStep) {
238 
239  uint16_t thisX0 = 0;
240  if (aStep != NULL) {
241  G4StepPoint* hitPoint = aStep->GetPreStepPoint();
242  G4LogicalVolume* lv = hitPoint->GetTouchable()->GetVolume(0)->GetLogicalVolume();
243 
244  if (useWeight) {
245  G4ThreeVector localPoint = setToLocal(hitPoint->GetPosition(),
246  hitPoint->GetTouchable());
247  double crlength = crystalLength(lv);
248  double radl = hitPoint->GetMaterial()->GetRadlen();
249  double detz = (float)(0.5*crlength + localPoint.z());
250  thisX0 = (uint16_t)floor(detz/radl);
251  }
252  }
253  return thisX0;
254 }
255 
256 uint32_t ECalSD::setDetUnitId(G4Step * aStep) {
257  if (numberingScheme == 0) {
258  return EBDetId(1,1)();
259  } else {
260  getBaseNumber(aStep);
262  }
263 }
264 
266  if (scheme != 0) {
267  edm::LogInfo("EcalSim") << "EcalSD: updates numbering scheme for "
268  << GetName() << "\n";
269  if (numberingScheme) delete numberingScheme;
270  numberingScheme = scheme;
271  }
272 }
273 
274 
275 void ECalSD::initMap(G4String sd, const DDCompactView & cpv) {
276 
277  G4String attribute = "ReadOutName";
279  DDValue ddv(attribute,sd,0);
280  filter.setCriteria(ddv,DDCompOp::equals);
281  DDFilteredView fv(cpv);
282  fv.addFilter(filter);
283  fv.firstChild();
284 
285  std::vector<G4LogicalVolume*> lvused;
286  const G4LogicalVolumeStore * lvs = G4LogicalVolumeStore::GetInstance();
287  std::vector<G4LogicalVolume *>::const_iterator lvcite;
288  std::map<std::string, G4LogicalVolume *> nameMap;
289  for (auto lvi = lvs->begin(), lve = lvs->end(); lvi != lve; ++lvi)
290  nameMap.insert(std::make_pair((*lvi)->GetName(), *lvi));
291 
292  bool dodet=true;
293  while (dodet) {
294  const std::string &matname = fv.logicalPart().material().name().name();
295  const std::string &lvname = fv.logicalPart().name().name();
296  G4LogicalVolume* lv = nameMap[lvname];
297  if (depth1Name != " ") {
298  if (strncmp(lvname.c_str(), depth1Name.c_str(), 4) == 0) {
299  if (!any(useDepth1, lv)) {
300  useDepth1.push_back(lv);
301 #ifdef DebugLog
302  LogDebug("EcalSim") << "ECalSD::initMap Logical Volume " << lvname
303  <<" in Depth 1 volume list";
304 #endif
305  }
306  G4LogicalVolume* lvr = nameMap[lvname + "_refl"];
307  if (lvr != 0 && !any(useDepth1, lvr)) {
308  useDepth1.push_back(lvr);
309 #ifdef DebugLog
310  LogDebug("EcalSim") << "ECalSD::initMap Logical Volume " << lvname << "_refl"
311  <<" in Depth 1 volume list";
312 #endif
313  }
314  }
315  }
316  if (depth2Name != " ") {
317  if (strncmp(lvname.c_str(), depth2Name.c_str(), 4) == 0) {
318  if (!any(useDepth2, lv)) {
319  useDepth2.push_back(lv);
320 #ifdef DebugLog
321  LogDebug("EcalSim") << "ECalSD::initMap Logical Volume " << lvname
322  <<" in Depth 2 volume list";
323 #endif
324  }
325  G4LogicalVolume* lvr = nameMap[lvname + "_refl"];
326  if (lvr != 0 && !any(useDepth2,lvr)) {
327  useDepth2.push_back(lvr);
328 #ifdef DebugLog
329  LogDebug("EcalSim") << "ECalSD::initMap Logical Volume " << lvname << "_refl"
330  <<" in Depth 2 volume list";
331 #endif
332  }
333  }
334  }
335  if (lv != 0) {
336  if (crystalMat.size() == matname.size() && !strcmp(crystalMat.c_str(), matname.c_str())) {
337  if (!any(lvused,lv)) {
338  lvused.push_back(lv);
339  const DDSolid & sol = fv.logicalPart().solid();
340  const std::vector<double> & paras = sol.parameters();
341 #ifdef DebugLog
342  LogDebug("EcalSim") << "ECalSD::initMap (for " << sd << "): Solid "
343  << lvname << " Shape " << sol.shape()
344  << " Parameter 0 = "<< paras[0]
345  << " Logical Volume " << lv;
346 #endif
347  if (sol.shape() == ddtrap) {
348  double dz = 2*paras[0];
349  xtalLMap.insert(std::pair<G4LogicalVolume*,double>(lv,dz));
350  lv = nameMap[lvname + "_refl"];
351  if (lv != 0)
352  xtalLMap.insert(std::pair<G4LogicalVolume*,double>(lv,dz));
353  }
354  }
355  } else {
356  if (!any(noWeight,lv)) {
357  noWeight.push_back(lv);
358 #ifdef DebugLog
359  LogDebug("EcalSim") << "ECalSD::initMap Logical Volume " << lvname
360  << " Material " << matname <<" in noWeight list";
361 #endif
362  }
363  lv = nameMap[lvname];
364  if (lv != 0 && !any(noWeight,lv)) {
365  noWeight.push_back(lv);
366 #ifdef DebugLog
367  LogDebug("EcalSim") << "ECalSD::initMap Logical Volume " << lvname
368  << " Material " << matname <<" in noWeight list";
369 #endif
370  }
371  }
372  }
373  dodet = fv.next();
374  }
375 #ifdef DebugLog
376  LogDebug("EcalSim") << "ECalSD: Length Table for " << attribute << " = "
377  << sd << ":";
378  std::map<G4LogicalVolume*,double>::const_iterator ite = xtalLMap.begin();
379  int i=0;
380  for (; ite != xtalLMap.end(); ite++, i++) {
381  G4String name = "Unknown";
382  if (ite->first != 0) name = (ite->first)->GetName();
383  LogDebug("EcalSim") << " " << i << " " << ite->first << " " << name
384  << " L = " << ite->second;
385  }
386 #endif
387 }
388 
389 double ECalSD::curve_LY(G4Step* aStep) {
390 
391  G4StepPoint* stepPoint = aStep->GetPreStepPoint();
392  G4LogicalVolume* lv = stepPoint->GetTouchable()->GetVolume(0)->GetLogicalVolume();
393 
394  double weight = 1.;
395  G4ThreeVector localPoint = setToLocal(stepPoint->GetPosition(),
396  stepPoint->GetTouchable());
397 
398  double crlength = crystalLength(lv);
399 
400  if(ageingWithSlopeLY){
401  //position along the crystal in mm from 0 to 230 (in EB)
402  double depth = 0.5 * crlength + localPoint.z();
403 
404  if (depth >= -0.1 || depth <= crlength+0.1)
405  weight = ageing.calcLightCollectionEfficiencyWeighted(currentID.unitID(), depth/crlength);
406  }
407  else{
408  double dapd = 0.5 * crlength - localPoint.z();
409  if (dapd >= -0.1 || dapd <= crlength+0.1) {
410  if (dapd <= 100.)
411  weight = 1.0 + slopeLY - dapd * 0.01 * slopeLY;
412  } else {
413  edm::LogWarning("EcalSim") << "ECalSD: light coll curve : wrong distance "
414  << "to APD " << dapd << " crlength = "
415  << crlength <<" crystal name = " <<lv->GetName()
416  << " z of localPoint = " << localPoint.z()
417  << " take weight = " << weight;
418  }
419  }
420 #ifdef DebugLog
421  LogDebug("EcalSim") << "ECalSD, light coll curve : " << dapd
422  << " crlength = " << crlength
423  << " crystal name = " << lv->GetName()
424  << " z of localPoint = " << localPoint.z()
425  << " take weight = " << weight;
426 #endif
427  return weight;
428 }
429 
430 double ECalSD::crystalLength(G4LogicalVolume* lv) {
431 
432  double length= 230.;
433  std::map<G4LogicalVolume*,double>::const_iterator ite = xtalLMap.find(lv);
434  if (ite != xtalLMap.end()) length = ite->second;
435  return length;
436 }
437 
438 void ECalSD::getBaseNumber(const G4Step* aStep) {
439 
441  const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
442  int theSize = touch->GetHistoryDepth()+1;
443  if ( theBaseNumber.getCapacity() < theSize ) theBaseNumber.setSize(theSize);
444  //Get name and copy numbers
445  if ( theSize > 1 ) {
446  for (int ii = 0; ii < theSize ; ii++) {
447  theBaseNumber.addLevel(touch->GetVolume(ii)->GetName(),touch->GetReplicaNumber(ii));
448 #ifdef DebugLog
449  LogDebug("EcalSim") << "ECalSD::getBaseNumber(): Adding level " << ii
450  << ": " << touch->GetVolume(ii)->GetName() << "["
451  << touch->GetReplicaNumber(ii) << "]";
452 #endif
453  }
454  }
455 
456 }
457 
458 double ECalSD::getBirkL3(G4Step* aStep) {
459 
460  double weight = 1.;
461  double charge = aStep->GetPreStepPoint()->GetCharge();
462 
463  if (charge != 0. && aStep->GetStepLength() > 0) {
464  G4Material* mat = aStep->GetPreStepPoint()->GetMaterial();
465  double density = mat->GetDensity();
466  double dedx = aStep->GetTotalEnergyDeposit()/aStep->GetStepLength();
467  double rkb = birk1/density;
468  if (dedx > 0) {
469  weight = 1. - birkSlope*log(rkb*dedx);
470  if (weight < birkCut) weight = birkCut;
471  else if (weight > 1.) weight = 1.;
472  }
473 #ifdef DebugLog
474  LogDebug("EcalSim") << "ECalSD::getBirkL3 in " << mat->GetName()
475  << " Charge " << charge << " dE/dx " << dedx
476  << " Birk Const " << rkb << " Weight = " << weight
477  << " dE " << aStep->GetTotalEnergyDeposit();
478 #endif
479  }
480  return weight;
481 
482 }
483 
484 std::vector<double> ECalSD::getDDDArray(const std::string & str,
485  const DDsvalues_type & sv) {
486 
487 #ifdef DebugLog
488  LogDebug("EcalSim") << "ECalSD:getDDDArray called for " << str;
489 #endif
490  DDValue value(str);
491  if (DDfetch(&sv,value)) {
492 #ifdef DebugLog
493  LogDebug("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,
504  const DDsvalues_type & sv) {
505 
506 #ifdef DebugLog
507  LogDebug("EcalSim") << "ECalSD:getStringArray called for " << str;
508 #endif
509  DDValue value(str);
510  if (DDfetch(&sv,value)) {
511 #ifdef DebugLog
512  LogDebug("EcalSim") << value;
513 #endif
514  const std::vector<std::string> & fvec = value.strings();
515  return fvec;
516  } else {
517  std::vector<std::string> fvec;
518  return fvec;
519  }
520 }
#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
tuple ret
prodAgent to be discontinued
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:139
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:503
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:438
void addFilter(const DDFilter &, DDLogOp op=DDLogOp::AND)
virtual double getEnergyDeposit(G4Step *)
Definition: ECalSD.cc:137
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
const DDSolid & solid(void) const
Returns a reference object of the solid being the shape of this LogicalPart.
int ii
Definition: cuy.py:588
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:134
bool DDfetch(const DDsvalues_type *, DDValue &)
helper for retrieving DDValues from DDsvalues_type *.
Definition: DDsvalues.cc:80
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
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:35
static TrackerG4SimHitNumberingScheme & numberingScheme(const DDCompactView &cpv, const GeometricDet &det)
virtual uint16_t getDepth(G4Step *)
Definition: ECalSD.cc:225
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:133
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:19
double crystalLength(G4LogicalVolume *)
Definition: ECalSD.cc:430
DDSolidShape shape(void) const
The type of the solid.
Definition: DDSolid.cc:144
virtual int getTrackID(G4Track *)
Definition: ECalSD.cc:204
double getAttenuation(G4Step *aStep, double birk1, double birk2, double birk3)
Definition: CaloSD.cc:465
void initMap(G4String, const DDCompactView &)
Definition: ECalSD.cc:275
G4Track * theTrack
Definition: CaloSD.h:118
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:62
double curve_LY(G4Step *)
Definition: ECalSD.cc:389
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
std::vector< double > getDDDArray(const std::string &, const DDsvalues_type &)
Definition: ECalSD.cc:484
CaloHitID currentID
Definition: CaloSD.h:117
double sd
virtual uint16_t getRadiationLength(G4Step *)
Definition: ECalSD.cc:237
double getBirkL3(G4Step *)
Definition: ECalSD.cc:458
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
double getResponseWt(G4Track *)
Definition: CaloSD.cc:607
double birk2
Definition: ECalSD.h:55
bool firstChild()
set the current node to the first child ...
virtual uint32_t setDetUnitId(G4Step *)
Definition: ECalSD.cc:256
EcalNumberingScheme * numberingScheme
Definition: ECalSD.h:52
EnergyResolutionVsLumi ageing
Definition: ECalSD.h:61
std::map< G4LogicalVolume *, double > xtalLMap
Definition: ECalSD.h:58
int weight
Definition: histoStyle.py:50
uint32_t unitID() const
Definition: CaloHitID.h:22
void setNumberingScheme(EcalNumberingScheme *)
Definition: ECalSD.cc:265
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:87
const DDMaterial & material(void) const
Returns a reference object of the material this LogicalPart is made of.
void setCriteria(const DDValue &nameVal, DDCompOp, DDLogOp l=DDLogOp::AND, bool asString=true, bool merged=true)
Definition: DDFilter.cc:245
std::string crystalMat
Definition: ECalSD.h:57
virtual uint32_t getUnitID(const EcalBaseNumber &baseNumber) const =0
double birk3
Definition: ECalSD.h:55
The DDGenericFilter is a runtime-parametrized Filter looking on DDSpecifcs.
Definition: DDFilter.h:32
G4ThreeVector setToLocal(const G4ThreeVector &, const G4VTouchable *)
Definition: CaloSD.cc:296