Geant4 10.7.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) override
 
void BuildPhysicsTable (const G4ParticleDefinition &aParticleType) override
 
void PreparePhysicsTable (const G4ParticleDefinition &part) override
 
void Initialise ()
 
G4double GetMeanFreePath (const G4Track &aTrack, G4double, G4ForceCondition *) override
 
G4double GetMeanLifeTime (const G4Track &aTrack, G4ForceCondition *) override
 
G4VParticleChangePostStepDoIt (const G4Track &aTrack, const G4Step &aStep) override
 
G4VParticleChangeAtRestDoIt (const G4Track &aTrack, const G4Step &aStep) override
 
G4double GetScintillationYieldByParticleType (const G4Track &aTrack, const G4Step &aStep)
 
G4double GetScintillationYieldByParticleType (const G4Track &aTrack, const G4Step &aStep, G4double &yield1, G4double &yield2, G4double &yield3)
 
void SetTrackSecondariesFirst (const G4bool state)
 
G4bool GetTrackSecondariesFirst () const
 
void SetFiniteRiseTime (const G4bool state)
 
G4bool GetFiniteRiseTime () const
 
void SetScintillationYieldFactor (const G4double yieldfactor)
 
G4double GetScintillationYieldFactor () const
 
void SetScintillationExcitationRatio (const G4double ratio)
 
G4double GetScintillationExcitationRatio () const
 
G4PhysicsTableGetFastIntegralTable () const
 
G4PhysicsTableGetSlowIntegralTable () const
 
G4PhysicsTableGetIntegralTable1 () const
 
G4PhysicsTableGetIntegralTable2 () const
 
G4PhysicsTableGetIntegralTable3 () const
 
void AddSaturation (G4EmSaturation *sat)
 
void RemoveSaturation ()
 
G4EmSaturationGetSaturation () const
 
void SetScintillationByParticleType (const G4bool)
 
G4bool GetScintillationByParticleType () const
 
void SetEnhancedTimeConstants (G4bool)
 
G4bool GetEnhancedTimeConstants () const
 
void SetScintillationTrackInfo (const G4bool trackType)
 
G4bool GetScintillationTrackInfo () const
 
void SetStackPhotons (const G4bool)
 
G4bool GetStackPhotons () const
 
G4int GetNumPhotons () const
 
void DumpPhysicsTable () const
 
- Public Member Functions inherited from G4VRestDiscreteProcess
 G4VRestDiscreteProcess (const G4String &, G4ProcessType aType=fNotDefined)
 
 G4VRestDiscreteProcess (G4VRestDiscreteProcess &)
 
virtual ~G4VRestDiscreteProcess ()
 
G4VRestDiscreteProcessoperator= (const G4VRestDiscreteProcess &)=delete
 
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 ()
 
G4VProcessoperator= (const G4VProcess &)=delete
 
G4bool operator== (const G4VProcess &right) const
 
G4bool 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
 
virtual void ProcessDescription (std::ostream &outfile) const
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 
virtual void SetMasterProcess (G4VProcess *masterP)
 
const G4VProcessGetMasterProcess () const
 
virtual void BuildWorkerPhysicsTable (const G4ParticleDefinition &part)
 
virtual void PrepareWorkerPhysicsTable (const G4ParticleDefinition &)
 

Protected Attributes

G4PhysicsTablefIntegralTable1
 
G4PhysicsTablefIntegralTable2
 
G4PhysicsTablefIntegralTable3
 
- Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager = nullptr
 
G4VParticleChangepParticleChange = nullptr
 
G4ParticleChange aParticleChange
 
G4double theNumberOfInteractionLengthLeft = -1.0
 
G4double currentInteractionLength = -1.0
 
G4double theInitialNumberOfInteractionLength = -1.0
 
G4String theProcessName
 
G4String thePhysicsTableFileName
 
G4ProcessType theProcessType = fNotDefined
 
G4int theProcessSubType = -1
 
G4double thePILfactor = 1.0
 
G4int verboseLevel = 0
 
G4bool enableAtRestDoIt = true
 
G4bool enableAlongStepDoIt = true
 
G4bool enablePostStepDoIt = true
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 
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 prevStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

Detailed Description

Definition at line 78 of file G4Scintillation.hh.

Constructor & Destructor Documentation

◆ G4Scintillation()

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

Definition at line 80 of file G4Scintillation.cc.

82 : G4VRestDiscreteProcess(processName, type)
83 , fIntegralTable1(nullptr)
84 , fIntegralTable2(nullptr)
85 , fIntegralTable3(nullptr)
86 , fNumPhotons(0)
87 , fEmSaturation(nullptr)
88{
90
91#ifdef G4DEBUG_SCINTILLATION
92 ScintTrackEDep = 0.;
93 ScintTrackYield = 0.;
94#endif
95
96 if(verboseLevel > 0)
97 {
98 G4cout << GetProcessName() << " is created " << G4endl;
99 }
100 Initialise();
101}
@ fScintillation
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4PhysicsTable * fIntegralTable3
G4PhysicsTable * fIntegralTable1
G4PhysicsTable * fIntegralTable2
G4int verboseLevel
Definition: G4VProcess.hh:356
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:406
const G4String & GetProcessName() const
Definition: G4VProcess.hh:382

◆ ~G4Scintillation()

G4Scintillation::~G4Scintillation ( )

Definition at line 104 of file G4Scintillation.cc.

105{
106 if(fIntegralTable1 != nullptr)
107 {
109 delete fIntegralTable1;
110 }
111 if(fIntegralTable2 != nullptr)
112 {
114 delete fIntegralTable2;
115 }
116 if(fIntegralTable3 != nullptr)
117 {
119 delete fIntegralTable3;
120 }
121}
void clearAndDestroy()

Member Function Documentation

◆ AddSaturation()

void G4Scintillation::AddSaturation ( G4EmSaturation sat)
inline

Definition at line 335 of file G4Scintillation.hh.

336{
337 fEmSaturation = sat;
338}

Referenced by G4OpticalPhysics::ConstructProcess().

◆ AtRestDoIt()

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

Reimplemented from G4VRestDiscreteProcess.

Definition at line 319 of file G4Scintillation.cc.

323{
324 return G4Scintillation::PostStepDoIt(aTrack, aStep);
325}
G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep) override

◆ BuildPhysicsTable()

void G4Scintillation::BuildPhysicsTable ( const G4ParticleDefinition aParticleType)
overridevirtual

Reimplemented from G4VProcess.

Definition at line 153 of file G4Scintillation.cc.

154{
155 if(fIntegralTable1 != nullptr)
156 {
158 delete fIntegralTable1;
159 fIntegralTable1 = nullptr;
160 }
161 if(fIntegralTable2 != nullptr)
162 {
164 delete fIntegralTable2;
165 fIntegralTable2 = nullptr;
166 }
167 if(fIntegralTable3 != nullptr)
168 {
170 delete fIntegralTable3;
171 fIntegralTable3 = nullptr;
172 }
173
174 const G4MaterialTable* materialTable = G4Material::GetMaterialTable();
175 size_t numOfMaterials = G4Material::GetNumberOfMaterials();
176
177 // create new physics table
178 if(!fIntegralTable1)
179 fIntegralTable1 = new G4PhysicsTable(numOfMaterials);
180 if(!fIntegralTable2)
181 fIntegralTable2 = new G4PhysicsTable(numOfMaterials);
182 if(!fIntegralTable3)
183 fIntegralTable3 = new G4PhysicsTable(numOfMaterials);
184
185 for(size_t i = 0; i < numOfMaterials; ++i)
186 {
190
191 // Retrieve vector of scintillation wavelength intensity for
192 // the material from the material's optical properties table.
194 ((*materialTable)[i])->GetMaterialPropertiesTable();
195
196 if(MPT)
197 {
198 // integral table 1 is either FASTCOMPONENT or SCINTILLATIONCOMPONENT1
200 if(!MPV)
202 if(MPV)
203 {
204 // Retrieve the first intensity point in vector
205 // of (photon energy, intensity) pairs
206 G4double currentIN = (*MPV)[0];
207 if(currentIN >= 0.0)
208 {
209 // Create first (photon energy, Scintillation Integral pair
210 G4double currentPM = MPV->Energy(0);
211 G4double currentCII = 0.0;
212 vector1->InsertValues(currentPM, currentCII);
213
214 // Set previous values to current ones prior to loop
215 G4double prevPM = currentPM;
216 G4double prevCII = currentCII;
217 G4double prevIN = currentIN;
218
219 // loop over all (photon energy, intensity)
220 // pairs stored for this material
221 for(size_t ii = 1; ii < MPV->GetVectorLength(); ++ii)
222 {
223 currentPM = MPV->Energy(ii);
224 currentIN = (*MPV)[ii];
225 currentCII =
226 prevCII + 0.5 * (currentPM - prevPM) * (prevIN + currentIN);
227
228 vector1->InsertValues(currentPM, currentCII);
229
230 prevPM = currentPM;
231 prevCII = currentCII;
232 prevIN = currentIN;
233 }
234 }
235 }
236
237 // integral table 2 is SCINTILLATIONCOMPONENT2
239 if(MPV)
240 {
241 // Retrieve the first intensity point in vector
242 // of (photon energy, intensity) pairs
243 G4double currentIN = (*MPV)[0];
244 if(currentIN >= 0.0)
245 {
246 // Create first (photon energy, Scintillation Integral pair
247 G4double currentPM = MPV->Energy(0);
248 G4double currentCII = 0.0;
249 vector2->InsertValues(currentPM, currentCII);
250
251 // Set previous values to current ones prior to loop
252 G4double prevPM = currentPM;
253 G4double prevCII = currentCII;
254 G4double prevIN = currentIN;
255
256 // loop over all (photon energy, intensity)
257 // pairs stored for this material
258 for(size_t ii = 1; ii < MPV->GetVectorLength(); ++ii)
259 {
260 currentPM = MPV->Energy(ii);
261 currentIN = (*MPV)[ii];
262 currentCII =
263 prevCII + 0.5 * (currentPM - prevPM) * (prevIN + currentIN);
264
265 vector2->InsertValues(currentPM, currentCII);
266
267 prevPM = currentPM;
268 prevCII = currentCII;
269 prevIN = currentIN;
270 }
271 }
272 }
273 // integral table 3 is either SLOWCOMPONENT or SCINTILLATIONCOMPONENT3
274 MPV = MPT->GetProperty(kSLOWCOMPONENT);
275 if(!MPV)
277 if(MPV)
278 {
279 // Retrieve the first intensity point in vector
280 // of (photon energy, intensity) pairs
281 G4double currentIN = (*MPV)[0];
282 if(currentIN >= 0.0)
283 {
284 // Create first (photon energy, Scintillation Integral pair
285 G4double currentPM = MPV->Energy(0);
286 G4double currentCII = 0.0;
287 vector3->InsertValues(currentPM, currentCII);
288
289 // Set previous values to current ones prior to loop
290 G4double prevPM = currentPM;
291 G4double prevCII = currentCII;
292 G4double prevIN = currentIN;
293
294 // loop over all (photon energy, intensity)
295 // pairs stored for this material
296 for(size_t ii = 1; ii < MPV->GetVectorLength(); ++ii)
297 {
298 currentPM = MPV->Energy(ii);
299 currentIN = (*MPV)[ii];
300 currentCII =
301 prevCII + 0.5 * (currentPM - prevPM) * (prevIN + currentIN);
302
303 vector3->InsertValues(currentPM, currentCII);
304
305 prevPM = currentPM;
306 prevCII = currentCII;
307 prevIN = currentIN;
308 }
309 }
310 }
311 }
312 fIntegralTable1->insertAt(i, vector1);
313 fIntegralTable2->insertAt(i, vector2);
314 fIntegralTable3->insertAt(i, vector3);
315 }
316}
@ kSCINTILLATIONCOMPONENT1
@ kSCINTILLATIONCOMPONENT2
@ kSCINTILLATIONCOMPONENT3
std::vector< G4Material * > G4MaterialTable
double G4double
Definition: G4Types.hh:83
G4MaterialPropertyVector * GetProperty(const char *key, G4bool warning=false)
static size_t GetNumberOfMaterials()
Definition: G4Material.cc:644
static G4MaterialTable * GetMaterialTable()
Definition: G4Material.cc:637
void InsertValues(G4double energy, G4double value)
void insertAt(std::size_t, G4PhysicsVector *)
G4double Energy(std::size_t index) const
std::size_t GetVectorLength() const

◆ DumpPhysicsTable()

void G4Scintillation::DumpPhysicsTable ( ) const

Definition at line 1068 of file G4Scintillation.cc.

1069{
1071 if(fIntegralTable1)
1072 {
1073 for(size_t i = 0; i < fIntegralTable1->entries(); ++i)
1074 {
1076 v->DumpValues();
1077 }
1078 }
1079 if(fIntegralTable2)
1080 {
1081 for(size_t i = 0; i < fIntegralTable2->entries(); ++i)
1082 {
1084 v->DumpValues();
1085 }
1086 }
1087 if(fIntegralTable3)
1088 {
1089 for(size_t i = 0; i < fIntegralTable3->entries(); ++i)
1090 {
1092 v->DumpValues();
1093 }
1094 }
1095}
std::size_t entries() const
void DumpValues(G4double unitE=1.0, G4double unitV=1.0) const

◆ GetEnhancedTimeConstants()

G4bool G4Scintillation::GetEnhancedTimeConstants ( ) const
inline

Definition at line 357 of file G4Scintillation.hh.

358{
359 return fEnhancedTimeConstants;
360}

◆ GetFastIntegralTable()

G4PhysicsTable * G4Scintillation::GetFastIntegralTable ( ) const
inline

Definition at line 315 of file G4Scintillation.hh.

316{
317 return fIntegralTable1;
318}

◆ GetFiniteRiseTime()

G4bool G4Scintillation::GetFiniteRiseTime ( ) const
inline

Definition at line 283 of file G4Scintillation.hh.

284{
285 return fFiniteRiseTime;
286}

◆ GetIntegralTable1()

G4PhysicsTable * G4Scintillation::GetIntegralTable1 ( ) const
inline

Definition at line 320 of file G4Scintillation.hh.

321{
322 return fIntegralTable1;
323}

◆ GetIntegralTable2()

G4PhysicsTable * G4Scintillation::GetIntegralTable2 ( ) const
inline

Definition at line 325 of file G4Scintillation.hh.

326{
327 return fIntegralTable2;
328}

◆ GetIntegralTable3()

G4PhysicsTable * G4Scintillation::GetIntegralTable3 ( ) const
inline

Definition at line 330 of file G4Scintillation.hh.

331{
332 return fIntegralTable3;
333}

◆ GetMeanFreePath()

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

Implements G4VRestDiscreteProcess.

Definition at line 703 of file G4Scintillation.cc.

705{
707 return DBL_MAX;
708}
G4double condition(const G4ErrorSymMatrix &m)
@ StronglyForced
#define DBL_MAX
Definition: templates.hh:62

◆ GetMeanLifeTime()

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

Implements G4VRestDiscreteProcess.

Definition at line 711 of file G4Scintillation.cc.

713{
714 *condition = Forced;
715 return DBL_MAX;
716}
@ Forced

◆ GetNumPhotons()

G4int G4Scintillation::GetNumPhotons ( ) const
inline

Definition at line 379 of file G4Scintillation.hh.

379{ return fNumPhotons; }

◆ GetSaturation()

G4EmSaturation * G4Scintillation::GetSaturation ( ) const
inline

Definition at line 342 of file G4Scintillation.hh.

343{
344 return fEmSaturation;
345}

◆ GetScintillationByParticleType()

G4bool G4Scintillation::GetScintillationByParticleType ( ) const
inline

Definition at line 347 of file G4Scintillation.hh.

348{
349 return fScintillationByParticleType;
350}

◆ GetScintillationExcitationRatio()

G4double G4Scintillation::GetScintillationExcitationRatio ( ) const
inline

Definition at line 305 of file G4Scintillation.hh.

306{
307 return fExcitationRatio;
308}

◆ GetScintillationTrackInfo()

G4bool G4Scintillation::GetScintillationTrackInfo ( ) const
inline

Definition at line 367 of file G4Scintillation.hh.

368{
369 return fScintillationTrackInfo;
370}

◆ GetScintillationYieldByParticleType() [1/2]

G4double G4Scintillation::GetScintillationYieldByParticleType ( const G4Track aTrack,
const G4Step aStep 
)

Definition at line 741 of file G4Scintillation.cc.

743{
744 // Get the G4MaterialPropertyVector containing the scintillation
745 // yield as a function of the energy deposited and particle type
746
748 G4MaterialPropertyVector* scintVector = nullptr;
751
752 // Protons
753 if(pDef == G4Proton::ProtonDefinition())
754 scintVector = MPT->GetProperty(kPROTONSCINTILLATIONYIELD);
755
756 // Deuterons
757 else if(pDef == G4Deuteron::DeuteronDefinition())
758 scintVector = MPT->GetProperty(kDEUTERONSCINTILLATIONYIELD);
759
760 // Tritons
761 else if(pDef == G4Triton::TritonDefinition())
762 scintVector = MPT->GetProperty(kTRITONSCINTILLATIONYIELD);
763
764 // Alphas
765 else if(pDef == G4Alpha::AlphaDefinition())
766 scintVector = MPT->GetProperty(kALPHASCINTILLATIONYIELD);
767
768 // Ions (particles derived from G4VIon and G4Ions) and recoil ions
769 // below the production cut from neutrons after hElastic
770 else if(pDef->GetParticleType() == "nucleus" ||
772 scintVector = MPT->GetProperty(kIONSCINTILLATIONYIELD);
773
774 // Electrons (must also
775 else if(pDef == G4Electron::ElectronDefinition() ||
776 pDef == G4Gamma::GammaDefinition())
777 scintVector = MPT->GetProperty(kELECTRONSCINTILLATIONYIELD);
778
779 // Default for particles not enumerated/listed above
780 // includes gamma to account for shell-binding energy
781 // attributed to gamma from standard photoelectric effect)
782 else
783 scintVector = MPT->GetProperty(kELECTRONSCINTILLATIONYIELD);
784
785 // Throw an exception if no scintillation yield vector is found
786 if(!scintVector)
787 {
789 ed << "\nG4Scintillation::PostStepDoIt(): "
790 << "Request for scintillation yield for energy deposit and particle\n"
791 << "type without correct entry in MaterialPropertiesTable.\n"
792 << "ScintillationByParticleType requires at minimum that \n"
793 << "ELECTRONSCINTILLATIONYIELD is set by the user\n"
794 << G4endl;
795 G4String comments = "Missing MaterialPropertiesTable entry - No correct "
796 "entry in MaterialPropertiesTable";
797 G4Exception("G4Scintillation::PostStepDoIt", "Scint01", FatalException, ed,
798 comments);
799 }
800
801 ///////////////////////////////////////
802 // Calculate the scintillation light //
803 ///////////////////////////////////////
804 // To account for potential nonlinearity and scintillation photon
805 // density along the track, light (L) is produced according to:
806 // L_currentStep = L(PreStepKE) - L(PreStepKE - EDep)
807
808 G4double ScintillationYield = 0.;
809 G4double StepEnergyDeposit = aStep.GetTotalEnergyDeposit();
810 G4double PreStepKineticEnergy = aStep.GetPreStepPoint()->GetKineticEnergy();
811
812 if(PreStepKineticEnergy <= scintVector->GetMaxEnergy())
813 {
814 G4double Yield1 = scintVector->Value(PreStepKineticEnergy);
815 G4double Yield2 =
816 scintVector->Value(PreStepKineticEnergy - StepEnergyDeposit);
817 ScintillationYield = Yield1 - Yield2;
818 }
819 else
820 {
822 ed << "\nG4Scintillation::GetScintillationYieldByParticleType(): Request\n"
823 << "for scintillation light yield above the available energy range\n"
824 << "specified in G4MaterialPropertiesTable. A linear interpolation\n"
825 << "will be performed to compute the scintillation light yield using\n"
826 << "(L_max / E_max) as the photon yield per unit energy." << G4endl;
827 G4String cmt = "\nScintillation yield may be unphysical!\n";
828 G4Exception("G4Scintillation::GetScintillationYieldByParticleType()",
829 "Scint03", JustWarning, ed, cmt);
830
831 // Units: [# scintillation photons]
832 ScintillationYield = scintVector->GetMaxValue() /
833 scintVector->GetMaxEnergy() * StepEnergyDeposit;
834 }
835
836#ifdef G4DEBUG_SCINTILLATION
837 // Increment track aggregators
838 ScintTrackYield += ScintillationYield;
839 ScintTrackEDep += StepEnergyDeposit;
840
841 G4cout << "\n--- G4Scintillation::GetScintillationYieldByParticleType() ---\n"
842 << "--\n"
843 << "-- Name = "
844 << aTrack.GetParticleDefinition()->GetParticleName() << "\n"
845 << "-- TrackID = " << aTrack.GetTrackID() << "\n"
846 << "-- ParentID = " << aTrack.GetParentID() << "\n"
847 << "-- Current KE = " << aTrack.GetKineticEnergy() / MeV
848 << " MeV\n"
849 << "-- Step EDep = " << aStep.GetTotalEnergyDeposit() / MeV
850 << " MeV\n"
851 << "-- Track EDep = " << ScintTrackEDep / MeV << " MeV\n"
852 << "-- Vertex KE = " << aTrack.GetVertexKineticEnergy() / MeV
853 << " MeV\n"
854 << "-- Step yield = " << ScintillationYield << " photons\n"
855 << "-- Track yield = " << ScintTrackYield << " photons\n"
856 << G4endl;
857
858 // The track has terminated within or has left the scintillator volume
859 if((aTrack.GetTrackStatus() == fStopButAlive) or
861 {
862 // Reset aggregators for the next track
863 ScintTrackEDep = 0.;
864 ScintTrackYield = 0.;
865 }
866#endif
867
868 return ScintillationYield;
869}
@ JustWarning
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
@ kELECTRONSCINTILLATIONYIELD
@ kALPHASCINTILLATIONYIELD
@ kPROTONSCINTILLATIONYIELD
@ kDEUTERONSCINTILLATIONYIELD
@ kIONSCINTILLATIONYIELD
@ kTRITONSCINTILLATIONYIELD
@ fGeomBoundary
Definition: G4StepStatus.hh:43
@ fStopButAlive
static G4Alpha * AlphaDefinition()
Definition: G4Alpha.cc:83
static G4Deuteron * DeuteronDefinition()
Definition: G4Deuteron.cc:88
G4ParticleDefinition * GetDefinition() const
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:88
static G4Gamma * GammaDefinition()
Definition: G4Gamma.cc:80
G4MaterialPropertiesTable * GetMaterialPropertiesTable() const
Definition: G4Material.hh:254
static G4Neutron * NeutronDefinition()
Definition: G4Neutron.cc:98
const G4String & GetParticleType() const
const G4String & GetParticleName() const
G4double Value(G4double theEnergy, std::size_t &lastidx) const
G4double GetMaxEnergy() const
static G4Proton * ProtonDefinition()
Definition: G4Proton.cc:87
G4StepStatus GetStepStatus() const
G4double GetKineticEnergy() const
G4StepPoint * GetPreStepPoint() const
G4double GetTotalEnergyDeposit() const
G4StepPoint * GetPostStepPoint() const
G4TrackStatus GetTrackStatus() const
G4double GetVertexKineticEnergy() const
G4int GetTrackID() const
const G4ParticleDefinition * GetParticleDefinition() const
G4Material * GetMaterial() const
const G4DynamicParticle * GetDynamicParticle() const
G4double GetKineticEnergy() const
G4int GetParentID() const
static G4Triton * TritonDefinition()
Definition: G4Triton.cc:89

Referenced by PostStepDoIt().

◆ GetScintillationYieldByParticleType() [2/2]

G4double G4Scintillation::GetScintillationYieldByParticleType ( const G4Track aTrack,
const G4Step aStep,
G4double yield1,
G4double yield2,
G4double yield3 
)

Definition at line 872 of file G4Scintillation.cc.

875{
876 // new in 10.7, allow multiple time constants with ScintByParticleType
877 // Get the G4MaterialPropertyVector containing the scintillation
878 // yield as a function of the energy deposited and particle type
879
881 G4MaterialPropertyVector* yieldVector = nullptr;
884
885 // Protons
886 if(pDef == G4Proton::ProtonDefinition())
887 {
888 yieldVector = MPT->GetProperty(kPROTONSCINTILLATIONYIELD);
891 : 1.;
894 : 0.;
897 : 0.;
898 }
899
900 // Deuterons
901 else if(pDef == G4Deuteron::DeuteronDefinition())
902 {
903 yieldVector = MPT->GetProperty(kDEUTERONSCINTILLATIONYIELD);
906 : 1.;
909 : 0.;
912 : 0.;
913 }
914
915 // Tritons
916 else if(pDef == G4Triton::TritonDefinition())
917 {
918 yieldVector = MPT->GetProperty(kTRITONSCINTILLATIONYIELD);
921 : 1.;
924 : 0.;
927 : 0.;
928 }
929
930 // Alphas
931 else if(pDef == G4Alpha::AlphaDefinition())
932 {
933 yieldVector = MPT->GetProperty(kALPHASCINTILLATIONYIELD);
936 : 1.;
939 : 0.;
942 : 0.;
943 }
944
945 // Ions (particles derived from G4VIon and G4Ions) and recoil ions
946 // below the production cut from neutrons after hElastic
947 else if(pDef->GetParticleType() == "nucleus" ||
949 {
950 yieldVector = MPT->GetProperty(kIONSCINTILLATIONYIELD);
953 : 1.;
956 : 0.;
959 : 0.;
960 }
961
962 // Electrons (must also account for shell-binding energy
963 // attributed to gamma from standard photoelectric effect)
964 // and, default for particles not enumerated/listed above
965 else
966 {
967 yieldVector = MPT->GetProperty(kELECTRONSCINTILLATIONYIELD);
970 : 1.;
973 : 0.;
976 : 0.;
977 }
978
979 // Throw an exception if no scintillation yield vector is found
980 if(!yieldVector)
981 {
983 ed << "\nG4Scintillation::PostStepDoIt(): "
984 << "Request for scintillation yield for energy deposit and particle\n"
985 << "type without correct entry in MaterialPropertiesTable.\n"
986 << "ScintillationByParticleType requires at minimum that \n"
987 << "ELECTRONSCINTILLATIONYIELD is set by the user\n"
988 << G4endl;
989 G4String comments = "Missing MaterialPropertiesTable entry - No correct "
990 "entry in MaterialPropertiesTable";
991 G4Exception("G4Scintillation::PostStepDoIt", "Scint01", FatalException, ed,
992 comments);
993 }
994
995 ///////////////////////////////////////
996 // Calculate the scintillation light //
997 ///////////////////////////////////////
998 // To account for potential nonlinearity and scintillation photon
999 // density along the track, light (L) is produced according to:
1000 // L_currentStep = L(PreStepKE) - L(PreStepKE - EDep)
1001
1002 G4double ScintillationYield = 0.;
1003 G4double StepEnergyDeposit = aStep.GetTotalEnergyDeposit();
1004 G4double PreStepKineticEnergy = aStep.GetPreStepPoint()->GetKineticEnergy();
1005
1006 if(PreStepKineticEnergy <= yieldVector->GetMaxEnergy())
1007 {
1008 // G4double Yield1 = yieldVector->Value(PreStepKineticEnergy);
1009 // G4double Yield2 = yieldVector->Value(PreStepKineticEnergy -
1010 // StepEnergyDeposit); ScintillationYield = Yield1 - Yield2;
1011 ScintillationYield =
1012 yieldVector->Value(PreStepKineticEnergy) -
1013 yieldVector->Value(PreStepKineticEnergy - StepEnergyDeposit);
1014 }
1015 else
1016 {
1018 ed << "\nG4Scintillation::GetScintillationYieldByParticleType(): Request\n"
1019 << "for scintillation light yield above the available energy range\n"
1020 << "specified in G4MaterialPropertiesTable. A linear interpolation\n"
1021 << "will be performed to compute the scintillation light yield using\n"
1022 << "(L_max / E_max) as the photon yield per unit energy." << G4endl;
1023 G4String cmt = "\nScintillation yield may be unphysical!\n";
1024 G4Exception("G4Scintillation::GetScintillationYieldByParticleType()",
1025 "Scint03", JustWarning, ed, cmt);
1026
1027 // Units: [# scintillation photons]
1028 ScintillationYield = yieldVector->GetMaxValue() /
1029 yieldVector->GetMaxEnergy() * StepEnergyDeposit;
1030 }
1031
1032#ifdef G4DEBUG_SCINTILLATION
1033 // Increment track aggregators
1034 ScintTrackYield += ScintillationYield;
1035 ScintTrackEDep += StepEnergyDeposit;
1036
1037 G4cout << "\n--- G4Scintillation::GetScintillationYieldByParticleType() ---\n"
1038 << "--\n"
1039 << "-- Name = "
1040 << aTrack.GetParticleDefinition()->GetParticleName() << "\n"
1041 << "-- TrackID = " << aTrack.GetTrackID() << "\n"
1042 << "-- ParentID = " << aTrack.GetParentID() << "\n"
1043 << "-- Current KE = " << aTrack.GetKineticEnergy() / MeV
1044 << " MeV\n"
1045 << "-- Step EDep = " << aStep.GetTotalEnergyDeposit() / MeV
1046 << " MeV\n"
1047 << "-- Track EDep = " << ScintTrackEDep / MeV << " MeV\n"
1048 << "-- Vertex KE = " << aTrack.GetVertexKineticEnergy() / MeV
1049 << " MeV\n"
1050 << "-- Step yield = " << ScintillationYield << " photons\n"
1051 << "-- Track yield = " << ScintTrackYield << " photons\n"
1052 << G4endl;
1053
1054 // The track has terminated within or has left the scintillator volume
1055 if((aTrack.GetTrackStatus() == fStopButAlive) or
1057 {
1058 // Reset aggregators for the next track
1059 ScintTrackEDep = 0.;
1060 ScintTrackYield = 0.;
1061 }
1062#endif
1063
1064 return ScintillationYield;
1065}
@ kTRITONSCINTILLATIONYIELD1
@ kDEUTERONSCINTILLATIONYIELD3
@ kIONSCINTILLATIONYIELD1
@ kDEUTERONSCINTILLATIONYIELD2
@ kTRITONSCINTILLATIONYIELD2
@ kALPHASCINTILLATIONYIELD2
@ kELECTRONSCINTILLATIONYIELD3
@ kALPHASCINTILLATIONYIELD1
@ kELECTRONSCINTILLATIONYIELD2
@ kIONSCINTILLATIONYIELD2
@ kIONSCINTILLATIONYIELD3
@ kPROTONSCINTILLATIONYIELD2
@ kDEUTERONSCINTILLATIONYIELD1
@ kTRITONSCINTILLATIONYIELD3
@ kPROTONSCINTILLATIONYIELD3
@ kELECTRONSCINTILLATIONYIELD1
@ kALPHASCINTILLATIONYIELD3
@ kPROTONSCINTILLATIONYIELD1
G4bool ConstPropertyExists(const G4String &key) const
G4double GetConstProperty(const G4String &key) const

◆ GetScintillationYieldFactor()

G4double G4Scintillation::GetScintillationYieldFactor ( ) const
inline

Definition at line 294 of file G4Scintillation.hh.

295{
296 return fYieldFactor;
297}

◆ GetSlowIntegralTable()

G4PhysicsTable * G4Scintillation::GetSlowIntegralTable ( ) const
inline

Definition at line 310 of file G4Scintillation.hh.

311{
312 return fIntegralTable3;
313}

◆ GetStackPhotons()

G4bool G4Scintillation::GetStackPhotons ( ) const
inline

Definition at line 377 of file G4Scintillation.hh.

377{ return fStackingFlag; }

◆ GetTrackSecondariesFirst()

G4bool G4Scintillation::GetTrackSecondariesFirst ( ) const
inline

Definition at line 273 of file G4Scintillation.hh.

274{
275 return fTrackSecondariesFirst;
276}

◆ Initialise()

void G4Scintillation::Initialise ( )

Definition at line 138 of file G4Scintillation.cc.

139{
150}
G4int GetScintVerboseLevel() const
G4bool GetScintStackPhotons() const
static G4OpticalParameters * Instance()
G4double GetScintExcitationRatio() const
G4bool GetScintEnhancedTimeConstants() const
G4bool GetScintByParticleType() const
G4bool GetScintFiniteRiseTime() const
G4double GetScintYieldFactor() const
G4bool GetScintTrackInfo() const
G4bool GetScintTrackSecondariesFirst() const
void SetScintillationYieldFactor(const G4double yieldfactor)
void SetTrackSecondariesFirst(const G4bool state)
void SetStackPhotons(const G4bool)
void SetScintillationExcitationRatio(const G4double ratio)
void SetScintillationTrackInfo(const G4bool trackType)
void SetEnhancedTimeConstants(G4bool)
void SetFiniteRiseTime(const G4bool state)
void SetScintillationByParticleType(const G4bool)
void SetVerboseLevel(G4int value)
Definition: G4VProcess.hh:412

Referenced by G4Scintillation(), and PreparePhysicsTable().

◆ IsApplicable()

G4bool G4Scintillation::IsApplicable ( const G4ParticleDefinition aParticleType)
overridevirtual

Reimplemented from G4VProcess.

Definition at line 124 of file G4Scintillation.cc.

125{
126 return (aParticleType.GetPDGCharge() == 0.0 || aParticleType.IsShortLived())
127 ? false
128 : true;
129}
G4double GetPDGCharge() const

Referenced by LBE::ConstructOp(), and G4OpticalPhysics::ConstructProcess().

◆ PostStepDoIt()

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

N_timeconstants != 1 and still scnt == 0

Reimplemented from G4VRestDiscreteProcess.

Definition at line 328 of file G4Scintillation.cc.

334{
336 fNumPhotons = 0;
337
338 const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
339 const G4Material* aMaterial = aTrack.GetMaterial();
340
341 G4StepPoint* pPreStepPoint = aStep.GetPreStepPoint();
342 G4StepPoint* pPostStepPoint = aStep.GetPostStepPoint();
343
344 G4ThreeVector x0 = pPreStepPoint->GetPosition();
345 G4ThreeVector p0 = aStep.GetDeltaPosition().unit();
346 G4double t0 = pPreStepPoint->GetGlobalTime();
347
348 G4double TotalEnergyDeposit = aStep.GetTotalEnergyDeposit();
349
351 if(!MPT)
352 return G4VRestDiscreteProcess::PostStepDoIt(aTrack, aStep);
353
354 G4int N_timeconstants = 1;
355
356 // only needed for old (two time constants) version
357 G4MaterialPropertyVector* Fast_Intensity = nullptr;
358 G4MaterialPropertyVector* Slow_Intensity = nullptr;
359
360 if(fEnhancedTimeConstants)
361 {
363 N_timeconstants = 3;
365 N_timeconstants = 2;
366 else if(!(MPT->GetProperty(kSCINTILLATIONCOMPONENT1)))
367 {
368 // no components were specified
369 return G4VRestDiscreteProcess::PostStepDoIt(aTrack, aStep);
370 }
371 }
372 else
373 { // OLD METHOD
374 Fast_Intensity = MPT->GetProperty(kFASTCOMPONENT);
375 Slow_Intensity = MPT->GetProperty(kSLOWCOMPONENT);
376 if(!Fast_Intensity && !Slow_Intensity)
377 return G4VRestDiscreteProcess::PostStepDoIt(aTrack, aStep);
378 if(Fast_Intensity && Slow_Intensity)
379 N_timeconstants = 2;
380 }
381
382 G4double ResolutionScale = MPT->GetConstProperty(kRESOLUTIONSCALE);
383 G4double MeanNumberOfPhotons;
384
385 G4double yield1 = 0.;
386 G4double yield2 = 0.;
387 G4double yield3 = 0.;
388 G4double sum_yields = 0.;
389
390 if(!fEnhancedTimeConstants)
391 {
392 // Scintillation depends on particle type, energy deposited
393 if(fScintillationByParticleType)
394 {
395 MeanNumberOfPhotons = GetScintillationYieldByParticleType(aTrack, aStep);
396 }
397 else
398 {
399 // The default linear scintillation process
400 // Units: [# scintillation photons / MeV]
401 MeanNumberOfPhotons =
402 MPT->GetConstProperty(kSCINTILLATIONYIELD) * fYieldFactor;
403 // Birk's correction via fEmSaturation and specifying scintillation by
404 // by particle type are physically mutually exclusive
405 if(fEmSaturation)
406 MeanNumberOfPhotons *=
407 (fEmSaturation->VisibleEnergyDepositionAtAStep(&aStep));
408 else
409 MeanNumberOfPhotons *= TotalEnergyDeposit;
410 }
411 }
412
413 else
414 {
415 if(fScintillationByParticleType)
416 {
417 MeanNumberOfPhotons = GetScintillationYieldByParticleType(
418 aTrack, aStep, yield1, yield2, yield3);
419 }
420 else
421 {
424 : 1.;
427 : 0.;
430 : 0.;
431 // The default linear scintillation process
432 // Units: [# scintillation photons / MeV]
433 MeanNumberOfPhotons =
434 MPT->GetConstProperty(kSCINTILLATIONYIELD) * fYieldFactor;
435 // Birk's correction via fEmSaturation and specifying scintillation by
436 // by particle type are physically mutually exclusive
437 if(fEmSaturation)
438 MeanNumberOfPhotons *=
439 (fEmSaturation->VisibleEnergyDepositionAtAStep(&aStep));
440 else
441 MeanNumberOfPhotons *= TotalEnergyDeposit;
442 }
443 sum_yields = yield1 + yield2 + yield3;
444 }
445
446 if(MeanNumberOfPhotons > 10.)
447 {
448 G4double sigma = ResolutionScale * std::sqrt(MeanNumberOfPhotons);
449 fNumPhotons = G4int(G4RandGauss::shoot(MeanNumberOfPhotons, sigma) + 0.5);
450 }
451 else
452 {
453 fNumPhotons = G4int(G4Poisson(MeanNumberOfPhotons));
454 }
455
456 if(fNumPhotons <= 0 || !fStackingFlag)
457 {
458 // return unchanged particle and no secondaries
460 return G4VRestDiscreteProcess::PostStepDoIt(aTrack, aStep);
461 }
462
464
465 if(fTrackSecondariesFirst)
466 {
467 if(aTrack.GetTrackStatus() == fAlive)
469 }
470
471 G4int materialIndex = aMaterial->GetIndex();
472
473 // Retrieve the Scintillation Integral for this material
474 // new G4PhysicsOrderedFreeVector allocated to hold CII's
475 size_t numPhot = fNumPhotons;
476 G4double scintTime = 0.;
477 G4double riseTime = 0.;
478 G4PhysicsOrderedFreeVector* scintIntegral = nullptr;
479 G4ScintillationType scintType = Slow;
480
481 for(G4int scnt = 0; scnt < N_timeconstants; ++scnt)
482 {
483 // Original method with FAST and SLOW
484 if(!fEnhancedTimeConstants)
485 {
486 if(scnt == 0)
487 {
488 if(N_timeconstants == 1)
489 {
490 if(Fast_Intensity)
491 {
492 scintTime = MPT->GetConstProperty(kFASTTIMECONSTANT);
493 if(fFiniteRiseTime)
494 {
496 }
497 scintType = Fast;
498 scintIntegral =
499 (G4PhysicsOrderedFreeVector*) ((*fIntegralTable1)(materialIndex));
500 }
501 if(Slow_Intensity)
502 {
503 scintTime = MPT->GetConstProperty(kSLOWTIMECONSTANT);
504 if(fFiniteRiseTime)
505 {
507 }
508 scintType = Slow;
509 scintIntegral =
510 (G4PhysicsOrderedFreeVector*) ((*fIntegralTable3)(materialIndex));
511 }
512 }
513 else
514 { /// N_timeconstants != 1 and still scnt == 0
515 G4double yieldRatio = MPT->GetConstProperty(kYIELDRATIO);
516 if(fExcitationRatio == 1.0 || fExcitationRatio == 0.0)
517 {
518 numPhot = G4int(std::min(yieldRatio, 1.0) * fNumPhotons);
519 }
520 else
521 {
522 numPhot = G4int(std::min(fExcitationRatio, 1.0) * fNumPhotons);
523 }
524 scintTime = MPT->GetConstProperty(kFASTTIMECONSTANT);
525 if(fFiniteRiseTime)
526 {
528 }
529 scintType = Fast;
530 scintIntegral =
531 (G4PhysicsOrderedFreeVector*) ((*fIntegralTable1)(materialIndex));
532 }
533 }
534 else
535 { // scnt != 0
536 numPhot = fNumPhotons - numPhot;
537 scintTime = MPT->GetConstProperty(kSLOWTIMECONSTANT);
538 if(fFiniteRiseTime)
539 {
541 }
542 scintType = Slow;
543 scintIntegral =
544 (G4PhysicsOrderedFreeVector*) ((*fIntegralTable3)(materialIndex));
545 }
546 }
547 else
548 { // fEnhancedTimeConstants == true
549 // in the new method, if there is 1 time constant it is #1, etc.
550 // Note: fExcitationRatio is not used
551 if(scnt == 0)
552 {
553 if(N_timeconstants == 1)
554 {
555 numPhot = fNumPhotons;
556 }
557 else
558 {
559 numPhot = yield1 / sum_yields * fNumPhotons;
560 }
562 if(fFiniteRiseTime)
563 {
565 }
566 scintType = Fast;
567 scintIntegral =
568 (G4PhysicsOrderedFreeVector*) ((*fIntegralTable1)(materialIndex));
569 }
570 else if(scnt == 1)
571 {
572 // to be consistent with old version (due to double->int conversion)
573 if(N_timeconstants == 2)
574 {
575 numPhot = fNumPhotons - numPhot;
576 }
577 else
578 {
579 numPhot = yield2 / sum_yields * fNumPhotons;
580 }
582 if(fFiniteRiseTime)
583 {
585 }
586 scintType = Medium;
587 scintIntegral =
588 (G4PhysicsOrderedFreeVector*) ((*fIntegralTable2)(materialIndex));
589 }
590 else if(scnt == 2)
591 {
592 numPhot = yield3 / sum_yields * fNumPhotons;
594 if(fFiniteRiseTime)
595 {
597 }
598 scintType = Slow;
599 scintIntegral =
600 (G4PhysicsOrderedFreeVector*) ((*fIntegralTable3)(materialIndex));
601 }
602 }
603
604 if(!scintIntegral)
605 continue;
606
607 G4double CIImax = scintIntegral->GetMaxValue();
608 for(size_t i = 0; i < numPhot; ++i)
609 {
610 // Determine photon energy
611 G4double CIIvalue = G4UniformRand() * CIImax;
612 G4double sampledEnergy = scintIntegral->GetEnergy(CIIvalue);
613
614 if(verboseLevel > 1)
615 {
616 G4cout << "sampledEnergy = " << sampledEnergy << G4endl;
617 G4cout << "CIIvalue = " << CIIvalue << G4endl;
618 }
619
620 // Generate random photon direction
621 G4double cost = 1. - 2. * G4UniformRand();
622 G4double sint = std::sqrt((1. - cost) * (1. + cost));
623 G4double phi = twopi * G4UniformRand();
624 G4double sinp = std::sin(phi);
625 G4double cosp = std::cos(phi);
626 G4ParticleMomentum photonMomentum(sint * cosp, sint * sinp, cost);
627
628 // Determine polarization of new photon
629 G4ThreeVector photonPolarization(cost * cosp, cost * sinp, -sint);
630 G4ThreeVector perp = photonMomentum.cross(photonPolarization);
631 phi = twopi * G4UniformRand();
632 sinp = std::sin(phi);
633 cosp = std::cos(phi);
634 photonPolarization = (cosp * photonPolarization + sinp * perp).unit();
635
636 // Generate a new photon:
637 G4DynamicParticle* scintPhoton =
638 new G4DynamicParticle(opticalphoton, photonMomentum);
639 scintPhoton->SetPolarization(photonPolarization);
640 scintPhoton->SetKineticEnergy(sampledEnergy);
641
642 // Generate new G4Track object:
643 G4double rand = G4UniformRand();
644 if(aParticle->GetDefinition()->GetPDGCharge() == 0)
645 {
646 rand = 1.0;
647 }
648
649 // emission time distribution
650 G4double delta = rand * aStep.GetStepLength();
651 G4double deltaTime =
652 delta /
653 (pPreStepPoint->GetVelocity() +
654 rand * (pPostStepPoint->GetVelocity() - pPreStepPoint->GetVelocity()) /
655 2.);
656 if(riseTime == 0.0)
657 {
658 deltaTime -= scintTime * std::log(G4UniformRand());
659 }
660 else
661 {
662 deltaTime += sample_time(riseTime, scintTime);
663 }
664
665 G4double secTime = t0 + deltaTime;
666 G4ThreeVector secPosition = x0 + rand * aStep.GetDeltaPosition();
667
668 G4Track* secTrack = new G4Track(scintPhoton, secTime, secPosition);
669 secTrack->SetTouchableHandle(
671 secTrack->SetParentID(aTrack.GetTrackID());
672 if(fScintillationTrackInfo)
673 secTrack->SetUserInformation(
674 new G4ScintillationTrackInformation(scintType));
676 }
677 }
678
679 if(verboseLevel > 1)
680 {
681 G4cout << "\n Exiting from G4Scintillation::DoIt -- NumberOfSecondaries = "
683 }
684
685 return G4VRestDiscreteProcess::PostStepDoIt(aTrack, aStep);
686}
@ kSCINTILLATIONTIMECONSTANT1
@ kSCINTILLATIONRISETIME2
@ kSCINTILLATIONRISETIME1
@ kFASTSCINTILLATIONRISETIME
@ kSLOWSCINTILLATIONRISETIME
@ kSCINTILLATIONRISETIME3
@ kSCINTILLATIONTIMECONSTANT3
@ kSCINTILLATIONTIMECONSTANT2
G4long G4Poisson(G4double mean)
Definition: G4Poisson.hh:50
@ fSuspend
@ fAlive
int G4int
Definition: G4Types.hh:85
#define G4UniformRand()
Definition: Randomize.hh:52
Hep3Vector unit() const
Hep3Vector cross(const Hep3Vector &) const
void SetPolarization(const G4ThreeVector &)
void SetKineticEnergy(G4double aEnergy)
G4double VisibleEnergyDepositionAtAStep(const G4Step *) const
size_t GetIndex() const
Definition: G4Material.hh:258
void AddSecondary(G4Track *aSecondary)
virtual void Initialize(const G4Track &)
G4double GetEnergy(G4double aValue)
G4double GetScintillationYieldByParticleType(const G4Track &aTrack, const G4Step &aStep)
G4double GetVelocity() const
G4double GetGlobalTime() const
const G4ThreeVector & GetPosition() const
const G4TouchableHandle & GetTouchableHandle() const
G4ThreeVector GetDeltaPosition() const
G4double GetStepLength() const
void SetTouchableHandle(const G4TouchableHandle &apValue)
void SetUserInformation(G4VUserTrackInformation *aValue) const
void SetParentID(const G4int aValue)
void ProposeTrackStatus(G4TrackStatus status)
G4int GetNumberOfSecondaries() const
void SetNumberOfSecondaries(G4int totSecondaries)
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:327
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)

Referenced by AtRestDoIt().

◆ PreparePhysicsTable()

void G4Scintillation::PreparePhysicsTable ( const G4ParticleDefinition part)
overridevirtual

Reimplemented from G4VProcess.

Definition at line 132 of file G4Scintillation.cc.

133{
134 Initialise();
135}

◆ RemoveSaturation()

void G4Scintillation::RemoveSaturation ( )
inline

Definition at line 340 of file G4Scintillation.hh.

340{ fEmSaturation = nullptr; }

Referenced by SetScintillationByParticleType().

◆ SetEnhancedTimeConstants()

void G4Scintillation::SetEnhancedTimeConstants ( G4bool  val)
inline

Definition at line 352 of file G4Scintillation.hh.

353{
354 fEnhancedTimeConstants = val;
355}

Referenced by Initialise().

◆ SetFiniteRiseTime()

void G4Scintillation::SetFiniteRiseTime ( const G4bool  state)
inline

Definition at line 278 of file G4Scintillation.hh.

279{
280 fFiniteRiseTime = state;
281}

Referenced by Initialise().

◆ SetScintillationByParticleType()

void G4Scintillation::SetScintillationByParticleType ( const G4bool  scintType)

Definition at line 689 of file G4Scintillation.cc.

690{
691 if(fEmSaturation && scintType)
692 {
693 G4Exception("G4Scintillation::SetScintillationByParticleType", "Scint02",
695 "Redefinition: Birks Saturation is replaced by "
696 "ScintillationByParticleType!");
698 }
699 fScintillationByParticleType = scintType;
700}

Referenced by Initialise().

◆ SetScintillationExcitationRatio()

void G4Scintillation::SetScintillationExcitationRatio ( const G4double  ratio)
inline

Definition at line 299 of file G4Scintillation.hh.

301{
302 fExcitationRatio = ratio;
303}

Referenced by LBE::ConstructOp(), and Initialise().

◆ SetScintillationTrackInfo()

void G4Scintillation::SetScintillationTrackInfo ( const G4bool  trackType)
inline

Definition at line 362 of file G4Scintillation.hh.

363{
364 fScintillationTrackInfo = trackType;
365}

Referenced by Initialise().

◆ SetScintillationYieldFactor()

void G4Scintillation::SetScintillationYieldFactor ( const G4double  yieldfactor)
inline

Definition at line 288 of file G4Scintillation.hh.

290{
291 fYieldFactor = yieldfactor;
292}

Referenced by LBE::ConstructOp(), and Initialise().

◆ SetStackPhotons()

void G4Scintillation::SetStackPhotons ( const G4bool  stackingFlag)
inline

Definition at line 372 of file G4Scintillation.hh.

373{
374 fStackingFlag = stackingFlag;
375}

Referenced by Initialise().

◆ SetTrackSecondariesFirst()

void G4Scintillation::SetTrackSecondariesFirst ( const G4bool  state)
inline

Definition at line 268 of file G4Scintillation.hh.

269{
270 fTrackSecondariesFirst = state;
271}

Referenced by LBE::ConstructOp(), and Initialise().

Member Data Documentation

◆ fIntegralTable1

G4PhysicsTable* G4Scintillation::fIntegralTable1
protected

◆ fIntegralTable2

G4PhysicsTable* G4Scintillation::fIntegralTable2
protected

◆ fIntegralTable3

G4PhysicsTable* G4Scintillation::fIntegralTable3
protected

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