CMS 3D CMS Logo

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