ECCE @ EIC Software
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4Scintillation.cc
Go to the documentation of this file. Or view the newest version in sPHENIX GitHub for file G4Scintillation.cc
1 // ********************************************************************
2 // * License and Disclaimer *
3 // * *
4 // * The Geant4 software is copyright of the Copyright Holders of *
5 // * the Geant4 Collaboration. It is provided under the terms and *
6 // * conditions of the Geant4 Software License, included in the file *
7 // * LICENSE and available at http://cern.ch/geant4/license . These *
8 // * include a list of copyright holders. *
9 // * *
10 // * Neither the authors of this software system, nor their employing *
11 // * institutes,nor the agencies providing financial support for this *
12 // * work make any representation or warranty, express or implied, *
13 // * regarding this software system or assume any liability for its *
14 // * use. Please see the license in the file LICENSE and URL above *
15 // * for the full disclaimer and the limitation of liability. *
16 // * *
17 // * This code implementation is the result of the scientific and *
18 // * technical work of the GEANT4 collaboration. *
19 // * By using, copying, modifying or distributing the software (or *
20 // * any work based on the software) you agree to acknowledge its *
21 // * use in resulting scientific publications, and indicate your *
22 // * acceptance of all terms of the Geant4 Software license. *
23 // ********************************************************************
24 //
25 //
26 //
28 // Scintillation Light Class Implementation
30 //
31 // File: G4Scintillation.cc
32 // Description: RestDiscrete Process - Generation of Scintillation Photons
33 // Version: 1.0
34 // Created: 1998-11-07
35 // Author: Peter Gumplinger
36 // Updated: 2010-10-20 Allow the scintillation yield to be a function
37 // of energy deposited by particle type
38 // Thanks to Zach Hartwig (Department of Nuclear
39 // Science and Engineeering - MIT)
40 // 2010-09-22 by Peter Gumplinger
41 // > scintillation rise time included, thanks to
42 // > Martin Goettlich/DESY
43 // 2005-08-17 by Peter Gumplinger
44 // > change variable name MeanNumPhotons -> MeanNumberOfPhotons
45 // 2005-07-28 by Peter Gumplinger
46 // > add G4ProcessType to constructor
47 // 2004-08-05 by Peter Gumplinger
48 // > changed StronglyForced back to Forced in GetMeanLifeTime
49 // 2002-11-21 by Peter Gumplinger
50 // > change to use G4Poisson for small MeanNumberOfPhotons
51 // 2002-11-07 by Peter Gumplinger
52 // > now allow for fast and slow scintillation component
53 // 2002-11-05 by Peter Gumplinger
54 // > now use scintillation constants from G4Material
55 // 2002-05-09 by Peter Gumplinger
56 // > use only the PostStepPoint location for the origin of
57 // scintillation photons when energy is lost to the medium
58 // by a neutral particle
59 // 2000-09-18 by Peter Gumplinger
60 // > change: aSecondaryPosition=x0+rand*aStep.GetDeltaPosition();
61 // aSecondaryTrack->SetTouchable(0);
62 // 2001-09-17, migration of Materials to pure STL (mma)
63 // 2003-06-03, V.Ivanchenko fix compilation warnings
64 //
65 // mail: gum@triumf.ca
66 //
68 
69 #include "G4ios.hh"
70 #include "globals.hh"
71 #include "G4PhysicalConstants.hh"
72 #include "G4SystemOfUnits.hh"
73 #include "G4ParticleTypes.hh"
74 #include "G4EmProcessSubType.hh"
76 
77 #include "G4Scintillation.hh"
78 
80  G4ProcessType type)
81  : G4VRestDiscreteProcess(processName, type) ,
82  fTrackSecondariesFirst(false),
83  fFiniteRiseTime(false),
84  fYieldFactor(1.0),
85  fExcitationRatio(1.0),
86  fScintillationByParticleType(false),
87  fScintillationTrackInfo(false),
88  fStackingFlag(true),
89  fNumPhotons(0),
90  fEmSaturation(nullptr)
91 {
93 
94 #ifdef G4DEBUG_SCINTILLATION
95  ScintTrackEDep = 0.;
96  ScintTrackYield = 0.;
97 #endif
98 
99  fFastIntegralTable = nullptr;
100  fSlowIntegralTable = nullptr;
101 
102  if (verboseLevel>0) {
103  G4cout << GetProcessName() << " is created " << G4endl;
104  }
105 }
106 
108 {
109  if (fFastIntegralTable != nullptr) {
111  delete fFastIntegralTable;
112  }
113  if (fSlowIntegralTable != nullptr) {
115  delete fSlowIntegralTable;
116  }
117 }
118 
120 {
121  return (aParticleType.GetPDGCharge() == 0.0 ||
122  aParticleType.IsShortLived()) ? false : true;
123 }
124 
126 {
127  if (fFastIntegralTable != nullptr) {
129  delete fFastIntegralTable;
130  fFastIntegralTable = nullptr;
131  }
132  if (fSlowIntegralTable != nullptr) {
134  delete fSlowIntegralTable;
135  fSlowIntegralTable = nullptr;
136  }
138 }
139 
140 
141 // AtRestDoIt
142 // ----------
143 //
145 G4Scintillation::AtRestDoIt(const G4Track& aTrack, const G4Step& aStep)
146 
147 // This routine simply calls the equivalent PostStepDoIt since all the
148 // necessary information resides in aStep.GetTotalEnergyDeposit()
149 
150 {
151  return G4Scintillation::PostStepDoIt(aTrack, aStep);
152 }
153 
154 // PostStepDoIt
155 // -------------
156 //
158 G4Scintillation::PostStepDoIt(const G4Track& aTrack, const G4Step& aStep)
159 
160 // This routine is called for each tracking step of a charged particle
161 // in a scintillator. A Poisson/Gauss-distributed number of photons is
162 // generated according to the scintillation yield formula, distributed
163 // evenly along the track segment and uniformly into 4pi.
164 
165 {
166  aParticleChange.Initialize(aTrack);
167  fNumPhotons = 0;
168 
169  const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
170  const G4Material* aMaterial = aTrack.GetMaterial();
171 
172  G4StepPoint* pPreStepPoint = aStep.GetPreStepPoint();
173  G4StepPoint* pPostStepPoint = aStep.GetPostStepPoint();
174 
175  G4ThreeVector x0 = pPreStepPoint->GetPosition();
176  G4ThreeVector p0 = aStep.GetDeltaPosition().unit();
177  G4double t0 = pPreStepPoint->GetGlobalTime();
178 
179  G4double TotalEnergyDeposit = aStep.GetTotalEnergyDeposit();
180 
181  G4MaterialPropertiesTable* aMaterialPropertiesTable =
182  aMaterial->GetMaterialPropertiesTable();
183  if (!aMaterialPropertiesTable)
184  return G4VRestDiscreteProcess::PostStepDoIt(aTrack, aStep);
185 
186  G4MaterialPropertyVector* Fast_Intensity =
187  aMaterialPropertiesTable->GetProperty(kFASTCOMPONENT);
188  G4MaterialPropertyVector* Slow_Intensity =
189  aMaterialPropertiesTable->GetProperty(kSLOWCOMPONENT);
190 
191  if (!Fast_Intensity && !Slow_Intensity )
192  return G4VRestDiscreteProcess::PostStepDoIt(aTrack, aStep);
193 
194  G4int nscnt = 1;
195  if (Fast_Intensity && Slow_Intensity) nscnt = 2;
196 
197  G4double ScintillationYield = 0.;
198 
199  // Scintillation depends on particle type, energy deposited
201 
202  ScintillationYield =
204 
205  // The default linear scintillation process
206  } else {
207 
208  ScintillationYield = aMaterialPropertiesTable->
209  GetConstProperty(kSCINTILLATIONYIELD);
210 
211  // Units: [# scintillation photons / MeV]
212  ScintillationYield *= fYieldFactor;
213  }
214 
215  G4double ResolutionScale = aMaterialPropertiesTable->
216  GetConstProperty(kRESOLUTIONSCALE);
217 
218  // Birks law saturation:
219 
220  //G4double constBirks = 0.0;
221 
222  //constBirks = aMaterial->GetIonisation()->GetBirksConstant();
223 
224  G4double MeanNumberOfPhotons;
225 
226  // Birk's correction via fEmSaturation and specifying scintillation by
227  // by particle type are physically mutually exclusive
228 
230  MeanNumberOfPhotons = ScintillationYield;
231  else if (fEmSaturation)
232  MeanNumberOfPhotons = ScintillationYield*
234  else
235  MeanNumberOfPhotons = ScintillationYield*TotalEnergyDeposit;
236 
237  if (MeanNumberOfPhotons > 10.)
238  {
239  G4double sigma = ResolutionScale * std::sqrt(MeanNumberOfPhotons);
240  fNumPhotons=G4int(G4RandGauss::shoot(MeanNumberOfPhotons,sigma)+0.5);
241  }
242  else
243  {
244  fNumPhotons = G4int(G4Poisson(MeanNumberOfPhotons));
245  }
246 
247  if ( fNumPhotons <= 0 || !fStackingFlag )
248  {
249  // return unchanged particle and no secondaries
250 
252 
253  return G4VRestDiscreteProcess::PostStepDoIt(aTrack, aStep);
254  }
255 
257 
259 
261  if (aTrack.GetTrackStatus() == fAlive )
263  }
264 
266 
267  G4int materialIndex = aMaterial->GetIndex();
268 
269  // Retrieve the Scintillation Integral for this material
270  // new G4PhysicsOrderedFreeVector allocated to hold CII's
271 
272  G4int Num = fNumPhotons;
273 
274  for (G4int scnt = 1; scnt <= nscnt; scnt++) {
275 
276  G4double ScintillationTime = 0.*ns;
277  G4double ScintillationRiseTime = 0.*ns;
278  G4PhysicsOrderedFreeVector* ScintillationIntegral = nullptr;
279  G4ScintillationType ScintillationType = Slow;
280 
281  if (scnt == 1) {
282  if (nscnt == 1) {
283  if (Fast_Intensity) {
284  ScintillationTime = aMaterialPropertiesTable->
285  GetConstProperty(kFASTTIMECONSTANT);
286  if (fFiniteRiseTime) {
287  ScintillationRiseTime = aMaterialPropertiesTable->
288  GetConstProperty(kFASTSCINTILLATIONRISETIME);
289  }
290  ScintillationType = Fast;
291  ScintillationIntegral =
293  ((*fFastIntegralTable)(materialIndex));
294  }
295  if (Slow_Intensity) {
296  ScintillationTime = aMaterialPropertiesTable->
297  GetConstProperty(kSLOWTIMECONSTANT);
298  if (fFiniteRiseTime) {
299  ScintillationRiseTime = aMaterialPropertiesTable->
300  GetConstProperty(kSLOWSCINTILLATIONRISETIME);
301  }
302  ScintillationType = Slow;
303  ScintillationIntegral =
305  ((*fSlowIntegralTable)(materialIndex));
306  }
307  }
308  else {
309  G4double yieldRatio = aMaterialPropertiesTable->
310  GetConstProperty(kYIELDRATIO);
311  if ( fExcitationRatio == 1.0 || fExcitationRatio == 0.0) {
312  Num = G4int (std::min(yieldRatio,1.0) * fNumPhotons);
313  }
314  else {
316  }
317  ScintillationTime = aMaterialPropertiesTable->
318  GetConstProperty(kFASTTIMECONSTANT);
319  if (fFiniteRiseTime) {
320  ScintillationRiseTime = aMaterialPropertiesTable->
321  GetConstProperty(kFASTSCINTILLATIONRISETIME);
322  }
323  ScintillationType = Fast;
324  ScintillationIntegral =
326  ((*fFastIntegralTable)(materialIndex));
327  }
328  }
329  else {
330  Num = fNumPhotons - Num;
331  ScintillationTime = aMaterialPropertiesTable->
332  GetConstProperty(kSLOWTIMECONSTANT);
333  if (fFiniteRiseTime) {
334  ScintillationRiseTime = aMaterialPropertiesTable->
335  GetConstProperty(kSLOWSCINTILLATIONRISETIME);
336  }
337  ScintillationType = Slow;
338  ScintillationIntegral =
340  ((*fSlowIntegralTable)(materialIndex));
341  }
342 
343  if (!ScintillationIntegral) continue;
344 
345  // Max Scintillation Integral
346 
347  G4double CIImax = ScintillationIntegral->GetMaxValue();
348 
349  for (G4int i = 0; i < Num; i++) {
350 
351  // Determine photon energy
352 
353  G4double CIIvalue = G4UniformRand()*CIImax;
354  G4double sampledEnergy =
355  ScintillationIntegral->GetEnergy(CIIvalue);
356 
357  if (verboseLevel>1) {
358  G4cout << "sampledEnergy = " << sampledEnergy << G4endl;
359  G4cout << "CIIvalue = " << CIIvalue << G4endl;
360  }
361 
362  // Generate random photon direction
363 
364  G4double cost = 1. - 2.*G4UniformRand();
365  G4double sint = std::sqrt((1.-cost)*(1.+cost));
366 
368  G4double sinp = std::sin(phi);
369  G4double cosp = std::cos(phi);
370 
371  G4double px = sint*cosp;
372  G4double py = sint*sinp;
373  G4double pz = cost;
374 
375  // Create photon momentum direction vector
376 
377  G4ParticleMomentum photonMomentum(px, py, pz);
378 
379  // Determine polarization of new photon
380 
381  G4double sx = cost*cosp;
382  G4double sy = cost*sinp;
383  G4double sz = -sint;
384 
385  G4ThreeVector photonPolarization(sx, sy, sz);
386 
387  G4ThreeVector perp = photonMomentum.cross(photonPolarization);
388 
389  phi = twopi*G4UniformRand();
390  sinp = std::sin(phi);
391  cosp = std::cos(phi);
392 
393  photonPolarization = cosp * photonPolarization + sinp * perp;
394 
395  photonPolarization = photonPolarization.unit();
396 
397  // Generate a new photon:
398 
399  G4DynamicParticle* aScintillationPhoton =
401  photonMomentum);
402  aScintillationPhoton->SetPolarization
403  (photonPolarization.x(),
404  photonPolarization.y(),
405  photonPolarization.z());
406 
407  aScintillationPhoton->SetKineticEnergy(sampledEnergy);
408 
409  // Generate new G4Track object:
410 
411  G4double rand;
412 
413  if (aParticle->GetDefinition()->GetPDGCharge() != 0) {
414  rand = G4UniformRand();
415  } else {
416  rand = 1.0;
417  }
418 
419  G4double delta = rand * aStep.GetStepLength();
420  G4double deltaTime = delta / (pPreStepPoint->GetVelocity()+
421  rand*(pPostStepPoint->GetVelocity()-
422  pPreStepPoint->GetVelocity())/2.);
423 
424  // emission time distribution
425  if (ScintillationRiseTime==0.0) {
426  deltaTime = deltaTime -
427  ScintillationTime * std::log( G4UniformRand() );
428  } else {
429  deltaTime = deltaTime +
430  sample_time(ScintillationRiseTime, ScintillationTime);
431  }
432 
433  G4double aSecondaryTime = t0 + deltaTime;
434 
435  G4ThreeVector aSecondaryPosition =
436  x0 + rand * aStep.GetDeltaPosition();
437 
438  G4Track* aSecondaryTrack = new G4Track(aScintillationPhoton,
439  aSecondaryTime,
440  aSecondaryPosition);
441 
442  aSecondaryTrack->SetTouchableHandle(
444  // aSecondaryTrack->SetTouchableHandle((G4VTouchable*)0);
445 
446  aSecondaryTrack->SetParentID(aTrack.GetTrackID());
447 
448  if (fScintillationTrackInfo) aSecondaryTrack->
449  SetUserInformation(new G4ScintillationTrackInformation(ScintillationType));
450 
451  aParticleChange.AddSecondary(aSecondaryTrack);
452 
453  }
454  }
455 
456  if (verboseLevel>0) {
457  G4cout << "\n Exiting from G4Scintillation::DoIt -- NumberOfSecondaries = "
459  }
460 
461  return G4VRestDiscreteProcess::PostStepDoIt(aTrack, aStep);
462 }
463 
464 // BuildThePhysicsTable for the scintillation process
465 // --------------------------------------------------
466 //
467 
469 {
470  if (fFastIntegralTable && fSlowIntegralTable) return;
471 
472  const G4MaterialTable* theMaterialTable =
474  G4int numOfMaterials = G4Material::GetNumberOfMaterials();
475 
476  // create new physics table
477 
479  new G4PhysicsTable(numOfMaterials);
481  new G4PhysicsTable(numOfMaterials);
482 
483  // loop for materials
484 
485  for (G4int i=0 ; i < numOfMaterials; i++)
486  {
487  G4PhysicsOrderedFreeVector* aPhysicsOrderedFreeVector =
489  G4PhysicsOrderedFreeVector* bPhysicsOrderedFreeVector =
491 
492  // Retrieve vector of scintillation wavelength intensity for
493  // the material from the material's optical properties table.
494 
495  G4Material* aMaterial = (*theMaterialTable)[i];
496 
497  G4MaterialPropertiesTable* aMaterialPropertiesTable =
498  aMaterial->GetMaterialPropertiesTable();
499 
500  if (aMaterialPropertiesTable) {
501 
502  G4MaterialPropertyVector* theFastLightVector =
503  aMaterialPropertiesTable->GetProperty(kFASTCOMPONENT);
504 
505  if (theFastLightVector) {
506 
507  // Retrieve the first intensity point in vector
508  // of (photon energy, intensity) pairs
509 
510  G4double currentIN = (*theFastLightVector)[0];
511 
512  if (currentIN >= 0.0) {
513 
514  // Create first (photon energy, Scintillation
515  // Integral pair
516 
517  G4double currentPM = theFastLightVector->Energy(0);
518 
519  G4double currentCII = 0.0;
520 
521  aPhysicsOrderedFreeVector->
522  InsertValues(currentPM , currentCII);
523 
524  // Set previous values to current ones prior to loop
525 
526  G4double prevPM = currentPM;
527  G4double prevCII = currentCII;
528  G4double prevIN = currentIN;
529 
530  // loop over all (photon energy, intensity)
531  // pairs stored for this material
532 
533  for (size_t ii = 1;
534  ii < theFastLightVector->GetVectorLength();
535  ++ii)
536  {
537  currentPM = theFastLightVector->Energy(ii);
538  currentIN = (*theFastLightVector)[ii];
539 
540  currentCII = 0.5 * (prevIN + currentIN);
541 
542  currentCII = prevCII +
543  (currentPM - prevPM) * currentCII;
544 
545  aPhysicsOrderedFreeVector->
546  InsertValues(currentPM, currentCII);
547 
548  prevPM = currentPM;
549  prevCII = currentCII;
550  prevIN = currentIN;
551  }
552 
553  }
554  }
555 
556  G4MaterialPropertyVector* theSlowLightVector =
557  aMaterialPropertiesTable->GetProperty(kSLOWCOMPONENT);
558 
559  if (theSlowLightVector) {
560 
561  // Retrieve the first intensity point in vector
562  // of (photon energy, intensity) pairs
563 
564  G4double currentIN = (*theSlowLightVector)[0];
565 
566  if (currentIN >= 0.0) {
567 
568  // Create first (photon energy, Scintillation
569  // Integral pair
570 
571  G4double currentPM = theSlowLightVector->Energy(0);
572 
573  G4double currentCII = 0.0;
574 
575  bPhysicsOrderedFreeVector->
576  InsertValues(currentPM , currentCII);
577 
578  // Set previous values to current ones prior to loop
579 
580  G4double prevPM = currentPM;
581  G4double prevCII = currentCII;
582  G4double prevIN = currentIN;
583 
584  // loop over all (photon energy, intensity)
585  // pairs stored for this material
586 
587  for (size_t ii = 1;
588  ii < theSlowLightVector->GetVectorLength();
589  ++ii)
590  {
591  currentPM = theSlowLightVector->Energy(ii);
592  currentIN = (*theSlowLightVector)[ii];
593 
594  currentCII = 0.5 * (prevIN + currentIN);
595 
596  currentCII = prevCII +
597  (currentPM - prevPM) * currentCII;
598 
599  bPhysicsOrderedFreeVector->
600  InsertValues(currentPM, currentCII);
601 
602  prevPM = currentPM;
603  prevCII = currentCII;
604  prevIN = currentIN;
605  }
606 
607  }
608  }
609  }
610 
611  // The scintillation integral(s) for a given material
612  // will be inserted in the table(s) according to the
613  // position of the material in the material table.
614 
615  fFastIntegralTable->insertAt(i,aPhysicsOrderedFreeVector);
616  fSlowIntegralTable->insertAt(i,bPhysicsOrderedFreeVector);
617 
618  }
619 }
620 
622 {
623  if (fEmSaturation) {
624  G4Exception("G4Scintillation::SetScintillationByParticleType", "Scint02",
625  JustWarning, "Redefinition: Birks Saturation is replaced by ScintillationByParticleType!");
627  }
628  fScintillationByParticleType = scintType;
629 }
630 
631 // GetMeanFreePath
632 // ---------------
633 //
634 
636  G4double ,
638 {
639  *condition = StronglyForced;
640 
641  return DBL_MAX;
642 
643 }
644 
645 // GetMeanLifeTime
646 // ---------------
647 //
648 
651 {
652  *condition = Forced;
653 
654  return DBL_MAX;
655 
656 }
657 
659 {
660 // tau1: rise time and tau2: decay time
661 
662  // Loop checking, 07-Aug-2015, Vladimir Ivanchenko
663  while(1) {
664  // two random numbers
665  G4double ran1 = G4UniformRand();
666  G4double ran2 = G4UniformRand();
667  //
668  // exponential distribution as envelope function: very efficient
669  //
670  G4double d = (tau1+tau2)/tau2;
671  // make sure the envelope function is
672  // always larger than the bi-exponential
673  G4double t = -1.0*tau2*std::log(1-ran1);
674  G4double gg = d*single_exp(t,tau2);
675  if (ran2 <= bi_exp(t,tau1,tau2)/gg) return t;
676  }
677  return -1.0;
678 }
679 
682 {
684  // Get the scintillation yield vector //
686 
688 
689  G4MaterialPropertyVector *Scint_Yield_Vector = nullptr;
690 
691  G4MaterialPropertiesTable *aMaterialPropertiesTable
693 
694  // Get the G4MaterialPropertyVector containing the scintillation
695  // yield as a function of the energy deposited and particle type
696 
697  // Protons
698  if(pDef==G4Proton::ProtonDefinition())
699  Scint_Yield_Vector = aMaterialPropertiesTable->
700  GetProperty(kPROTONSCINTILLATIONYIELD);
701 
702  // Deuterons
703  else if(pDef==G4Deuteron::DeuteronDefinition())
704  Scint_Yield_Vector = aMaterialPropertiesTable->
705  GetProperty(kDEUTERONSCINTILLATIONYIELD);
706 
707  // Tritons
708  else if(pDef==G4Triton::TritonDefinition())
709  Scint_Yield_Vector = aMaterialPropertiesTable->
710  GetProperty(kTRITONSCINTILLATIONYIELD);
711 
712  // Alphas
713  else if(pDef==G4Alpha::AlphaDefinition())
714  Scint_Yield_Vector = aMaterialPropertiesTable->
715  GetProperty(kALPHASCINTILLATIONYIELD);
716 
717  // Ions (particles derived from G4VIon and G4Ions) and recoil ions
718  // below the production cut from neutrons after hElastic
719  else if(pDef->GetParticleType()== "nucleus" ||
721  Scint_Yield_Vector = aMaterialPropertiesTable->
722  GetProperty(kIONSCINTILLATIONYIELD);
723 
724  // Electrons (must also account for shell-binding energy
725  // attributed to gamma from standard photoelectric effect)
726  else if(pDef==G4Electron::ElectronDefinition() ||
727  pDef==G4Gamma::GammaDefinition())
728  Scint_Yield_Vector = aMaterialPropertiesTable->
729  GetProperty(kELECTRONSCINTILLATIONYIELD);
730 
731  // Default for particles not enumerated/listed above
732  else
733  Scint_Yield_Vector = aMaterialPropertiesTable->
734  GetProperty(kELECTRONSCINTILLATIONYIELD);
735 
736  // If the user has specified none of the above particles then the
737  // default is the electron scintillation yield
738  if(!Scint_Yield_Vector)
739  Scint_Yield_Vector = aMaterialPropertiesTable->
740  GetProperty(kELECTRONSCINTILLATIONYIELD);
741 
742  // Throw an exception if no scintillation yield vector is found
743  if (!Scint_Yield_Vector) {
745  ed << "\nG4Scintillation::PostStepDoIt(): "
746  << "Request for scintillation yield for energy deposit and particle\n"
747  << "type without correct entry in MaterialPropertiesTable.\n"
748  << "ScintillationByParticleType requires at minimum that \n"
749  << "ELECTRONSCINTILLATIONYIELD is set by the user\n"
750  << G4endl;
751  G4String comments = "Missing MaterialPropertiesTable entry - No correct entry in MaterialPropertiesTable";
752  G4Exception("G4Scintillation::PostStepDoIt","Scint01",
753  FatalException,ed,comments);
754  }
755 
757  // Calculate the scintillation light //
759  // To account for potential nonlinearity and scintillation photon
760  // density along the track, light (L) is produced according to:
761  //
762  // L_currentStep = L(PreStepKE) - L(PreStepKE - EDep)
763 
764  G4double ScintillationYield = 0.;
765 
766  G4double StepEnergyDeposit = aStep.GetTotalEnergyDeposit();
767  G4double PreStepKineticEnergy = aStep.GetPreStepPoint()->GetKineticEnergy();
768 
769  if(PreStepKineticEnergy <= Scint_Yield_Vector->GetMaxEnergy()){
770  G4double Yield1 = Scint_Yield_Vector->Value(PreStepKineticEnergy);
771  G4double Yield2 = Scint_Yield_Vector->
772  Value(PreStepKineticEnergy - StepEnergyDeposit);
773  ScintillationYield = Yield1 - Yield2;
774  } else {
776  ed << "\nG4Scintillation::GetScintillationYieldByParticleType(): Request\n"
777  << "for scintillation light yield above the available energy range\n"
778  << "specifed in G4MaterialPropertiesTable. A linear interpolation\n"
779  << "will be performed to compute the scintillation light yield using\n"
780  << "(L_max / E_max) as the photon yield per unit energy."
781  << G4endl;
782  G4String cmt = "\nScintillation yield may be unphysical!\n";
783  G4Exception("G4Scintillation::GetScintillationYieldByParticleType()",
784  "Scint03", JustWarning, ed, cmt);
785 
786  G4double LinearYield = Scint_Yield_Vector->GetMaxValue()
787  / Scint_Yield_Vector->GetMaxEnergy();
788 
789  // Units: [# scintillation photons]
790  ScintillationYield = LinearYield * StepEnergyDeposit;
791  }
792 
793 #ifdef G4DEBUG_SCINTILLATION
794 
795  // Increment track aggregators
796  ScintTrackYield += ScintillationYield;
797  ScintTrackEDep += StepEnergyDeposit;
798 
799  G4cout << "\n--- G4Scintillation::GetScintillationYieldByParticleType() ---\n"
800  << "--\n"
801  << "-- Name = " << aTrack.GetParticleDefinition()->GetParticleName() << "\n"
802  << "-- TrackID = " << aTrack.GetTrackID() << "\n"
803  << "-- ParentID = " << aTrack.GetParentID() << "\n"
804  << "-- Current KE = " << aTrack.GetKineticEnergy()/MeV << " MeV\n"
805  << "-- Step EDep = " << aStep.GetTotalEnergyDeposit()/MeV << " MeV\n"
806  << "-- Track EDep = " << ScintTrackEDep/MeV << " MeV\n"
807  << "-- Vertex KE = " << aTrack.GetVertexKineticEnergy()/MeV << " MeV\n"
808  << "-- Step yield = " << ScintillationYield << " photons\n"
809  << "-- Track yield = " << ScintTrackYield << " photons\n"
810  << G4endl;
811 
812  // The track has terminated within or has left the scintillator volume
813  if( (aTrack.GetTrackStatus() == fStopButAlive) or
814  (aStep.GetPostStepPoint()->GetStepStatus() == fGeomBoundary) ){
815 
816  // Reset aggregators for the next track
817  ScintTrackEDep = 0.;
818  ScintTrackYield = 0.;
819  }
820 
821 #endif
822 
823  return ScintillationYield;
824 }
825 
827 {
828  if (fFastIntegralTable) {
829  G4int PhysicsTableSize = fFastIntegralTable->entries();
831 
832  for (G4int i = 0 ; i < PhysicsTableSize ; i++ )
833  {
835  v->DumpValues();
836  }
837  }
838 
839  if (fSlowIntegralTable) {
840  G4int PhysicsTableSize = fSlowIntegralTable->entries();
842 
843  for (G4int i = 0 ; i < PhysicsTableSize ; i++ )
844  {
846  v->DumpValues();
847  }
848  }
849 }