CMS 3D CMS Logo

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