Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4Scintillation Class Reference

#include <G4Scintillation.hh>

+ Inheritance diagram for G4Scintillation:

Public Member Functions

 G4Scintillation (const G4String &processName="Scintillation", G4ProcessType type=fElectromagnetic)
 
 ~G4Scintillation ()
 
G4bool IsApplicable (const G4ParticleDefinition &aParticleType)
 
G4double GetMeanFreePath (const G4Track &aTrack, G4double, G4ForceCondition *)
 
G4double GetMeanLifeTime (const G4Track &aTrack, G4ForceCondition *)
 
G4VParticleChangePostStepDoIt (const G4Track &aTrack, const G4Step &aStep)
 
G4VParticleChangeAtRestDoIt (const G4Track &aTrack, const G4Step &aStep)
 
void SetTrackSecondariesFirst (const G4bool state)
 
void SetFiniteRiseTime (const G4bool state)
 
G4bool GetTrackSecondariesFirst () const
 
G4bool GetFiniteRiseTime () const
 
void SetScintillationYieldFactor (const G4double yieldfactor)
 
G4double GetScintillationYieldFactor () const
 
void SetScintillationExcitationRatio (const G4double excitationratio)
 
G4double GetScintillationExcitationRatio () const
 
G4PhysicsTableGetFastIntegralTable () const
 
G4PhysicsTableGetSlowIntegralTable () const
 
void AddSaturation (G4EmSaturation *sat)
 
void RemoveSaturation ()
 
G4EmSaturationGetSaturation () const
 
void SetScintillationByParticleType (const G4bool)
 
G4bool GetScintillationByParticleType () const
 
void DumpPhysicsTable () const
 
- Public Member Functions inherited from G4VRestDiscreteProcess
 G4VRestDiscreteProcess (const G4String &, G4ProcessType aType=fNotDefined)
 
 G4VRestDiscreteProcess (G4VRestDiscreteProcess &)
 
virtual ~G4VRestDiscreteProcess ()
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &)
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &)
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &, const G4Step &)
 
- Public Member Functions inherited from G4VProcess
 G4VProcess (const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
 
 G4VProcess (const G4VProcess &right)
 
virtual ~G4VProcess ()
 
G4int operator== (const G4VProcess &right) const
 
G4int operator!= (const G4VProcess &right) const
 
virtual G4VParticleChangePostStepDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)=0
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &track, G4ForceCondition *condition)=0
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)=0
 
G4double GetCurrentInteractionLength () const
 
void SetPILfactor (G4double value)
 
G4double GetPILfactor () const
 
G4double AlongStepGPIL (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
 
G4double AtRestGPIL (const G4Track &track, G4ForceCondition *condition)
 
G4double PostStepGPIL (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4bool IsApplicable (const G4ParticleDefinition &)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void PreparePhysicsTable (const G4ParticleDefinition &)
 
virtual G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
virtual G4bool RetrievePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
const G4StringGetPhysicsTableFileName (const G4ParticleDefinition *, const G4String &directory, const G4String &tableName, G4bool ascii=false)
 
const G4StringGetProcessName () const
 
G4ProcessType GetProcessType () const
 
void SetProcessType (G4ProcessType)
 
G4int GetProcessSubType () const
 
void SetProcessSubType (G4int)
 
virtual void StartTracking (G4Track *)
 
virtual void EndTracking ()
 
virtual void SetProcessManager (const G4ProcessManager *)
 
virtual const G4ProcessManagerGetProcessManager ()
 
virtual void ResetNumberOfInteractionLengthLeft ()
 
G4double GetNumberOfInteractionLengthLeft () const
 
G4double GetTotalNumberOfInteractionLengthTraversed () const
 
G4bool isAtRestDoItIsEnabled () const
 
G4bool isAlongStepDoItIsEnabled () const
 
G4bool isPostStepDoItIsEnabled () const
 
virtual void DumpInfo () const
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 

Protected Member Functions

void BuildThePhysicsTable ()
 
virtual G4double GetMeanFreePath (const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)=0
 
virtual G4double GetMeanLifeTime (const G4Track &aTrack, G4ForceCondition *condition)=0
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double previousStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

Protected Attributes

G4PhysicsTabletheSlowIntegralTable
 
G4PhysicsTabletheFastIntegralTable
 
G4bool fTrackSecondariesFirst
 
G4bool fFiniteRiseTime
 
G4double YieldFactor
 
G4double ExcitationRatio
 
G4bool scintillationByParticleType
 
- Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager
 
G4VParticleChangepParticleChange
 
G4ParticleChange aParticleChange
 
G4double theNumberOfInteractionLengthLeft
 
G4double currentInteractionLength
 
G4double theInitialNumberOfInteractionLength
 
G4String theProcessName
 
G4String thePhysicsTableFileName
 
G4ProcessType theProcessType
 
G4int theProcessSubType
 
G4double thePILfactor
 
G4bool enableAtRestDoIt
 
G4bool enableAlongStepDoIt
 
G4bool enablePostStepDoIt
 
G4int verboseLevel
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 

Detailed Description

Definition at line 88 of file G4Scintillation.hh.

Constructor & Destructor Documentation

◆ G4Scintillation()

G4Scintillation::G4Scintillation ( const G4String processName = "Scintillation",
G4ProcessType  type = fElectromagnetic 
)

Definition at line 96 of file G4Scintillation.cc.

98 : G4VRestDiscreteProcess(processName, type)
99{
101
103 fFiniteRiseTime = false;
104
105 YieldFactor = 1.0;
106 ExcitationRatio = 1.0;
107
109
112
113 if (verboseLevel>0) {
114 G4cout << GetProcessName() << " is created " << G4endl;
115 }
116
118
119 emSaturation = NULL;
120}
@ fScintillation
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
G4bool scintillationByParticleType
G4double ExcitationRatio
G4bool fTrackSecondariesFirst
G4PhysicsTable * theSlowIntegralTable
G4PhysicsTable * theFastIntegralTable
G4int verboseLevel
Definition: G4VProcess.hh:368
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:403
const G4String & GetProcessName() const
Definition: G4VProcess.hh:379

◆ ~G4Scintillation()

G4Scintillation::~G4Scintillation ( )

Definition at line 126 of file G4Scintillation.cc.

127{
128 if (theFastIntegralTable != NULL) {
131 }
132 if (theSlowIntegralTable != NULL) {
135 }
136}
void clearAndDestroy()

Member Function Documentation

◆ AddSaturation()

void G4Scintillation::AddSaturation ( G4EmSaturation sat)
inline

Definition at line 183 of file G4Scintillation.hh.

183{ emSaturation = sat; }

Referenced by G4OpticalPhysics::AddScintillationSaturation(), and G4OpticalPhysics::ConstructProcess().

◆ AtRestDoIt()

G4VParticleChange * G4Scintillation::AtRestDoIt ( const G4Track aTrack,
const G4Step aStep 
)
virtual

Reimplemented from G4VRestDiscreteProcess.

Definition at line 146 of file G4Scintillation.cc.

151{
152 return G4Scintillation::PostStepDoIt(aTrack, aStep);
153}
G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep)

◆ BuildThePhysicsTable()

void G4Scintillation::BuildThePhysicsTable ( )
protected

Definition at line 537 of file G4Scintillation.cc.

538{
540
541 const G4MaterialTable* theMaterialTable =
543 G4int numOfMaterials = G4Material::GetNumberOfMaterials();
544
545 // create new physics table
546
549
550 // loop for materials
551
552 for (G4int i=0 ; i < numOfMaterials; i++)
553 {
554 G4PhysicsOrderedFreeVector* aPhysicsOrderedFreeVector =
556 G4PhysicsOrderedFreeVector* bPhysicsOrderedFreeVector =
558
559 // Retrieve vector of scintillation wavelength intensity for
560 // the material from the material's optical properties table.
561
562 G4Material* aMaterial = (*theMaterialTable)[i];
563
564 G4MaterialPropertiesTable* aMaterialPropertiesTable =
565 aMaterial->GetMaterialPropertiesTable();
566
567 if (aMaterialPropertiesTable) {
568
569 G4MaterialPropertyVector* theFastLightVector =
570 aMaterialPropertiesTable->GetProperty("FASTCOMPONENT");
571
572 if (theFastLightVector) {
573
574 // Retrieve the first intensity point in vector
575 // of (photon energy, intensity) pairs
576
577 G4double currentIN = (*theFastLightVector)[0];
578
579 if (currentIN >= 0.0) {
580
581 // Create first (photon energy, Scintillation
582 // Integral pair
583
584 G4double currentPM = theFastLightVector->Energy(0);
585
586 G4double currentCII = 0.0;
587
588 aPhysicsOrderedFreeVector->
589 InsertValues(currentPM , currentCII);
590
591 // Set previous values to current ones prior to loop
592
593 G4double prevPM = currentPM;
594 G4double prevCII = currentCII;
595 G4double prevIN = currentIN;
596
597 // loop over all (photon energy, intensity)
598 // pairs stored for this material
599
600 for (size_t ii = 1;
601 ii < theFastLightVector->GetVectorLength();
602 ++ii)
603 {
604 currentPM = theFastLightVector->Energy(ii);
605 currentIN = (*theFastLightVector)[ii];
606
607 currentCII = 0.5 * (prevIN + currentIN);
608
609 currentCII = prevCII +
610 (currentPM - prevPM) * currentCII;
611
612 aPhysicsOrderedFreeVector->
613 InsertValues(currentPM, currentCII);
614
615 prevPM = currentPM;
616 prevCII = currentCII;
617 prevIN = currentIN;
618 }
619
620 }
621 }
622
623 G4MaterialPropertyVector* theSlowLightVector =
624 aMaterialPropertiesTable->GetProperty("SLOWCOMPONENT");
625
626 if (theSlowLightVector) {
627
628 // Retrieve the first intensity point in vector
629 // of (photon energy, intensity) pairs
630
631 G4double currentIN = (*theSlowLightVector)[0];
632
633 if (currentIN >= 0.0) {
634
635 // Create first (photon energy, Scintillation
636 // Integral pair
637
638 G4double currentPM = theSlowLightVector->Energy(0);
639
640 G4double currentCII = 0.0;
641
642 bPhysicsOrderedFreeVector->
643 InsertValues(currentPM , currentCII);
644
645 // Set previous values to current ones prior to loop
646
647 G4double prevPM = currentPM;
648 G4double prevCII = currentCII;
649 G4double prevIN = currentIN;
650
651 // loop over all (photon energy, intensity)
652 // pairs stored for this material
653
654 for (size_t ii = 1;
655 ii < theSlowLightVector->GetVectorLength();
656 ++ii)
657 {
658 currentPM = theSlowLightVector->Energy(ii);
659 currentIN = (*theSlowLightVector)[ii];
660
661 currentCII = 0.5 * (prevIN + currentIN);
662
663 currentCII = prevCII +
664 (currentPM - prevPM) * currentCII;
665
666 bPhysicsOrderedFreeVector->
667 InsertValues(currentPM, currentCII);
668
669 prevPM = currentPM;
670 prevCII = currentCII;
671 prevIN = currentIN;
672 }
673
674 }
675 }
676 }
677
678 // The scintillation integral(s) for a given material
679 // will be inserted in the table(s) according to the
680 // position of the material in the material table.
681
682 theFastIntegralTable->insertAt(i,aPhysicsOrderedFreeVector);
683 theSlowIntegralTable->insertAt(i,bPhysicsOrderedFreeVector);
684
685 }
686}
std::vector< G4Material * > G4MaterialTable
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
G4MaterialPropertyVector * GetProperty(const char *key)
static const G4MaterialTable * GetMaterialTable()
Definition: G4Material.cc:562
G4MaterialPropertiesTable * GetMaterialPropertiesTable() const
Definition: G4Material.hh:251
static size_t GetNumberOfMaterials()
Definition: G4Material.cc:569
void insertAt(size_t, G4PhysicsVector *)
size_t GetVectorLength() const
G4double Energy(size_t index) const

Referenced by G4Scintillation().

◆ DumpPhysicsTable()

void G4Scintillation::DumpPhysicsTable ( ) const
inline

Definition at line 312 of file G4Scintillation.hh.

313{
315 G4int PhysicsTableSize = theFastIntegralTable->entries();
317
318 for (G4int i = 0 ; i < PhysicsTableSize ; i++ )
319 {
321 v->DumpValues();
322 }
323 }
324
326 G4int PhysicsTableSize = theSlowIntegralTable->entries();
328
329 for (G4int i = 0 ; i < PhysicsTableSize ; i++ )
330 {
332 v->DumpValues();
333 }
334 }
335}
size_t entries() const

◆ GetFastIntegralTable()

G4PhysicsTable * G4Scintillation::GetFastIntegralTable ( ) const
inline

Definition at line 306 of file G4Scintillation.hh.

307{
309}

◆ GetFiniteRiseTime()

G4bool G4Scintillation::GetFiniteRiseTime ( ) const
inline

Definition at line 270 of file G4Scintillation.hh.

271{
272 return fFiniteRiseTime;
273}

◆ GetMeanFreePath()

G4double G4Scintillation::GetMeanFreePath ( const G4Track aTrack,
G4double  ,
G4ForceCondition condition 
)
virtual

Implements G4VRestDiscreteProcess.

Definition at line 705 of file G4Scintillation.cc.

708{
710
711 return DBL_MAX;
712
713}
G4double condition(const G4ErrorSymMatrix &m)
@ StronglyForced
#define DBL_MAX
Definition: templates.hh:83

◆ GetMeanLifeTime()

G4double G4Scintillation::GetMeanLifeTime ( const G4Track aTrack,
G4ForceCondition condition 
)
virtual

Implements G4VRestDiscreteProcess.

Definition at line 719 of file G4Scintillation.cc.

721{
722 *condition = Forced;
723
724 return DBL_MAX;
725
726}
@ Forced

◆ GetSaturation()

G4EmSaturation * G4Scintillation::GetSaturation ( ) const
inline

Definition at line 189 of file G4Scintillation.hh.

189{ return emSaturation; }

◆ GetScintillationByParticleType()

G4bool G4Scintillation::GetScintillationByParticleType ( ) const
inline

Definition at line 196 of file G4Scintillation.hh.

◆ GetScintillationExcitationRatio()

G4double G4Scintillation::GetScintillationExcitationRatio ( ) const
inline

Definition at line 294 of file G4Scintillation.hh.

295{
296 return ExcitationRatio;
297}

◆ GetScintillationYieldFactor()

G4double G4Scintillation::GetScintillationYieldFactor ( ) const
inline

Definition at line 282 of file G4Scintillation.hh.

283{
284 return YieldFactor;
285}

◆ GetSlowIntegralTable()

G4PhysicsTable * G4Scintillation::GetSlowIntegralTable ( ) const
inline

Definition at line 300 of file G4Scintillation.hh.

301{
303}

◆ GetTrackSecondariesFirst()

G4bool G4Scintillation::GetTrackSecondariesFirst ( ) const
inline

Definition at line 264 of file G4Scintillation.hh.

265{
267}

◆ IsApplicable()

G4bool G4Scintillation::IsApplicable ( const G4ParticleDefinition aParticleType)
inlinevirtual

Reimplemented from G4VProcess.

Definition at line 243 of file G4Scintillation.hh.

244{
245 if (aParticleType.GetParticleName() == "opticalphoton") return false;
246 if (aParticleType.IsShortLived()) return false;
247
248 return true;
249}
const G4String & GetParticleName() const

Referenced by G4OpticalPhysics::ConstructProcess().

◆ PostStepDoIt()

G4VParticleChange * G4Scintillation::PostStepDoIt ( const G4Track aTrack,
const G4Step aStep 
)
virtual

Reimplemented from G4VRestDiscreteProcess.

Definition at line 159 of file G4Scintillation.cc.

166{
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("FASTCOMPONENT");
188 G4MaterialPropertyVector* Slow_Intensity =
189 aMaterialPropertiesTable->GetProperty("SLOWCOMPONENT");
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
200 // The scintillation response is a function of the energy
201 // deposited by particle types.
202
203 // Get the definition of the current particle
204 G4ParticleDefinition *pDef = aParticle->GetDefinition();
205 G4MaterialPropertyVector *Scint_Yield_Vector = NULL;
206
207 // Obtain the G4MaterialPropertyVectory containing the
208 // scintillation light yield as a function of the deposited
209 // energy for the current particle type
210
211 // Protons
212 if(pDef==G4Proton::ProtonDefinition())
213 Scint_Yield_Vector = aMaterialPropertiesTable->
214 GetProperty("PROTONSCINTILLATIONYIELD");
215
216 // Deuterons
217 else if(pDef==G4Deuteron::DeuteronDefinition())
218 Scint_Yield_Vector = aMaterialPropertiesTable->
219 GetProperty("DEUTERONSCINTILLATIONYIELD");
220
221 // Tritons
222 else if(pDef==G4Triton::TritonDefinition())
223 Scint_Yield_Vector = aMaterialPropertiesTable->
224 GetProperty("TRITONSCINTILLATIONYIELD");
225
226 // Alphas
227 else if(pDef==G4Alpha::AlphaDefinition())
228 Scint_Yield_Vector = aMaterialPropertiesTable->
229 GetProperty("ALPHASCINTILLATIONYIELD");
230
231 // Ions (particles derived from G4VIon and G4Ions)
232 // and recoil ions below tracking cut from neutrons after hElastic
233 else if(pDef->GetParticleType()== "nucleus" ||
235 Scint_Yield_Vector = aMaterialPropertiesTable->
236 GetProperty("IONSCINTILLATIONYIELD");
237
238 // Electrons (must also account for shell-binding energy
239 // attributed to gamma from standard PhotoElectricEffect)
240 else if(pDef==G4Electron::ElectronDefinition() ||
242 Scint_Yield_Vector = aMaterialPropertiesTable->
243 GetProperty("ELECTRONSCINTILLATIONYIELD");
244
245 // Default for particles not enumerated/listed above
246 else
247 Scint_Yield_Vector = aMaterialPropertiesTable->
248 GetProperty("ELECTRONSCINTILLATIONYIELD");
249
250 // If the user has not specified yields for (p,d,t,a,carbon)
251 // then these unspecified particles will default to the
252 // electron's scintillation yield
253 if(!Scint_Yield_Vector){
254 Scint_Yield_Vector = aMaterialPropertiesTable->
255 GetProperty("ELECTRONSCINTILLATIONYIELD");
256 }
257
258 // Throw an exception if no scintillation yield is found
259 if (!Scint_Yield_Vector) {
261 ed << "\nG4Scintillation::PostStepDoIt(): "
262 << "Request for scintillation yield for energy deposit and particle type without correct entry in MaterialPropertiesTable\n"
263 << "ScintillationByParticleType requires at minimum that ELECTRONSCINTILLATIONYIELD is set by the user\n"
264 << G4endl;
265 G4String comments = "Missing MaterialPropertiesTable entry - No correct entry in MaterialPropertiesTable";
266 G4Exception("G4Scintillation::PostStepDoIt","Scint01",
267 FatalException,ed,comments);
268 return G4VRestDiscreteProcess::PostStepDoIt(aTrack, aStep);
269 }
270
271 if (verboseLevel>1) {
272 G4cout << "\n"
273 << "Particle = " << pDef->GetParticleName() << "\n"
274 << "Energy Dep. = " << TotalEnergyDeposit/MeV << "\n"
275 << "Yield = "
276 << Scint_Yield_Vector->Value(TotalEnergyDeposit)
277 << "\n" << G4endl;
278 }
279
280 // Obtain the scintillation yield using the total energy
281 // deposited by the particle in this step.
282
283 // Units: [# scintillation photons]
284 ScintillationYield = Scint_Yield_Vector->
285 Value(TotalEnergyDeposit);
286 } else {
287 // The default linear scintillation process
288 ScintillationYield = aMaterialPropertiesTable->
289 GetConstProperty("SCINTILLATIONYIELD");
290
291 // Units: [# scintillation photons / MeV]
292 ScintillationYield *= YieldFactor;
293 }
294
295 G4double ResolutionScale = aMaterialPropertiesTable->
296 GetConstProperty("RESOLUTIONSCALE");
297
298 // Birks law saturation:
299
300 //G4double constBirks = 0.0;
301
302 //constBirks = aMaterial->GetIonisation()->GetBirksConstant();
303
304 G4double MeanNumberOfPhotons;
305
306 // Birk's correction via emSaturation and specifying scintillation by
307 // by particle type are physically mutually exclusive
308
310 MeanNumberOfPhotons = ScintillationYield;
311 else if (emSaturation)
312 MeanNumberOfPhotons = ScintillationYield*
313 (emSaturation->VisibleEnergyDeposition(&aStep));
314 else
315 MeanNumberOfPhotons = ScintillationYield*TotalEnergyDeposit;
316
317 G4int NumPhotons;
318
319 if (MeanNumberOfPhotons > 10.)
320 {
321 G4double sigma = ResolutionScale * std::sqrt(MeanNumberOfPhotons);
322 NumPhotons = G4int(G4RandGauss::shoot(MeanNumberOfPhotons,sigma)+0.5);
323 }
324 else
325 {
326 NumPhotons = G4int(G4Poisson(MeanNumberOfPhotons));
327 }
328
329 if (NumPhotons <= 0)
330 {
331 // return unchanged particle and no secondaries
332
334
335 return G4VRestDiscreteProcess::PostStepDoIt(aTrack, aStep);
336 }
337
338 ////////////////////////////////////////////////////////////////
339
341
343 if (aTrack.GetTrackStatus() == fAlive )
345 }
346
347 ////////////////////////////////////////////////////////////////
348
349 G4int materialIndex = aMaterial->GetIndex();
350
351 // Retrieve the Scintillation Integral for this material
352 // new G4PhysicsOrderedFreeVector allocated to hold CII's
353
354 G4int Num = NumPhotons;
355
356 for (G4int scnt = 1; scnt <= nscnt; scnt++) {
357
358 G4double ScintillationTime = 0.*ns;
359 G4double ScintillationRiseTime = 0.*ns;
360 G4PhysicsOrderedFreeVector* ScintillationIntegral = NULL;
361
362 if (scnt == 1) {
363 if (nscnt == 1) {
364 if(Fast_Intensity){
365 ScintillationTime = aMaterialPropertiesTable->
366 GetConstProperty("FASTTIMECONSTANT");
367 if (fFiniteRiseTime) {
368 ScintillationRiseTime = aMaterialPropertiesTable->
369 GetConstProperty("FASTSCINTILLATIONRISETIME");
370 }
371 ScintillationIntegral =
373 }
374 if(Slow_Intensity){
375 ScintillationTime = aMaterialPropertiesTable->
376 GetConstProperty("SLOWTIMECONSTANT");
377 if (fFiniteRiseTime) {
378 ScintillationRiseTime = aMaterialPropertiesTable->
379 GetConstProperty("SLOWSCINTILLATIONRISETIME");
380 }
381 ScintillationIntegral =
383 }
384 }
385 else {
386 G4double YieldRatio = aMaterialPropertiesTable->
387 GetConstProperty("YIELDRATIO");
388 if ( ExcitationRatio == 1.0 ) {
389 Num = G4int (std::min(YieldRatio,1.0) * NumPhotons);
390 }
391 else {
392 Num = G4int (std::min(ExcitationRatio,1.0) * NumPhotons);
393 }
394 ScintillationTime = aMaterialPropertiesTable->
395 GetConstProperty("FASTTIMECONSTANT");
396 if (fFiniteRiseTime) {
397 ScintillationRiseTime = aMaterialPropertiesTable->
398 GetConstProperty("FASTSCINTILLATIONRISETIME");
399 }
400 ScintillationIntegral =
402 }
403 }
404 else {
405 Num = NumPhotons - Num;
406 ScintillationTime = aMaterialPropertiesTable->
407 GetConstProperty("SLOWTIMECONSTANT");
408 if (fFiniteRiseTime) {
409 ScintillationRiseTime = aMaterialPropertiesTable->
410 GetConstProperty("SLOWSCINTILLATIONRISETIME");
411 }
412 ScintillationIntegral =
414 }
415
416 if (!ScintillationIntegral) continue;
417
418 // Max Scintillation Integral
419
420 G4double CIImax = ScintillationIntegral->GetMaxValue();
421
422 for (G4int i = 0; i < Num; i++) {
423
424 // Determine photon energy
425
426 G4double CIIvalue = G4UniformRand()*CIImax;
427 G4double sampledEnergy =
428 ScintillationIntegral->GetEnergy(CIIvalue);
429
430 if (verboseLevel>1) {
431 G4cout << "sampledEnergy = " << sampledEnergy << G4endl;
432 G4cout << "CIIvalue = " << CIIvalue << G4endl;
433 }
434
435 // Generate random photon direction
436
437 G4double cost = 1. - 2.*G4UniformRand();
438 G4double sint = std::sqrt((1.-cost)*(1.+cost));
439
440 G4double phi = twopi*G4UniformRand();
441 G4double sinp = std::sin(phi);
442 G4double cosp = std::cos(phi);
443
444 G4double px = sint*cosp;
445 G4double py = sint*sinp;
446 G4double pz = cost;
447
448 // Create photon momentum direction vector
449
450 G4ParticleMomentum photonMomentum(px, py, pz);
451
452 // Determine polarization of new photon
453
454 G4double sx = cost*cosp;
455 G4double sy = cost*sinp;
456 G4double sz = -sint;
457
458 G4ThreeVector photonPolarization(sx, sy, sz);
459
460 G4ThreeVector perp = photonMomentum.cross(photonPolarization);
461
462 phi = twopi*G4UniformRand();
463 sinp = std::sin(phi);
464 cosp = std::cos(phi);
465
466 photonPolarization = cosp * photonPolarization + sinp * perp;
467
468 photonPolarization = photonPolarization.unit();
469
470 // Generate a new photon:
471
472 G4DynamicParticle* aScintillationPhoton =
474 photonMomentum);
475 aScintillationPhoton->SetPolarization
476 (photonPolarization.x(),
477 photonPolarization.y(),
478 photonPolarization.z());
479
480 aScintillationPhoton->SetKineticEnergy(sampledEnergy);
481
482 // Generate new G4Track object:
483
484 G4double rand;
485
486 if (aParticle->GetDefinition()->GetPDGCharge() != 0) {
487 rand = G4UniformRand();
488 } else {
489 rand = 1.0;
490 }
491
492 G4double delta = rand * aStep.GetStepLength();
493 G4double deltaTime = delta /
494 ((pPreStepPoint->GetVelocity()+
495 pPostStepPoint->GetVelocity())/2.);
496
497 // emission time distribution
498 if (ScintillationRiseTime==0.0) {
499 deltaTime = deltaTime -
500 ScintillationTime * std::log( G4UniformRand() );
501 } else {
502 deltaTime = deltaTime +
503 sample_time(ScintillationRiseTime, ScintillationTime);
504 }
505
506 G4double aSecondaryTime = t0 + deltaTime;
507
508 G4ThreeVector aSecondaryPosition =
509 x0 + rand * aStep.GetDeltaPosition();
510
511 G4Track* aSecondaryTrack =
512 new G4Track(aScintillationPhoton,aSecondaryTime,aSecondaryPosition);
513
514 aSecondaryTrack->SetTouchableHandle(
516 // aSecondaryTrack->SetTouchableHandle((G4VTouchable*)0);
517
518 aSecondaryTrack->SetParentID(aTrack.GetTrackID());
519
520 aParticleChange.AddSecondary(aSecondaryTrack);
521
522 }
523 }
524
525 if (verboseLevel>0) {
526 G4cout << "\n Exiting from G4Scintillation::DoIt -- NumberOfSecondaries = "
528 }
529
530 return G4VRestDiscreteProcess::PostStepDoIt(aTrack, aStep);
531}
@ FatalException
G4long G4Poisson(G4double mean)
Definition: G4Poisson.hh:50
@ fSuspend
@ fAlive
#define G4UniformRand()
Definition: Randomize.hh:53
Hep3Vector unit() const
Hep3Vector cross(const Hep3Vector &) const
static G4Alpha * AlphaDefinition()
Definition: G4Alpha.cc:84
static G4Deuteron * DeuteronDefinition()
Definition: G4Deuteron.cc:89
void SetPolarization(G4double polX, G4double polY, G4double polZ)
G4ParticleDefinition * GetDefinition() const
void SetKineticEnergy(G4double aEnergy)
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:89
G4double VisibleEnergyDeposition(const G4ParticleDefinition *, const G4MaterialCutsCouple *, G4double length, G4double edepTotal, G4double edepNIEL=0.0)
static G4Gamma * GammaDefinition()
Definition: G4Gamma.cc:81
size_t GetIndex() const
Definition: G4Material.hh:261
static G4Neutron * NeutronDefinition()
Definition: G4Neutron.cc:99
static G4OpticalPhoton * OpticalPhoton()
void AddSecondary(G4Track *aSecondary)
virtual void Initialize(const G4Track &)
const G4String & GetParticleType() const
G4double GetPDGCharge() const
G4double GetEnergy(G4double aValue)
G4double Value(G4double theEnergy)
static G4Proton * ProtonDefinition()
Definition: G4Proton.cc:88
G4double GetVelocity() const
G4double GetGlobalTime() const
const G4ThreeVector & GetPosition() const
const G4TouchableHandle & GetTouchableHandle() const
G4ThreeVector GetDeltaPosition() const
G4StepPoint * GetPreStepPoint() const
G4double GetStepLength() const
G4double GetTotalEnergyDeposit() const
G4StepPoint * GetPostStepPoint() const
G4TrackStatus GetTrackStatus() const
G4int GetTrackID() const
void SetTouchableHandle(const G4TouchableHandle &apValue)
G4Material * GetMaterial() const
const G4DynamicParticle * GetDynamicParticle() const
void SetParentID(const G4int aValue)
static G4Triton * TritonDefinition()
Definition: G4Triton.cc:90
void ProposeTrackStatus(G4TrackStatus status)
G4int GetNumberOfSecondaries() const
void SetNumberOfSecondaries(G4int totSecondaries)
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:289
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
#define ns
Definition: xmlparse.cc:597

Referenced by AtRestDoIt().

◆ RemoveSaturation()

void G4Scintillation::RemoveSaturation ( )
inline

Definition at line 186 of file G4Scintillation.hh.

186{ emSaturation = NULL; }

Referenced by SetScintillationByParticleType().

◆ SetFiniteRiseTime()

void G4Scintillation::SetFiniteRiseTime ( const G4bool  state)
inline

Definition at line 258 of file G4Scintillation.hh.

259{
260 fFiniteRiseTime = state;
261}

Referenced by G4OpticalPhysics::ConstructProcess(), and G4OpticalPhysics::SetFiniteRiseTime().

◆ SetScintillationByParticleType()

void G4Scintillation::SetScintillationByParticleType ( const G4bool  scintType)

Definition at line 691 of file G4Scintillation.cc.

692{
693 if (emSaturation) {
694 G4Exception("G4Scintillation::SetScintillationByParticleType", "Scint02",
695 JustWarning, "Redefinition: Birks Saturation is replaced by ScintillationByParticleType!");
697 }
698 scintillationByParticleType = scintType;
699}
@ JustWarning

Referenced by G4OpticalPhysics::ConstructProcess(), and G4OpticalPhysics::SetScintillationByParticleType().

◆ SetScintillationExcitationRatio()

void G4Scintillation::SetScintillationExcitationRatio ( const G4double  excitationratio)
inline

Definition at line 288 of file G4Scintillation.hh.

289{
290 ExcitationRatio = excitationratio;
291}

Referenced by G4OpticalPhysics::ConstructProcess(), and G4OpticalPhysics::SetScintillationExcitationRatio().

◆ SetScintillationYieldFactor()

void G4Scintillation::SetScintillationYieldFactor ( const G4double  yieldfactor)
inline

Definition at line 276 of file G4Scintillation.hh.

277{
278 YieldFactor = yieldfactor;
279}

Referenced by G4OpticalPhysics::ConstructProcess(), and G4OpticalPhysics::SetScintillationYieldFactor().

◆ SetTrackSecondariesFirst()

void G4Scintillation::SetTrackSecondariesFirst ( const G4bool  state)
inline

Definition at line 252 of file G4Scintillation.hh.

253{
255}

Referenced by G4OpticalPhysics::SetTrackSecondariesFirst().

Member Data Documentation

◆ ExcitationRatio

G4double G4Scintillation::ExcitationRatio
protected

◆ fFiniteRiseTime

G4bool G4Scintillation::fFiniteRiseTime
protected

◆ fTrackSecondariesFirst

G4bool G4Scintillation::fTrackSecondariesFirst
protected

◆ scintillationByParticleType

G4bool G4Scintillation::scintillationByParticleType
protected

◆ theFastIntegralTable

G4PhysicsTable* G4Scintillation::theFastIntegralTable
protected

◆ theSlowIntegralTable

G4PhysicsTable* G4Scintillation::theSlowIntegralTable
protected

◆ YieldFactor

G4double G4Scintillation::YieldFactor
protected

The documentation for this class was generated from the following files: