Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
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 // $Id: G4Scintillation.cc 108423 2018-02-13 11:18:13Z gcosmo $
27 //
29 // Scintillation Light Class Implementation
31 //
32 // File: G4Scintillation.cc
33 // Description: RestDiscrete Process - Generation of Scintillation Photons
34 // Version: 1.0
35 // Created: 1998-11-07
36 // Author: Peter Gumplinger
37 // Updated: 2010-10-20 Allow the scintillation yield to be a function
38 // of energy deposited by particle type
39 // Thanks to Zach Hartwig (Department of Nuclear
40 // Science and Engineeering - MIT)
41 // 2010-09-22 by Peter Gumplinger
42 // > scintillation rise time included, thanks to
43 // > Martin Goettlich/DESY
44 // 2005-08-17 by Peter Gumplinger
45 // > change variable name MeanNumPhotons -> MeanNumberOfPhotons
46 // 2005-07-28 by Peter Gumplinger
47 // > add G4ProcessType to constructor
48 // 2004-08-05 by Peter Gumplinger
49 // > changed StronglyForced back to Forced in GetMeanLifeTime
50 // 2002-11-21 by Peter Gumplinger
51 // > change to use G4Poisson for small MeanNumberOfPhotons
52 // 2002-11-07 by Peter Gumplinger
53 // > now allow for fast and slow scintillation component
54 // 2002-11-05 by Peter Gumplinger
55 // > now use scintillation constants from G4Material
56 // 2002-05-09 by Peter Gumplinger
57 // > use only the PostStepPoint location for the origin of
58 // scintillation photons when energy is lost to the medium
59 // by a neutral particle
60 // 2000-09-18 by Peter Gumplinger
61 // > change: aSecondaryPosition=x0+rand*aStep.GetDeltaPosition();
62 // aSecondaryTrack->SetTouchable(0);
63 // 2001-09-17, migration of Materials to pure STL (mma)
64 // 2003-06-03, V.Ivanchenko fix compilation warnings
65 //
66 // mail: gum@triumf.ca
67 //
69 
70 #include "G4ios.hh"
71 #include "globals.hh"
72 #include "G4PhysicalConstants.hh"
73 #include "G4SystemOfUnits.hh"
74 #include "G4ParticleTypes.hh"
75 #include "G4EmProcessSubType.hh"
77 
78 #include "G4Scintillation.hh"
79 
81 // Class Implementation
83 
85  // static data members
87 
88 
89 //G4bool G4Scintillation::fTrackSecondariesFirst = false;
90 //G4bool G4Scintillation::fFiniteRiseTime = false;
91 //G4double G4Scintillation::fYieldFactor = 1.0;
92 //G4double G4Scintillation::fExcitationRatio = 1.0;
93 //G4bool G4Scintillation::fScintillationByParticleType = false;
94 //G4bool G4Scintillation::fScintillationTrackInfo = false;
95 //G4bool G4Scintillation::fStackingFlag = true;
96 //G4EmSaturation* G4Scintillation::fEmSaturation = NULL;
97 
99  // Operators
101 
102 // G4Scintillation::operator=(const G4Scintillation &right)
103 // {
104 // }
105 
107  // Constructors
109 
111  G4ProcessType type)
112  : G4VRestDiscreteProcess(processName, type) ,
113  fTrackSecondariesFirst(false),
114  fFiniteRiseTime(false),
115  fYieldFactor(1.0),
116  fExcitationRatio(1.0),
117  fScintillationByParticleType(false),
118  fScintillationTrackInfo(false),
119  fStackingFlag(true),
120  fNumPhotons(0),
121  fEmSaturation(nullptr)
122 {
124 
125 #ifdef G4DEBUG_SCINTILLATION
126  ScintTrackEDep = 0.;
127  ScintTrackYield = 0.;
128 #endif
129 
130  fFastIntegralTable = nullptr;
131  fSlowIntegralTable = nullptr;
132 
133  if (verboseLevel>0) {
134  G4cout << GetProcessName() << " is created " << G4endl;
135  }
136 }
137 
139  // Destructors
141 
143 {
144  if (fFastIntegralTable != nullptr) {
146  delete fFastIntegralTable;
147  }
148  if (fSlowIntegralTable != nullptr) {
150  delete fSlowIntegralTable;
151  }
152 }
153 
155  // Methods
157 
159 {
160  if (aParticleType.GetParticleName() == "opticalphoton") return false;
161  if (aParticleType.IsShortLived()) return false;
162 
163  return true;
164 }
165 
167 {
168  if (fFastIntegralTable != nullptr) {
170  delete fFastIntegralTable;
171  fFastIntegralTable = nullptr;
172  }
173  if (fSlowIntegralTable != nullptr) {
175  delete fSlowIntegralTable;
176  fSlowIntegralTable = nullptr;
177  }
179 }
180 
181 
182 // AtRestDoIt
183 // ----------
184 //
186 G4Scintillation::AtRestDoIt(const G4Track& aTrack, const G4Step& aStep)
187 
188 // This routine simply calls the equivalent PostStepDoIt since all the
189 // necessary information resides in aStep.GetTotalEnergyDeposit()
190 
191 {
192  return G4Scintillation::PostStepDoIt(aTrack, aStep);
193 }
194 
195 // PostStepDoIt
196 // -------------
197 //
199 G4Scintillation::PostStepDoIt(const G4Track& aTrack, const G4Step& aStep)
200 
201 // This routine is called for each tracking step of a charged particle
202 // in a scintillator. A Poisson/Gauss-distributed number of photons is
203 // generated according to the scintillation yield formula, distributed
204 // evenly along the track segment and uniformly into 4pi.
205 
206 {
207  aParticleChange.Initialize(aTrack);
208 
209  const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
210  const G4Material* aMaterial = aTrack.GetMaterial();
211 
212  G4StepPoint* pPreStepPoint = aStep.GetPreStepPoint();
213  G4StepPoint* pPostStepPoint = aStep.GetPostStepPoint();
214 
215  G4ThreeVector x0 = pPreStepPoint->GetPosition();
216  G4ThreeVector p0 = aStep.GetDeltaPosition().unit();
217  G4double t0 = pPreStepPoint->GetGlobalTime();
218 
219  G4double TotalEnergyDeposit = aStep.GetTotalEnergyDeposit();
220 
221  G4MaterialPropertiesTable* aMaterialPropertiesTable =
222  aMaterial->GetMaterialPropertiesTable();
223  if (!aMaterialPropertiesTable)
224  return G4VRestDiscreteProcess::PostStepDoIt(aTrack, aStep);
225 
226  G4MaterialPropertyVector* Fast_Intensity =
227  aMaterialPropertiesTable->GetProperty(kFASTCOMPONENT);
228  G4MaterialPropertyVector* Slow_Intensity =
229  aMaterialPropertiesTable->GetProperty(kSLOWCOMPONENT);
230 
231  if (!Fast_Intensity && !Slow_Intensity )
232  return G4VRestDiscreteProcess::PostStepDoIt(aTrack, aStep);
233 
234  G4int nscnt = 1;
235  if (Fast_Intensity && Slow_Intensity) nscnt = 2;
236 
237  G4double ScintillationYield = 0.;
238 
239  // Scintillation depends on particle type, energy deposited
241 
242  ScintillationYield =
244 
245  // The default linear scintillation process
246  } else {
247 
248  ScintillationYield = aMaterialPropertiesTable->
249  GetConstProperty(kSCINTILLATIONYIELD);
250 
251  // Units: [# scintillation photons / MeV]
252  ScintillationYield *= fYieldFactor;
253  }
254 
255  G4double ResolutionScale = aMaterialPropertiesTable->
256  GetConstProperty(kRESOLUTIONSCALE);
257 
258  // Birks law saturation:
259 
260  //G4double constBirks = 0.0;
261 
262  //constBirks = aMaterial->GetIonisation()->GetBirksConstant();
263 
264  G4double MeanNumberOfPhotons;
265 
266  // Birk's correction via fEmSaturation and specifying scintillation by
267  // by particle type are physically mutually exclusive
268 
270  MeanNumberOfPhotons = ScintillationYield;
271  else if (fEmSaturation)
272  MeanNumberOfPhotons = ScintillationYield*
274  else
275  MeanNumberOfPhotons = ScintillationYield*TotalEnergyDeposit;
276 
277  if (MeanNumberOfPhotons > 10.)
278  {
279  G4double sigma = ResolutionScale * std::sqrt(MeanNumberOfPhotons);
280  fNumPhotons=G4int(G4RandGauss::shoot(MeanNumberOfPhotons,sigma)+0.5);
281  }
282  else
283  {
284  fNumPhotons = G4int(G4Poisson(MeanNumberOfPhotons));
285  }
286 
287  if ( fNumPhotons <= 0 || !fStackingFlag )
288  {
289  // return unchanged particle and no secondaries
290 
292 
293  return G4VRestDiscreteProcess::PostStepDoIt(aTrack, aStep);
294  }
295 
297 
299 
301  if (aTrack.GetTrackStatus() == fAlive )
303  }
304 
306 
307  G4int materialIndex = aMaterial->GetIndex();
308 
309  // Retrieve the Scintillation Integral for this material
310  // new G4PhysicsOrderedFreeVector allocated to hold CII's
311 
312  G4int Num = fNumPhotons;
313 
314  for (G4int scnt = 1; scnt <= nscnt; scnt++) {
315 
316  G4double ScintillationTime = 0.*ns;
317  G4double ScintillationRiseTime = 0.*ns;
318  G4PhysicsOrderedFreeVector* ScintillationIntegral = nullptr;
319  G4ScintillationType ScintillationType = Slow;
320 
321  if (scnt == 1) {
322  if (nscnt == 1) {
323  if (Fast_Intensity) {
324  ScintillationTime = aMaterialPropertiesTable->
325  GetConstProperty(kFASTTIMECONSTANT);
326  if (fFiniteRiseTime) {
327  ScintillationRiseTime = aMaterialPropertiesTable->
328  GetConstProperty(kFASTSCINTILLATIONRISETIME);
329  }
330  ScintillationType = Fast;
331  ScintillationIntegral =
333  ((*fFastIntegralTable)(materialIndex));
334  }
335  if (Slow_Intensity) {
336  ScintillationTime = aMaterialPropertiesTable->
337  GetConstProperty(kSLOWTIMECONSTANT);
338  if (fFiniteRiseTime) {
339  ScintillationRiseTime = aMaterialPropertiesTable->
340  GetConstProperty(kSLOWSCINTILLATIONRISETIME);
341  }
342  ScintillationType = Slow;
343  ScintillationIntegral =
345  ((*fSlowIntegralTable)(materialIndex));
346  }
347  }
348  else {
349  G4double yieldRatio = aMaterialPropertiesTable->
350  GetConstProperty(kYIELDRATIO);
351  if ( fExcitationRatio == 1.0 || fExcitationRatio == 0.0) {
352  Num = G4int (std::min(yieldRatio,1.0) * fNumPhotons);
353  }
354  else {
356  }
357  ScintillationTime = aMaterialPropertiesTable->
358  GetConstProperty(kFASTTIMECONSTANT);
359  if (fFiniteRiseTime) {
360  ScintillationRiseTime = aMaterialPropertiesTable->
361  GetConstProperty(kFASTSCINTILLATIONRISETIME);
362  }
363  ScintillationType = Fast;
364  ScintillationIntegral =
366  ((*fFastIntegralTable)(materialIndex));
367  }
368  }
369  else {
370  Num = fNumPhotons - Num;
371  ScintillationTime = aMaterialPropertiesTable->
372  GetConstProperty(kSLOWTIMECONSTANT);
373  if (fFiniteRiseTime) {
374  ScintillationRiseTime = aMaterialPropertiesTable->
375  GetConstProperty(kSLOWSCINTILLATIONRISETIME);
376  }
377  ScintillationType = Slow;
378  ScintillationIntegral =
380  ((*fSlowIntegralTable)(materialIndex));
381  }
382 
383  if (!ScintillationIntegral) continue;
384 
385  // Max Scintillation Integral
386 
387  G4double CIImax = ScintillationIntegral->GetMaxValue();
388 
389  for (G4int i = 0; i < Num; i++) {
390 
391  // Determine photon energy
392 
393  G4double CIIvalue = G4UniformRand()*CIImax;
394  G4double sampledEnergy =
395  ScintillationIntegral->GetEnergy(CIIvalue);
396 
397  if (verboseLevel>1) {
398  G4cout << "sampledEnergy = " << sampledEnergy << G4endl;
399  G4cout << "CIIvalue = " << CIIvalue << G4endl;
400  }
401 
402  // Generate random photon direction
403 
404  G4double cost = 1. - 2.*G4UniformRand();
405  G4double sint = std::sqrt((1.-cost)*(1.+cost));
406 
407  G4double phi = twopi*G4UniformRand();
408  G4double sinp = std::sin(phi);
409  G4double cosp = std::cos(phi);
410 
411  G4double px = sint*cosp;
412  G4double py = sint*sinp;
413  G4double pz = cost;
414 
415  // Create photon momentum direction vector
416 
417  G4ParticleMomentum photonMomentum(px, py, pz);
418 
419  // Determine polarization of new photon
420 
421  G4double sx = cost*cosp;
422  G4double sy = cost*sinp;
423  G4double sz = -sint;
424 
425  G4ThreeVector photonPolarization(sx, sy, sz);
426 
427  G4ThreeVector perp = photonMomentum.cross(photonPolarization);
428 
429  phi = twopi*G4UniformRand();
430  sinp = std::sin(phi);
431  cosp = std::cos(phi);
432 
433  photonPolarization = cosp * photonPolarization + sinp * perp;
434 
435  photonPolarization = photonPolarization.unit();
436 
437  // Generate a new photon:
438 
439  G4DynamicParticle* aScintillationPhoton =
441  photonMomentum);
442  aScintillationPhoton->SetPolarization
443  (photonPolarization.x(),
444  photonPolarization.y(),
445  photonPolarization.z());
446 
447  aScintillationPhoton->SetKineticEnergy(sampledEnergy);
448 
449  // Generate new G4Track object:
450 
451  G4double rand;
452 
453  if (aParticle->GetDefinition()->GetPDGCharge() != 0) {
454  rand = G4UniformRand();
455  } else {
456  rand = 1.0;
457  }
458 
459  G4double delta = rand * aStep.GetStepLength();
460  G4double deltaTime = delta / (pPreStepPoint->GetVelocity()+
461  rand*(pPostStepPoint->GetVelocity()-
462  pPreStepPoint->GetVelocity())/2.);
463 
464  // emission time distribution
465  if (ScintillationRiseTime==0.0) {
466  deltaTime = deltaTime -
467  ScintillationTime * std::log( G4UniformRand() );
468  } else {
469  deltaTime = deltaTime +
470  sample_time(ScintillationRiseTime, ScintillationTime);
471  }
472 
473  G4double aSecondaryTime = t0 + deltaTime;
474 
475  G4ThreeVector aSecondaryPosition =
476  x0 + rand * aStep.GetDeltaPosition();
477 
478  G4Track* aSecondaryTrack = new G4Track(aScintillationPhoton,
479  aSecondaryTime,
480  aSecondaryPosition);
481 
482  aSecondaryTrack->SetTouchableHandle(
484  // aSecondaryTrack->SetTouchableHandle((G4VTouchable*)0);
485 
486  aSecondaryTrack->SetParentID(aTrack.GetTrackID());
487 
488  if (fScintillationTrackInfo) aSecondaryTrack->
489  SetUserInformation(new G4ScintillationTrackInformation(ScintillationType));
490 
491  aParticleChange.AddSecondary(aSecondaryTrack);
492 
493  }
494  }
495 
496  if (verboseLevel>0) {
497  G4cout << "\n Exiting from G4Scintillation::DoIt -- NumberOfSecondaries = "
499  }
500 
501  return G4VRestDiscreteProcess::PostStepDoIt(aTrack, aStep);
502 }
503 
504 // BuildThePhysicsTable for the scintillation process
505 // --------------------------------------------------
506 //
507 
509 {
510  if (fFastIntegralTable && fSlowIntegralTable) return;
511 
512  const G4MaterialTable* theMaterialTable =
514  G4int numOfMaterials = G4Material::GetNumberOfMaterials();
515 
516  // create new physics table
517 
519  new G4PhysicsTable(numOfMaterials);
521  new G4PhysicsTable(numOfMaterials);
522 
523  // loop for materials
524 
525  for (G4int i=0 ; i < numOfMaterials; i++)
526  {
527  G4PhysicsOrderedFreeVector* aPhysicsOrderedFreeVector =
529  G4PhysicsOrderedFreeVector* bPhysicsOrderedFreeVector =
531 
532  // Retrieve vector of scintillation wavelength intensity for
533  // the material from the material's optical properties table.
534 
535  G4Material* aMaterial = (*theMaterialTable)[i];
536 
537  G4MaterialPropertiesTable* aMaterialPropertiesTable =
538  aMaterial->GetMaterialPropertiesTable();
539 
540  if (aMaterialPropertiesTable) {
541 
542  G4MaterialPropertyVector* theFastLightVector =
543  aMaterialPropertiesTable->GetProperty(kFASTCOMPONENT);
544 
545  if (theFastLightVector) {
546 
547  // Retrieve the first intensity point in vector
548  // of (photon energy, intensity) pairs
549 
550  G4double currentIN = (*theFastLightVector)[0];
551 
552  if (currentIN >= 0.0) {
553 
554  // Create first (photon energy, Scintillation
555  // Integral pair
556 
557  G4double currentPM = theFastLightVector->Energy(0);
558 
559  G4double currentCII = 0.0;
560 
561  aPhysicsOrderedFreeVector->
562  InsertValues(currentPM , currentCII);
563 
564  // Set previous values to current ones prior to loop
565 
566  G4double prevPM = currentPM;
567  G4double prevCII = currentCII;
568  G4double prevIN = currentIN;
569 
570  // loop over all (photon energy, intensity)
571  // pairs stored for this material
572 
573  for (size_t ii = 1;
574  ii < theFastLightVector->GetVectorLength();
575  ++ii)
576  {
577  currentPM = theFastLightVector->Energy(ii);
578  currentIN = (*theFastLightVector)[ii];
579 
580  currentCII = 0.5 * (prevIN + currentIN);
581 
582  currentCII = prevCII +
583  (currentPM - prevPM) * currentCII;
584 
585  aPhysicsOrderedFreeVector->
586  InsertValues(currentPM, currentCII);
587 
588  prevPM = currentPM;
589  prevCII = currentCII;
590  prevIN = currentIN;
591  }
592 
593  }
594  }
595 
596  G4MaterialPropertyVector* theSlowLightVector =
597  aMaterialPropertiesTable->GetProperty(kSLOWCOMPONENT);
598 
599  if (theSlowLightVector) {
600 
601  // Retrieve the first intensity point in vector
602  // of (photon energy, intensity) pairs
603 
604  G4double currentIN = (*theSlowLightVector)[0];
605 
606  if (currentIN >= 0.0) {
607 
608  // Create first (photon energy, Scintillation
609  // Integral pair
610 
611  G4double currentPM = theSlowLightVector->Energy(0);
612 
613  G4double currentCII = 0.0;
614 
615  bPhysicsOrderedFreeVector->
616  InsertValues(currentPM , currentCII);
617 
618  // Set previous values to current ones prior to loop
619 
620  G4double prevPM = currentPM;
621  G4double prevCII = currentCII;
622  G4double prevIN = currentIN;
623 
624  // loop over all (photon energy, intensity)
625  // pairs stored for this material
626 
627  for (size_t ii = 1;
628  ii < theSlowLightVector->GetVectorLength();
629  ++ii)
630  {
631  currentPM = theSlowLightVector->Energy(ii);
632  currentIN = (*theSlowLightVector)[ii];
633 
634  currentCII = 0.5 * (prevIN + currentIN);
635 
636  currentCII = prevCII +
637  (currentPM - prevPM) * currentCII;
638 
639  bPhysicsOrderedFreeVector->
640  InsertValues(currentPM, currentCII);
641 
642  prevPM = currentPM;
643  prevCII = currentCII;
644  prevIN = currentIN;
645  }
646 
647  }
648  }
649  }
650 
651  // The scintillation integral(s) for a given material
652  // will be inserted in the table(s) according to the
653  // position of the material in the material table.
654 
655  fFastIntegralTable->insertAt(i,aPhysicsOrderedFreeVector);
656  fSlowIntegralTable->insertAt(i,bPhysicsOrderedFreeVector);
657 
658  }
659 }
660 
662 {
663  if (fEmSaturation) {
664  G4Exception("G4Scintillation::SetScintillationByParticleType", "Scint02",
665  JustWarning, "Redefinition: Birks Saturation is replaced by ScintillationByParticleType!");
667  }
668  fScintillationByParticleType = scintType;
669 }
670 
671 // GetMeanFreePath
672 // ---------------
673 //
674 
676  G4double ,
678 {
679  *condition = StronglyForced;
680 
681  return DBL_MAX;
682 
683 }
684 
685 // GetMeanLifeTime
686 // ---------------
687 //
688 
691 {
692  *condition = Forced;
693 
694  return DBL_MAX;
695 
696 }
697 
699 {
700 // tau1: rise time and tau2: decay time
701 
702  // Loop checking, 07-Aug-2015, Vladimir Ivanchenko
703  while(1) {
704  // two random numbers
705  G4double ran1 = G4UniformRand();
706  G4double ran2 = G4UniformRand();
707  //
708  // exponential distribution as envelope function: very efficient
709  //
710  G4double d = (tau1+tau2)/tau2;
711  // make sure the envelope function is
712  // always larger than the bi-exponential
713  G4double t = -1.0*tau2*std::log(1-ran1);
714  G4double gg = d*single_exp(t,tau2);
715  if (ran2 <= bi_exp(t,tau1,tau2)/gg) return t;
716  }
717  return -1.0;
718 }
719 
722 {
724  // Get the scintillation yield vector //
726 
728 
729  G4MaterialPropertyVector *Scint_Yield_Vector = nullptr;
730 
731  G4MaterialPropertiesTable *aMaterialPropertiesTable
733 
734  // Get the G4MaterialPropertyVector containing the scintillation
735  // yield as a function of the energy deposited and particle type
736 
737  // Protons
738  if(pDef==G4Proton::ProtonDefinition())
739  Scint_Yield_Vector = aMaterialPropertiesTable->
740  GetProperty(kPROTONSCINTILLATIONYIELD);
741 
742  // Deuterons
743  else if(pDef==G4Deuteron::DeuteronDefinition())
744  Scint_Yield_Vector = aMaterialPropertiesTable->
745  GetProperty(kDEUTERONSCINTILLATIONYIELD);
746 
747  // Tritons
748  else if(pDef==G4Triton::TritonDefinition())
749  Scint_Yield_Vector = aMaterialPropertiesTable->
750  GetProperty(kTRITONSCINTILLATIONYIELD);
751 
752  // Alphas
753  else if(pDef==G4Alpha::AlphaDefinition())
754  Scint_Yield_Vector = aMaterialPropertiesTable->
755  GetProperty(kALPHASCINTILLATIONYIELD);
756 
757  // Ions (particles derived from G4VIon and G4Ions) and recoil ions
758  // below the production cut from neutrons after hElastic
759  else if(pDef->GetParticleType()== "nucleus" ||
761  Scint_Yield_Vector = aMaterialPropertiesTable->
762  GetProperty(kIONSCINTILLATIONYIELD);
763 
764  // Electrons (must also account for shell-binding energy
765  // attributed to gamma from standard photoelectric effect)
766  else if(pDef==G4Electron::ElectronDefinition() ||
767  pDef==G4Gamma::GammaDefinition())
768  Scint_Yield_Vector = aMaterialPropertiesTable->
769  GetProperty(kELECTRONSCINTILLATIONYIELD);
770 
771  // Default for particles not enumerated/listed above
772  else
773  Scint_Yield_Vector = aMaterialPropertiesTable->
774  GetProperty(kELECTRONSCINTILLATIONYIELD);
775 
776  // If the user has specified none of the above particles then the
777  // default is the electron scintillation yield
778  if(!Scint_Yield_Vector)
779  Scint_Yield_Vector = aMaterialPropertiesTable->
780  GetProperty(kELECTRONSCINTILLATIONYIELD);
781 
782  // Throw an exception if no scintillation yield vector is found
783  if (!Scint_Yield_Vector) {
785  ed << "\nG4Scintillation::PostStepDoIt(): "
786  << "Request for scintillation yield for energy deposit and particle\n"
787  << "type without correct entry in MaterialPropertiesTable.\n"
788  << "ScintillationByParticleType requires at minimum that \n"
789  << "ELECTRONSCINTILLATIONYIELD is set by the user\n"
790  << G4endl;
791  G4String comments = "Missing MaterialPropertiesTable entry - No correct entry in MaterialPropertiesTable";
792  G4Exception("G4Scintillation::PostStepDoIt","Scint01",
793  FatalException,ed,comments);
794  }
795 
797  // Calculate the scintillation light //
799  // To account for potential nonlinearity and scintillation photon
800  // density along the track, light (L) is produced according to:
801  //
802  // L_currentStep = L(PreStepKE) - L(PreStepKE - EDep)
803 
804  G4double ScintillationYield = 0.;
805 
806  G4double StepEnergyDeposit = aStep.GetTotalEnergyDeposit();
807  G4double PreStepKineticEnergy = aStep.GetPreStepPoint()->GetKineticEnergy();
808 
809  if(PreStepKineticEnergy <= Scint_Yield_Vector->GetMaxEnergy()){
810  G4double Yield1 = Scint_Yield_Vector->Value(PreStepKineticEnergy);
811  G4double Yield2 = Scint_Yield_Vector->
812  Value(PreStepKineticEnergy - StepEnergyDeposit);
813  ScintillationYield = Yield1 - Yield2;
814  } else {
816  ed << "\nG4Scintillation::GetScintillationYieldByParticleType(): Request\n"
817  << "for scintillation light yield above the available energy range\n"
818  << "specifed in G4MaterialPropertiesTable. A linear interpolation\n"
819  << "will be performed to compute the scintillation light yield using\n"
820  << "(L_max / E_max) as the photon yield per unit energy."
821  << G4endl;
822  G4String cmt = "\nScintillation yield may be unphysical!\n";
823  G4Exception("G4Scintillation::GetScintillationYieldByParticleType()",
824  "Scint03", JustWarning, ed, cmt);
825 
826  G4double LinearYield = Scint_Yield_Vector->GetMaxValue()
827  / Scint_Yield_Vector->GetMaxEnergy();
828 
829  // Units: [# scintillation photons]
830  ScintillationYield = LinearYield * StepEnergyDeposit;
831  }
832 
833 #ifdef G4DEBUG_SCINTILLATION
834 
835  // Increment track aggregators
836  ScintTrackYield += ScintillationYield;
837  ScintTrackEDep += StepEnergyDeposit;
838 
839  G4cout << "\n--- G4Scintillation::GetScintillationYieldByParticleType() ---\n"
840  << "--\n"
841  << "-- Name = " << aTrack.GetParticleDefinition()->GetParticleName() << "\n"
842  << "-- TrackID = " << aTrack.GetTrackID() << "\n"
843  << "-- ParentID = " << aTrack.GetParentID() << "\n"
844  << "-- Current KE = " << aTrack.GetKineticEnergy()/MeV << " MeV\n"
845  << "-- Step EDep = " << aStep.GetTotalEnergyDeposit()/MeV << " MeV\n"
846  << "-- Track EDep = " << ScintTrackEDep/MeV << " MeV\n"
847  << "-- Vertex KE = " << aTrack.GetVertexKineticEnergy()/MeV << " MeV\n"
848  << "-- Step yield = " << ScintillationYield << " photons\n"
849  << "-- Track yield = " << ScintTrackYield << " photons\n"
850  << G4endl;
851 
852  // The track has terminated within or has left the scintillator volume
853  if( (aTrack.GetTrackStatus() == fStopButAlive) or
854  (aStep.GetPostStepPoint()->GetStepStatus() == fGeomBoundary) ){
855 
856  // Reset aggregators for the next track
857  ScintTrackEDep = 0.;
858  ScintTrackYield = 0.;
859  }
860 
861 #endif
862 
863  return ScintillationYield;
864 }
865 
867 {
868  if (fFastIntegralTable) {
869  G4int PhysicsTableSize = fFastIntegralTable->entries();
871 
872  for (G4int i = 0 ; i < PhysicsTableSize ; i++ )
873  {
875  v->DumpValues();
876  }
877  }
878 
879  if (fSlowIntegralTable) {
880  G4int PhysicsTableSize = fSlowIntegralTable->entries();
882 
883  for (G4int i = 0 ; i < PhysicsTableSize ; i++ )
884  {
886  v->DumpValues();
887  }
888  }
889 }
G4double Energy(size_t index) const
G4double GetKineticEnergy() const
G4double fExcitationRatio
static G4Triton * TritonDefinition()
Definition: G4Triton.cc:90
G4Scintillation(const G4String &processName="Scintillation", G4ProcessType type=fElectromagnetic)
G4double GetScintillationYieldByParticleType(const G4Track &aTrack, const G4Step &aStep)
G4PhysicsTable * fSlowIntegralTable
void SetPolarization(G4double polX, G4double polY, G4double polZ)
G4int GetNumberOfSecondaries() const
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:45
G4EmSaturation * fEmSaturation
static size_t GetNumberOfMaterials()
Definition: G4Material.cc:600
static constexpr double MeV
Definition: G4SIunits.hh:214
static G4Deuteron * DeuteronDefinition()
Definition: G4Deuteron.cc:89
static G4MaterialTable * GetMaterialTable()
Definition: G4Material.cc:593
void SetKineticEnergy(G4double aEnergy)
G4bool fTrackSecondariesFirst
G4StepPoint * GetPreStepPoint() const
size_t GetIndex() const
Definition: G4Material.hh:262
#define G4endl
Definition: G4ios.hh:61
G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep) override
G4double GetEnergy(G4double aValue)
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
G4bool fScintillationByParticleType
const G4String & GetParticleName() const
G4int GetTrackID() const
const G4String & GetParticleType() const
G4double GetPDGCharge() const
double z() const
G4bool fScintillationTrackInfo
const G4ParticleDefinition * GetParticleDefinition() const
void insertAt(size_t, G4PhysicsVector *)
G4double Value(G4double theEnergy, size_t &lastidx) const
G4double condition(const G4ErrorSymMatrix &m)
const G4TouchableHandle & GetTouchableHandle() const
void AddSecondary(G4Track *aSecondary)
G4double GetVelocity() const
G4double GetVertexKineticEnergy() const
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:289
void DumpValues(G4double unitE=1.0, G4double unitV=1.0) const
G4PhysicsTable * fFastIntegralTable
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
G4StepStatus GetStepStatus() const
G4ParticleDefinition * GetDefinition() const
static G4Proton * ProtonDefinition()
Definition: G4Proton.cc:88
virtual void Initialize(const G4Track &)
G4ProcessType
G4double GetStepLength() const
const G4ThreeVector & GetPosition() const
#define G4UniformRand()
Definition: Randomize.hh:53
static constexpr double twopi
Definition: G4SIunits.hh:76
static G4OpticalPhoton * OpticalPhoton()
void SetTouchableHandle(const G4TouchableHandle &apValue)
Float_t d
const G4String & GetProcessName() const
Definition: G4VProcess.hh:411
G4double GetTotalEnergyDeposit() const
G4double GetGlobalTime() const
Definition: G4Step.hh:76
G4StepPoint * GetPostStepPoint() const
G4Material * GetMaterial() const
G4double GetKineticEnergy() const
ThreeVector shoot(const G4int Ap, const G4int Af)
G4double single_exp(G4double t, G4double tau2)
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:89
Hep3Vector cross(const Hep3Vector &) const
Hep3Vector unit() const
std::vector< G4Material * > G4MaterialTable
G4VParticleChange * AtRestDoIt(const G4Track &aTrack, const G4Step &aStep) override
static G4Alpha * AlphaDefinition()
Definition: G4Alpha.cc:84
G4double VisibleEnergyDepositionAtAStep(const G4Step *) const
void SetScintillationByParticleType(const G4bool)
static G4Gamma * GammaDefinition()
Definition: G4Gamma.cc:81
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.hh:65
G4double GetMaxEnergy() const
int G4int
Definition: G4Types.hh:78
size_t entries() const
G4int verboseLevel
Definition: G4VProcess.hh:371
void SetParentID(const G4int aValue)
void BuildPhysicsTable(const G4ParticleDefinition &aParticleType) override
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:435
G4ForceCondition
static G4Neutron * NeutronDefinition()
Definition: G4Neutron.cc:99
G4double GetMeanFreePath(const G4Track &aTrack, G4double, G4ForceCondition *) override
G4GLOB_DLL std::ostream G4cout
double x() const
G4double bi_exp(G4double t, G4double tau1, G4double tau2)
G4MaterialPropertyVector * GetProperty(const char *key, G4bool warning=false)
G4bool IsApplicable(const G4ParticleDefinition &aParticleType) override
void clearAndDestroy()
G4long G4Poisson(G4double mean)
Definition: G4Poisson.hh:51
G4ThreeVector GetDeltaPosition() const
G4double GetMeanLifeTime(const G4Track &aTrack, G4ForceCondition *) override
double y() const
G4TrackStatus GetTrackStatus() const
void DumpPhysicsTable() const
#define DBL_MAX
Definition: templates.hh:83
void ProposeTrackStatus(G4TrackStatus status)
G4double sample_time(G4double tau1, G4double tau2)
#define ns
Definition: xmlparse.cc:614
const G4DynamicParticle * GetDynamicParticle() const
G4MaterialPropertiesTable * GetMaterialPropertiesTable() const
Definition: G4Material.hh:252
size_t GetVectorLength() const
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
G4int GetParentID() const
void SetNumberOfSecondaries(G4int totSecondaries)