CMS 3D CMS Logo

DreamSD.cc
Go to the documentation of this file.
1 
8 
9 #include "G4LogicalVolumeStore.hh"
10 #include "G4LogicalVolume.hh"
11 #include "G4Step.hh"
12 #include "G4Track.hh"
13 #include "G4VProcess.hh"
14 #include "G4Poisson.hh"
15 
16 // Histogramming
19 #include <TTree.h>
20 
21 // Cherenkov
24 
25 #include "G4SystemOfUnits.hh"
26 #include "G4PhysicalConstants.hh"
27 
28 //________________________________________________________________________________________
30  const SensitiveDetectorCatalog & clg,
31  edm::ParameterSet const & p, const SimTrackManager* manager) :
32  CaloSD(name, cpv, clg, p, manager) {
33 
35  useBirk= m_EC.getParameter<bool>("UseBirkLaw");
36  doCherenkov_ = m_EC.getParameter<bool>("doCherenkov");
37  birk1 = m_EC.getParameter<double>("BirkC1")*(g/(MeV*cm2));
38  birk2 = m_EC.getParameter<double>("BirkC2");
39  birk3 = m_EC.getParameter<double>("BirkC3");
40  slopeLY= m_EC.getParameter<double>("SlopeLightYield");
41  readBothSide_ = m_EC.getUntrackedParameter<bool>("ReadBothSide", false);
42 
43  chAngleIntegrals_.reset(nullptr);
44 
45  edm::LogInfo("EcalSim") << "Constructing a DreamSD with name " << GetName() << "\n"
46  << "DreamSD:: Use of Birks law is set to "
47  << useBirk << " with three constants kB = "
48  << birk1 << ", C1 = " << birk2 << ", C2 = "
49  << birk3 << "\n"
50  << " Slope for Light yield is set to "
51  << slopeLY << "\n"
52  << " Parameterization of Cherenkov is set to "
53  << doCherenkov_ << " and readout both sides is "
54  << readBothSide_;
55 
56  initMap(name,cpv);
57 
58  // Init histogramming
60 
61  if ( !tfile.isAvailable() )
62  throw cms::Exception("BadConfig") << "TFileService unavailable: "
63  << "please add it to config file";
64 
65  ntuple_ = tfile->make<TTree>("tree","Cherenkov photons");
66  if (doCherenkov_) {
67  ntuple_->Branch("nphotons",&nphotons_,"nphotons/I");
68  ntuple_->Branch("px",px_,"px[nphotons]/F");
69  ntuple_->Branch("py",py_,"py[nphotons]/F");
70  ntuple_->Branch("pz",pz_,"pz[nphotons]/F");
71  ntuple_->Branch("x",x_,"x[nphotons]/F");
72  ntuple_->Branch("y",y_,"y[nphotons]/F");
73  ntuple_->Branch("z",z_,"z[nphotons]/F");
74  }
75 
76 }
77 
78 //________________________________________________________________________________________
79 double DreamSD::getEnergyDeposit(const G4Step* aStep) {
80 
81  // take into account light collection curve for crystals
82  double weight = curve_LY(aStep, side);
83  if (useBirk) weight *= getAttenuation(aStep, birk1, birk2, birk3);
84  double edep = aStep->GetTotalEnergyDeposit() * weight;
85 
86  // Get Cerenkov contribution
87  if ( doCherenkov_ ) { edep += cherenkovDeposit_( aStep ); }
88  LogDebug("EcalSim") << "DreamSD:: " << aStep->GetPreStepPoint()->GetPhysicalVolume()->GetName()
89  << " Side " << side
90  <<" Light Collection Efficiency " << weight
91  << " Weighted Energy Deposit " << edep/MeV
92  << " MeV";
93 
94  return edep;
95 }
96 
97 //________________________________________________________________________________________
99 
100  // Get the material and set properties if needed
101  DimensionMap::const_iterator ite = xtalLMap.begin();
102  const G4LogicalVolume* lv = (ite->first);
103  G4Material* material = lv->GetMaterial();
104  edm::LogInfo("EcalSim") << "DreamSD::initRun: Initializes for material "
105  << material->GetName() << " in " << lv->GetName();
106  materialPropertiesTable = material->GetMaterialPropertiesTable();
107  if ( !materialPropertiesTable ) {
108  if ( !setPbWO2MaterialProperties_( material ) ) {
109  edm::LogWarning("EcalSim") << "Couldn't retrieve material properties table\n"
110  << " Material = " << material->GetName();
111  }
112  materialPropertiesTable = material->GetMaterialPropertiesTable();
113  }
114 }
115 
116 
117 //________________________________________________________________________________________
118 uint32_t DreamSD::setDetUnitId(const G4Step * aStep) {
119  const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
120  uint32_t id = (touch->GetReplicaNumber(1))*10 + (touch->GetReplicaNumber(0));
121  side = readBothSide_ ? -1 : 1;
122  if(side < 0) { ++id; }
123  LogDebug("EcalSim") << "DreamSD:: ID " << id;
124  return id;
125 }
126 
127 
128 //________________________________________________________________________________________
129 void DreamSD::initMap(const std::string& sd, const DDCompactView & cpv) {
130 
131  G4String attribute = "ReadOutName";
133  DDFilteredView fv(cpv,filter);
134  fv.firstChild();
135 
136  const G4LogicalVolumeStore * lvs = G4LogicalVolumeStore::GetInstance();
137  std::vector<G4LogicalVolume *>::const_iterator lvcite;
138  bool dodet=true;
139  while (dodet) {
140  const DDSolid & sol = fv.logicalPart().solid();
141  std::vector<double> paras(sol.parameters());
142  G4String name = sol.name().name();
143  G4LogicalVolume* lv=nullptr;
144  for (lvcite = lvs->begin(); lvcite != lvs->end(); lvcite++)
145  if ((*lvcite)->GetName() == name) {
146  lv = (*lvcite);
147  break;
148  }
149  LogDebug("EcalSim") << "DreamSD::initMap (for " << sd << "): Solid "
150  << name << " Shape " << sol.shape() <<" Parameter 0 = "
151  << paras[0] << " Logical Volume " << lv;
152  double length = 0, width = 0;
153  // Set length to be the largest size, width the smallest
154  std::sort( paras.begin(), paras.end() );
155  length = 2.0*paras.back();
156  width = 2.0*paras.front();
157  xtalLMap.insert( std::pair<G4LogicalVolume*,Doubles>(lv,Doubles(length,width)) );
158  dodet = fv.next();
159  }
160  LogDebug("EcalSim") << "DreamSD: Length Table for " << attribute << " = "
161  << sd << ":";
162  DimensionMap::const_iterator ite = xtalLMap.begin();
163  int i=0;
164  for (; ite != xtalLMap.end(); ite++, i++) {
165  G4String name = "Unknown";
166  if (ite->first != nullptr) name = (ite->first)->GetName();
167  LogDebug("EcalSim") << " " << i << " " << ite->first << " " << name
168  << " L = " << ite->second.first
169  << " W = " << ite->second.second;
170  }
171 }
172 
173 //________________________________________________________________________________________
174 double DreamSD::curve_LY(const G4Step* aStep, int flag) {
175 
176  auto const stepPoint = aStep->GetPreStepPoint();
177  auto const lv = stepPoint->GetTouchable()->GetVolume(0)->GetLogicalVolume();
178  G4String nameVolume= lv->GetName();
179 
180  double weight = 1.;
181  G4ThreeVector localPoint = setToLocal(stepPoint->GetPosition(),
182  stepPoint->GetTouchable());
183  double crlength = crystalLength(lv);
184  double localz = localPoint.x();
185  double dapd = 0.5 * crlength - flag*localz; // Distance from closest APD
186  if (dapd >= -0.1 || dapd <= crlength+0.1) {
187  if (dapd <= 100.)
188  weight = 1.0 + slopeLY - dapd * 0.01 * slopeLY;
189  } else {
190  edm::LogWarning("EcalSim") << "DreamSD: light coll curve : wrong distance "
191  << "to APD " << dapd << " crlength = "
192  << crlength << " crystal name = " << nameVolume
193  << " z of localPoint = " << localz
194  << " take weight = " << weight;
195  }
196  LogDebug("EcalSim") << "DreamSD, light coll curve : " << dapd
197  << " crlength = " << crlength
198  << " crystal name = " << nameVolume
199  << " z of localPoint = " << localz
200  << " take weight = " << weight;
201  return weight;
202 }
203 
204 //________________________________________________________________________________________
205 const double DreamSD::crystalLength(G4LogicalVolume* lv) const {
206 
207  double length= -1.;
208  DimensionMap::const_iterator ite = xtalLMap.find(lv);
209  if (ite != xtalLMap.end()) length = ite->second.first;
210  return length;
211 
212 }
213 
214 //________________________________________________________________________________________
215 const double DreamSD::crystalWidth(G4LogicalVolume* lv) const {
216 
217  double width= -1.;
218  DimensionMap::const_iterator ite = xtalLMap.find(lv);
219  if (ite != xtalLMap.end()) width = ite->second.second;
220  return width;
221 
222 }
223 
224 
225 //________________________________________________________________________________________
226 // Calculate total cherenkov deposit
227 // Inspired by Geant4's Cherenkov implementation
228 double DreamSD::cherenkovDeposit_(const G4Step* aStep ) {
229 
230  double cherenkovEnergy = 0;
231  if (!materialPropertiesTable) return cherenkovEnergy;
232  G4Material* material = aStep->GetTrack()->GetMaterial();
233 
234  // Retrieve refractive index
235  G4MaterialPropertyVector* Rindex = materialPropertiesTable->GetProperty("RINDEX");
236  if ( Rindex == nullptr ) {
237  edm::LogWarning("EcalSim") << "Couldn't retrieve refractive index";
238  return cherenkovEnergy;
239  }
240 
241  // V.Ivanchenko - temporary close log output for 9.5
242  // Material refraction properties
243  int Rlength = Rindex->GetVectorLength() - 1;
244  double Pmin = Rindex->Energy(0);
245  double Pmax = Rindex->Energy(Rlength);
246  LogDebug("EcalSim") << "Material properties: " << "\n"
247  << " Pmin = " << Pmin
248  << " Pmax = " << Pmax;
249 
250  // Get particle properties
251  const G4StepPoint* pPreStepPoint = aStep->GetPreStepPoint();
252  const G4StepPoint* pPostStepPoint = aStep->GetPostStepPoint();
253  const G4ThreeVector& x0 = pPreStepPoint->GetPosition();
254  G4ThreeVector p0 = aStep->GetDeltaPosition().unit();
255  const G4DynamicParticle* aParticle = aStep->GetTrack()->GetDynamicParticle();
256  const double charge = aParticle->GetDefinition()->GetPDGCharge();
257  // beta is averaged over step
258  double beta = 0.5*( pPreStepPoint->GetBeta() + pPostStepPoint->GetBeta() );
259  double BetaInverse = 1.0/beta;
260 
261  LogDebug("EcalSim") << "Particle properties: " << "\n"
262  << " charge = " << charge
263  << " beta = " << beta;
264 
265  // Now get number of photons generated in this step
266  double meanNumberOfPhotons = getAverageNumberOfPhotons_( charge, beta, material, Rindex );
267  if ( meanNumberOfPhotons <= 0.0 ) { // Don't do anything
268  LogDebug("EcalSim") << "Mean number of photons is zero: " << meanNumberOfPhotons
269  << ", stopping here";
270  return cherenkovEnergy;
271  }
272 
273  // number of photons is in unit of Geant4...
274  meanNumberOfPhotons *= aStep->GetStepLength();
275 
276  // Now get a poisson distribution
277  int numPhotons = static_cast<int>( G4Poisson(meanNumberOfPhotons) );
278  //edm::LogVerbatim("EcalSim") << "Number of photons = " << numPhotons;
279  if ( numPhotons <= 0 ) {
280  LogDebug("EcalSim") << "Poission number of photons is zero: " << numPhotons
281  << ", stopping here";
282  return cherenkovEnergy;
283  }
284 
285  // Material refraction properties
286  double dp = Pmax - Pmin;
287  double maxCos = BetaInverse / (*Rindex)[Rlength];
288  double maxSin2 = (1.0 - maxCos) * (1.0 + maxCos);
289 
290  // Finally: get contribution of each photon
291  for ( int iPhoton = 0; iPhoton<numPhotons; ++iPhoton ) {
292 
293  // Determine photon momentum
294  double randomNumber;
295  double sampledMomentum, sampledRI;
296  double cosTheta, sin2Theta;
297 
298  // sample a momentum (not sure why this is needed!)
299  do {
300  randomNumber = G4UniformRand();
301  sampledMomentum = Pmin + randomNumber * dp;
302  sampledRI = Rindex->Value(sampledMomentum);
303  cosTheta = BetaInverse / sampledRI;
304 
305  sin2Theta = (1.0 - cosTheta)*(1.0 + cosTheta);
306  randomNumber = G4UniformRand();
307 
308  } while (randomNumber*maxSin2 > sin2Theta);
309 
310  // Generate random position of photon on cone surface
311  // defined by Theta
312  randomNumber = G4UniformRand();
313 
314  double phi = twopi*randomNumber;
315  double sinPhi = sin(phi);
316  double cosPhi = cos(phi);
317 
318  // Create photon momentum direction vector
319  // The momentum direction is still w.r.t. the coordinate system where the primary
320  // particle direction is aligned with the z axis
321  double sinTheta = sqrt(sin2Theta);
322  double px = sinTheta*cosPhi;
323  double py = sinTheta*sinPhi;
324  double pz = cosTheta;
325  G4ThreeVector photonDirection(px, py, pz);
326 
327  // Rotate momentum direction back to global (crystal) reference system
328  photonDirection.rotateUz(p0);
329 
330  // Create photon position and momentum
331  randomNumber = G4UniformRand();
332  G4ThreeVector photonPosition = x0 + randomNumber * aStep->GetDeltaPosition();
333  G4ThreeVector photonMomentum = sampledMomentum*photonDirection;
334 
335  // Collect energy on APD
336  cherenkovEnergy += getPhotonEnergyDeposit_( photonMomentum, photonPosition, aStep );
337 
338  // Ntuple variables
339  nphotons_ = numPhotons;
340  px_[iPhoton] = photonMomentum.x();
341  py_[iPhoton] = photonMomentum.y();
342  pz_[iPhoton] = photonMomentum.z();
343  x_[iPhoton] = photonPosition.x();
344  y_[iPhoton] = photonPosition.y();
345  z_[iPhoton] = photonPosition.z();
346  }
347 
348  // Fill ntuple
349  ntuple_->Fill();
350 
351 
352  return cherenkovEnergy;
353 
354 }
355 
356 
357 //________________________________________________________________________________________
358 // Returns number of photons produced per GEANT-unit (millimeter) in the current medium.
359 // From G4Cerenkov.cc
361  const double beta,
362  const G4Material* aMaterial,
363  const G4MaterialPropertyVector* Rindex )
364 {
365  const G4double rFact = 369.81/(eV * cm);
366 
367  if( beta <= 0.0 ) return 0.0;
368 
369  double BetaInverse = 1./beta;
370 
371  // Vectors used in computation of Cerenkov Angle Integral:
372  // - Refraction Indices for the current material
373  // - new G4PhysicsOrderedFreeVector allocated to hold CAI's
374 
375  // Min and Max photon momenta
376  int Rlength = Rindex->GetVectorLength() - 1;
377  double Pmin = Rindex->Energy(0);
378  double Pmax = Rindex->Energy(Rlength);
379 
380  // Min and Max Refraction Indices
381  double nMin = (*Rindex)[0];
382  double nMax = (*Rindex)[Rlength];
383 
384  // Max Cerenkov Angle Integral
385  double CAImax = chAngleIntegrals_.get()->GetMaxValue();
386 
387  double dp = 0., ge = 0., CAImin = 0.;
388 
389  // If n(Pmax) < 1/Beta -- no photons generated
390  if ( nMax < BetaInverse) { }
391 
392  // otherwise if n(Pmin) >= 1/Beta -- photons generated
393  else if (nMin > BetaInverse) {
394  dp = Pmax - Pmin;
395  ge = CAImax;
396  }
397  // If n(Pmin) < 1/Beta, and n(Pmax) >= 1/Beta, then
398  // we need to find a P such that the value of n(P) == 1/Beta.
399  // Interpolation is performed by the GetPhotonEnergy() and
400  // GetProperty() methods of the G4MaterialPropertiesTable and
401  // the GetValue() method of G4PhysicsVector.
402  else {
403  Pmin = Rindex->Value(BetaInverse);
404  dp = Pmax - Pmin;
405  // need boolean for current implementation of G4PhysicsVector
406  // ==> being phased out
407  double CAImin = chAngleIntegrals_->Value(Pmin);
408  ge = CAImax - CAImin;
409 
410  }
411 
412  // Calculate number of photons
413  double numPhotons = rFact * charge/eplus * charge/eplus *
414  (dp - ge * BetaInverse*BetaInverse);
415 
416  LogDebug("EcalSim") << "@SUB=getAverageNumberOfPhotons"
417  << "CAImin = " << CAImin << "\n"
418  << "CAImax = " << CAImax << "\n"
419  << "dp = " << dp << ", ge = " << ge << "\n"
420  << "numPhotons = " << numPhotons;
421 
422 
423 
424  return numPhotons;
425 
426 }
427 
428 
429 //________________________________________________________________________________________
430 // Set lead tungstate material properties on the fly.
431 // Values from Ts42 detector construction
432 bool DreamSD::setPbWO2MaterialProperties_( G4Material* aMaterial ) {
433 
434  std::string pbWO2Name("E_PbWO4");
435  if ( pbWO2Name != aMaterial->GetName() ) { // Wrong material!
436  edm::LogWarning("EcalSim") << "This is not the right material: "
437  << "expecting " << pbWO2Name
438  << ", got " << aMaterial->GetName();
439  return false;
440  }
441 
442  G4MaterialPropertiesTable* table = new G4MaterialPropertiesTable();
443 
444  // Refractive index as a function of photon momentum
445  // FIXME: Should somehow put that in the configuration
446  const int nEntries = 14;
447  double PhotonEnergy[nEntries] = { 1.7712*eV, 1.8368*eV, 1.90745*eV, 1.98375*eV, 2.0664*eV,
448  2.15625*eV, 2.25426*eV, 2.3616*eV, 2.47968*eV, 2.61019*eV,
449  2.75521*eV, 2.91728*eV, 3.09961*eV, 3.30625*eV };
450  double RefractiveIndex[nEntries] = { 2.17728, 2.18025, 2.18357, 2.18753, 2.19285,
451  2.19813, 2.20441, 2.21337, 2.22328, 2.23619,
452  2.25203, 2.27381, 2.30282, 2.34666 };
453 
454  table->AddProperty( "RINDEX", PhotonEnergy, RefractiveIndex, nEntries );
455  aMaterial->SetMaterialPropertiesTable(table); // FIXME: could this leak? What does G4 do?
456 
457  // Calculate Cherenkov angle integrals:
458  // This is an ad-hoc solution (we hold it in the class, not in the material)
459  chAngleIntegrals_.reset(new G4PhysicsOrderedFreeVector());
460 
461  int index = 0;
462  double currentRI = RefractiveIndex[index];
463  double currentPM = PhotonEnergy[index];
464  double currentCAI = 0.0;
465  chAngleIntegrals_.get()->InsertValues(currentPM, currentCAI);
466  double prevPM = currentPM;
467  double prevCAI = currentCAI;
468  double prevRI = currentRI;
469  while ( ++index < nEntries ) {
470  currentRI = RefractiveIndex[index];
471  currentPM = PhotonEnergy[index];
472  currentCAI = 0.5*(1.0/(prevRI*prevRI) + 1.0/(currentRI*currentRI));
473  currentCAI = prevCAI + (currentPM - prevPM) * currentCAI;
474 
475  chAngleIntegrals_.get()->InsertValues(currentPM, currentCAI);
476 
477  prevPM = currentPM;
478  prevCAI = currentCAI;
479  prevRI = currentRI;
480  }
481 
482  LogDebug("EcalSim") << "Material properties set for " << aMaterial->GetName();
483 
484  return true;
485 
486 }
487 
488 
489 //________________________________________________________________________________________
490 // Calculate energy deposit of a photon on APD
491 // - simple tracing to APD position (straight line);
492 // - configurable reflection probability if not straight to APD;
493 // - APD response function
494 double DreamSD::getPhotonEnergyDeposit_( const G4ThreeVector& p,
495  const G4ThreeVector& x,
496  const G4Step* aStep )
497 {
498 
499  double energy = 0;
500 
501  // Crystal dimensions
502 
503  //edm::LogVerbatim("EcalSim") << p << x;
504 
505  // 1. Check if this photon goes straight to the APD:
506  // - assume that APD is at x=xtalLength/2.0
507  // - extrapolate from x=x0 to x=xtalLength/2.0 using momentum in x-y
508 
509  G4StepPoint* stepPoint = aStep->GetPreStepPoint();
510  G4LogicalVolume* lv = stepPoint->GetTouchable()->GetVolume(0)->GetLogicalVolume();
511  G4String nameVolume= lv->GetName();
512 
513  double crlength = crystalLength(lv);
514  double crwidth = crystalWidth(lv);
515  double dapd = 0.5 * crlength - x.x(); // Distance from closest APD
516  double y = p.y()/p.x()*dapd;
517 
518  LogDebug("EcalSim") << "Distance to APD: " << dapd
519  << " - y at APD: " << y;
520 
521  // Not straight: compute probability
522  if ( fabs(y)>crwidth*0.5 ) {
523 
524  }
525 
526  // 2. Retrieve efficiency for this wavelength (in nm, from MeV)
527  double waveLength = p.mag()*1.239e8;
528 
529 
530  energy = p.mag()*PMTResponse::getEfficiency(waveLength);
531 
532  LogDebug("EcalSim") << "Wavelength: " << waveLength << " - Energy: " << energy;
533 
534  return energy;
535 
536 }
#define LogDebug(id)
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
const std::vector< double > & parameters(void) const
Give the parameters of the solid.
Definition: DDSolid.cc:144
const DDLogicalPart & logicalPart() const
The logical-part of the current node in the filtered-view.
const N & name() const
Definition: DDBase.h:74
bool readBothSide_
Definition: DreamSD.h:59
double getPhotonEnergyDeposit_(const G4ParticleMomentum &p, const G4ThreeVector &x, const G4Step *aStep)
Returns energy deposit for a given photon.
Definition: DreamSD.cc:494
std::pair< double, double > Doubles
Definition: DreamSD.h:37
Definition: CaloSD.h:37
float y_[MAXPHOTONS]
Definition: DreamSD.h:73
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
uint32_t setDetUnitId(const G4Step *) override
Definition: DreamSD.cc:118
Definition: weight.py:1
double birk1
Definition: DreamSD.h:60
bool setPbWO2MaterialProperties_(G4Material *aMaterial)
Sets material properties at run-time...
Definition: DreamSD.cc:432
T * make(const Args &...args) const
make new ROOT object
Definition: TFileService.h:64
double birk3
Definition: DreamSD.h:60
const DDSolid & solid(void) const
Returns a reference object of the solid being the shape of this LogicalPart.
void initRun() override
Definition: DreamSD.cc:98
float px_[MAXPHOTONS]
Definition: DreamSD.h:72
Compact representation of the geometrical detector hierarchy.
Definition: DDCompactView.h:80
TTree * ntuple_
Definition: DreamSD.h:70
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
A DDSolid represents the shape of a part.
Definition: DDSolid.h:39
const double MeV
float z_[MAXPHOTONS]
Definition: DreamSD.h:73
bool next()
set current node to the next node in the filtered tree
double cherenkovDeposit_(const G4Step *aStep)
Returns the total energy due to Cherenkov radiation.
Definition: DreamSD.cc:228
T sqrt(T t)
Definition: SSEVec.h:18
int nphotons_
Definition: DreamSD.h:71
float pz_[MAXPHOTONS]
Definition: DreamSD.h:72
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
int side
Definition: DreamSD.h:64
bool isAvailable() const
Definition: Service.h:46
double getAverageNumberOfPhotons_(const double charge, const double beta, const G4Material *aMaterial, const G4MaterialPropertyVector *rIndex)
Returns average number of photons created by track.
Definition: DreamSD.cc:360
DDSolidShape shape(void) const
The type of the solid.
Definition: DDSolid.cc:138
double birk2
Definition: DreamSD.h:60
G4MaterialPropertiesTable * materialPropertiesTable
Definition: DreamSD.h:68
double curve_LY(const G4Step *, int)
Definition: DreamSD.cc:174
const double crystalWidth(G4LogicalVolume *) const
Definition: DreamSD.cc:215
G4ThreeVector setToLocal(const G4ThreeVector &, const G4VTouchable *) const
Definition: CaloSD.cc:307
DimensionMap xtalLMap
Definition: DreamSD.h:62
double sd
float x_[MAXPHOTONS]
Definition: DreamSD.h:73
std::unique_ptr< G4PhysicsOrderedFreeVector > chAngleIntegrals_
Table of Cherenkov angle integrals vs photon momentum.
Definition: DreamSD.h:67
double getEnergyDeposit(const G4Step *) override
Definition: DreamSD.cc:79
double getAttenuation(const G4Step *aStep, double birk1, double birk2, double birk3) const
Definition: CaloSD.cc:462
bool doCherenkov_
Definition: DreamSD.h:59
double slopeLY
Definition: DreamSD.h:61
DreamSD(const std::string &, const DDCompactView &, const SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *)
Definition: DreamSD.cc:29
bool firstChild()
set the current node to the first child ...
float py_[MAXPHOTONS]
Definition: DreamSD.h:72
const double crystalLength(G4LogicalVolume *) const
Definition: DreamSD.cc:205
const std::string & name() const
Returns the name.
Definition: DDName.cc:53
void initMap(const std::string &, const DDCompactView &)
Definition: DreamSD.cc:129
static const double getEfficiency(const double &waveLengthNm)
Return efficiency for given photon wavelength (in nm)
Definition: PMTResponse.cc:6
bool useBirk
Definition: DreamSD.h:59