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