Geant4 11.1.1
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 ()
 
 G4Scintillation (const G4Scintillation &right)=delete
 
G4Scintillationoperator= (const G4Scintillation &right)=delete
 
G4bool IsApplicable (const G4ParticleDefinition &aParticleType) override
 
void ProcessDescription (std::ostream &) const override
 
void DumpInfo () const 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 &yield1, G4double &yield2, G4double &yield3)
 
void SetTrackSecondariesFirst (const G4bool state)
 
G4bool GetTrackSecondariesFirst () const
 
void SetFiniteRiseTime (const G4bool state)
 
G4bool GetFiniteRiseTime () const
 
G4PhysicsTableGetIntegralTable1 () const
 
G4PhysicsTableGetIntegralTable2 () const
 
G4PhysicsTableGetIntegralTable3 () const
 
void AddSaturation (G4EmSaturation *sat)
 
void RemoveSaturation ()
 
G4EmSaturationGetSaturation () const
 
void SetScintillationByParticleType (const G4bool)
 
G4bool GetScintillationByParticleType () const
 
void SetScintillationTrackInfo (const G4bool trackType)
 
G4bool GetScintillationTrackInfo () const
 
void SetStackPhotons (const G4bool)
 
G4bool GetStackPhotons () const
 
G4int GetNumPhotons () const
 
void DumpPhysicsTable () const
 
void SetVerboseLevel (G4int)
 
- 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 const G4VProcessGetCreatorProcess () const
 
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 &)
 

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 ()
 
- 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
 

Detailed Description

Definition at line 67 of file G4Scintillation.hh.

Constructor & Destructor Documentation

◆ G4Scintillation() [1/2]

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

Definition at line 88 of file G4Scintillation.cc.

90 : G4VRestDiscreteProcess(processName, type)
91 , fIntegralTable1(nullptr)
92 , fIntegralTable2(nullptr)
93 , fIntegralTable3(nullptr)
94 , fEmSaturation(nullptr)
95 , fNumPhotons(0)
96{
97 secID = G4PhysicsModelCatalog::GetModelID("model_Scintillation");
99
100#ifdef G4DEBUG_SCINTILLATION
101 ScintTrackEDep = 0.;
102 ScintTrackYield = 0.;
103#endif
104
105 if(verboseLevel > 0)
106 {
107 G4cout << GetProcessName() << " is created " << G4endl;
108 }
109 Initialise();
110}
@ fScintillation
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
static G4int GetModelID(const G4int modelIndex)
G4int verboseLevel
Definition: G4VProcess.hh:360
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:410
const G4String & GetProcessName() const
Definition: G4VProcess.hh:386

◆ ~G4Scintillation()

G4Scintillation::~G4Scintillation ( )

Definition at line 113 of file G4Scintillation.cc.

114{
115 if(fIntegralTable1 != nullptr)
116 {
117 fIntegralTable1->clearAndDestroy();
118 delete fIntegralTable1;
119 }
120 if(fIntegralTable2 != nullptr)
121 {
122 fIntegralTable2->clearAndDestroy();
123 delete fIntegralTable2;
124 }
125 if(fIntegralTable3 != nullptr)
126 {
127 fIntegralTable3->clearAndDestroy();
128 delete fIntegralTable3;
129 }
130}
void clearAndDestroy()

◆ G4Scintillation() [2/2]

G4Scintillation::G4Scintillation ( const G4Scintillation right)
delete

Member Function Documentation

◆ AddSaturation()

void G4Scintillation::AddSaturation ( G4EmSaturation sat)
inline

Definition at line 246 of file G4Scintillation.hh.

247{
248 fEmSaturation = sat;
249}

Referenced by G4OpticalPhysics::ConstructProcess().

◆ AtRestDoIt()

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

Reimplemented from G4VRestDiscreteProcess.

Definition at line 338 of file G4Scintillation.cc.

342{
343 return G4Scintillation::PostStepDoIt(aTrack, aStep);
344}
G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep) override

◆ BuildPhysicsTable()

void G4Scintillation::BuildPhysicsTable ( const G4ParticleDefinition aParticleType)
overridevirtual

Reimplemented from G4VProcess.

Definition at line 178 of file G4Scintillation.cc.

179{
180 if(fIntegralTable1 != nullptr)
181 {
182 fIntegralTable1->clearAndDestroy();
183 delete fIntegralTable1;
184 fIntegralTable1 = nullptr;
185 }
186 if(fIntegralTable2 != nullptr)
187 {
188 fIntegralTable2->clearAndDestroy();
189 delete fIntegralTable2;
190 fIntegralTable2 = nullptr;
191 }
192 if(fIntegralTable3 != nullptr)
193 {
194 fIntegralTable3->clearAndDestroy();
195 delete fIntegralTable3;
196 fIntegralTable3 = nullptr;
197 }
198
199 const G4MaterialTable* materialTable = G4Material::GetMaterialTable();
200 std::size_t numOfMaterials = G4Material::GetNumberOfMaterials();
201
202 // create new physics table
203 if(!fIntegralTable1)
204 fIntegralTable1 = new G4PhysicsTable(numOfMaterials);
205 if(!fIntegralTable2)
206 fIntegralTable2 = new G4PhysicsTable(numOfMaterials);
207 if(!fIntegralTable3)
208 fIntegralTable3 = new G4PhysicsTable(numOfMaterials);
209
210 for(std::size_t i = 0; i < numOfMaterials; ++i)
211 {
212 auto vector1 = new G4PhysicsFreeVector();
213 auto vector2 = new G4PhysicsFreeVector();
214 auto vector3 = new G4PhysicsFreeVector();
215
216 // Retrieve vector of scintillation wavelength intensity for
217 // the material from the material's optical properties table.
219 ((*materialTable)[i])->GetMaterialPropertiesTable();
220
221 if(MPT)
222 {
225 if(MPV)
226 {
227 // Retrieve the first intensity point in vector
228 // of (photon energy, intensity) pairs
229 G4double currentIN = (*MPV)[0];
230 if(currentIN >= 0.0)
231 {
232 // Create first (photon energy, Scintillation Integral pair
233 G4double currentPM = MPV->Energy(0);
234 G4double currentCII = 0.0;
235 vector1->InsertValues(currentPM, currentCII);
236
237 // Set previous values to current ones prior to loop
238 G4double prevPM = currentPM;
239 G4double prevCII = currentCII;
240 G4double prevIN = currentIN;
241
242 // loop over all (photon energy, intensity)
243 // pairs stored for this material
244 for(std::size_t ii = 1; ii < MPV->GetVectorLength(); ++ii)
245 {
246 currentPM = MPV->Energy(ii);
247 currentIN = (*MPV)[ii];
248 currentCII =
249 prevCII + 0.5 * (currentPM - prevPM) * (prevIN + currentIN);
250
251 vector1->InsertValues(currentPM, currentCII);
252
253 prevPM = currentPM;
254 prevCII = currentCII;
255 prevIN = currentIN;
256 }
257 }
258 }
259
261 if(MPV)
262 {
263 // Retrieve the first intensity point in vector
264 // of (photon energy, intensity) pairs
265 G4double currentIN = (*MPV)[0];
266 if(currentIN >= 0.0)
267 {
268 // Create first (photon energy, Scintillation Integral pair
269 G4double currentPM = MPV->Energy(0);
270 G4double currentCII = 0.0;
271 vector2->InsertValues(currentPM, currentCII);
272
273 // Set previous values to current ones prior to loop
274 G4double prevPM = currentPM;
275 G4double prevCII = currentCII;
276 G4double prevIN = currentIN;
277
278 // loop over all (photon energy, intensity)
279 // pairs stored for this material
280 for(std::size_t ii = 1; ii < MPV->GetVectorLength(); ++ii)
281 {
282 currentPM = MPV->Energy(ii);
283 currentIN = (*MPV)[ii];
284 currentCII =
285 prevCII + 0.5 * (currentPM - prevPM) * (prevIN + currentIN);
286
287 vector2->InsertValues(currentPM, currentCII);
288
289 prevPM = currentPM;
290 prevCII = currentCII;
291 prevIN = currentIN;
292 }
293 }
294 }
296 if(MPV)
297 {
298 // Retrieve the first intensity point in vector
299 // of (photon energy, intensity) pairs
300 G4double currentIN = (*MPV)[0];
301 if(currentIN >= 0.0)
302 {
303 // Create first (photon energy, Scintillation Integral pair
304 G4double currentPM = MPV->Energy(0);
305 G4double currentCII = 0.0;
306 vector3->InsertValues(currentPM, currentCII);
307
308 // Set previous values to current ones prior to loop
309 G4double prevPM = currentPM;
310 G4double prevCII = currentCII;
311 G4double prevIN = currentIN;
312
313 // loop over all (photon energy, intensity)
314 // pairs stored for this material
315 for(std::size_t ii = 1; ii < MPV->GetVectorLength(); ++ii)
316 {
317 currentPM = MPV->Energy(ii);
318 currentIN = (*MPV)[ii];
319 currentCII =
320 prevCII + 0.5 * (currentPM - prevPM) * (prevIN + currentIN);
321
322 vector3->InsertValues(currentPM, currentCII);
323
324 prevPM = currentPM;
325 prevCII = currentCII;
326 prevIN = currentIN;
327 }
328 }
329 }
330 }
331 fIntegralTable1->insertAt(i, vector1);
332 fIntegralTable2->insertAt(i, vector2);
333 fIntegralTable3->insertAt(i, vector3);
334 }
335}
@ kSCINTILLATIONCOMPONENT1
@ kSCINTILLATIONCOMPONENT2
@ kSCINTILLATIONCOMPONENT3
std::vector< G4Material * > G4MaterialTable
double G4double
Definition: G4Types.hh:83
G4MaterialPropertyVector * GetProperty(const char *key) const
static size_t GetNumberOfMaterials()
Definition: G4Material.cc:684
static G4MaterialTable * GetMaterialTable()
Definition: G4Material.cc:677
void insertAt(std::size_t, G4PhysicsVector *)
G4double Energy(const std::size_t index) const
std::size_t GetVectorLength() const

◆ DumpInfo()

void G4Scintillation::DumpInfo ( ) const
inlineoverridevirtual

Reimplemented from G4VProcess.

Definition at line 86 of file G4Scintillation.hh.

void ProcessDescription(std::ostream &) const override

◆ DumpPhysicsTable()

void G4Scintillation::DumpPhysicsTable ( ) const

Definition at line 831 of file G4Scintillation.cc.

832{
833 if(fIntegralTable1)
834 {
835 for(std::size_t i = 0; i < fIntegralTable1->entries(); ++i)
836 {
837 ((G4PhysicsFreeVector*) (*fIntegralTable1)[i])->DumpValues();
838 }
839 }
840 if(fIntegralTable2)
841 {
842 for(std::size_t i = 0; i < fIntegralTable2->entries(); ++i)
843 {
844 ((G4PhysicsFreeVector*) (*fIntegralTable2)[i])->DumpValues();
845 }
846 }
847 if(fIntegralTable3)
848 {
849 for(std::size_t i = 0; i < fIntegralTable3->entries(); ++i)
850 {
851 ((G4PhysicsFreeVector*) (*fIntegralTable3)[i])->DumpValues();
852 }
853 }
854}
std::size_t entries() const

◆ GetFiniteRiseTime()

G4bool G4Scintillation::GetFiniteRiseTime ( ) const
inline

Definition at line 226 of file G4Scintillation.hh.

227{
228 return fFiniteRiseTime;
229}

◆ GetIntegralTable1()

G4PhysicsTable * G4Scintillation::GetIntegralTable1 ( ) const
inline

Definition at line 231 of file G4Scintillation.hh.

232{
233 return fIntegralTable1;
234}

◆ GetIntegralTable2()

G4PhysicsTable * G4Scintillation::GetIntegralTable2 ( ) const
inline

Definition at line 236 of file G4Scintillation.hh.

237{
238 return fIntegralTable2;
239}

◆ GetIntegralTable3()

G4PhysicsTable * G4Scintillation::GetIntegralTable3 ( ) const
inline

Definition at line 241 of file G4Scintillation.hh.

242{
243 return fIntegralTable3;
244}

◆ GetMeanFreePath()

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

Implements G4VRestDiscreteProcess.

Definition at line 597 of file G4Scintillation.cc.

599{
601 return DBL_MAX;
602}
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 605 of file G4Scintillation.cc.

607{
608 *condition = Forced;
609 return DBL_MAX;
610}
@ Forced

◆ GetNumPhotons()

G4int G4Scintillation::GetNumPhotons ( ) const
inline

Definition at line 270 of file G4Scintillation.hh.

270{ return fNumPhotons; }

◆ GetSaturation()

G4EmSaturation * G4Scintillation::GetSaturation ( ) const
inline

Definition at line 253 of file G4Scintillation.hh.

254{
255 return fEmSaturation;
256}

◆ GetScintillationByParticleType()

G4bool G4Scintillation::GetScintillationByParticleType ( ) const
inline

Definition at line 258 of file G4Scintillation.hh.

259{
260 return fScintillationByParticleType;
261}

◆ GetScintillationTrackInfo()

G4bool G4Scintillation::GetScintillationTrackInfo ( ) const
inline

Definition at line 263 of file G4Scintillation.hh.

264{
265 return fScintillationTrackInfo;
266}

◆ GetScintillationYieldByParticleType()

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

Definition at line 635 of file G4Scintillation.cc.

638{
639 // new in 10.7, allow multiple time constants with ScintByParticleType
640 // Get the G4MaterialPropertyVector containing the scintillation
641 // yield as a function of the energy deposited and particle type
642
644 G4MaterialPropertyVector* yieldVector = nullptr;
647
648 // Protons
649 if(pDef == G4Proton::ProtonDefinition())
650 {
651 yieldVector = MPT->GetProperty(kPROTONSCINTILLATIONYIELD);
654 : 1.;
657 : 0.;
660 : 0.;
661 }
662
663 // Deuterons
664 else if(pDef == G4Deuteron::DeuteronDefinition())
665 {
666 yieldVector = MPT->GetProperty(kDEUTERONSCINTILLATIONYIELD);
669 : 1.;
672 : 0.;
675 : 0.;
676 }
677
678 // Tritons
679 else if(pDef == G4Triton::TritonDefinition())
680 {
681 yieldVector = MPT->GetProperty(kTRITONSCINTILLATIONYIELD);
684 : 1.;
687 : 0.;
690 : 0.;
691 }
692
693 // Alphas
694 else if(pDef == G4Alpha::AlphaDefinition())
695 {
696 yieldVector = MPT->GetProperty(kALPHASCINTILLATIONYIELD);
699 : 1.;
702 : 0.;
705 : 0.;
706 }
707
708 // Ions (particles derived from G4VIon and G4Ions) and recoil ions
709 // below the production cut from neutrons after hElastic
710 else if(pDef->GetParticleType() == "nucleus" ||
712 {
713 yieldVector = MPT->GetProperty(kIONSCINTILLATIONYIELD);
716 : 1.;
719 : 0.;
722 : 0.;
723 }
724
725 // Electrons (must also account for shell-binding energy
726 // attributed to gamma from standard photoelectric effect)
727 // and, default for particles not enumerated/listed above
728 else
729 {
730 yieldVector = MPT->GetProperty(kELECTRONSCINTILLATIONYIELD);
733 : 1.;
736 : 0.;
739 : 0.;
740 }
741
742 // Throw an exception if no scintillation yield vector is found
743 if(!yieldVector)
744 {
746 ed << "\nG4Scintillation::PostStepDoIt(): "
747 << "Request for scintillation yield for energy deposit and particle\n"
748 << "type without correct entry in MaterialPropertiesTable.\n"
749 << "ScintillationByParticleType requires at minimum that \n"
750 << "ELECTRONSCINTILLATIONYIELD is set by the user\n"
751 << G4endl;
752 G4String comments = "Missing MaterialPropertiesTable entry - No correct "
753 "entry in MaterialPropertiesTable";
754 G4Exception("G4Scintillation::PostStepDoIt", "Scint01", FatalException, ed,
755 comments);
756 }
757
758 ///////////////////////////////////////
759 // Calculate the scintillation light //
760 ///////////////////////////////////////
761 // To account for potential nonlinearity and scintillation photon
762 // density along the track, light (L) is produced according to:
763 // L_currentStep = L(PreStepKE) - L(PreStepKE - EDep)
764
765 G4double ScintillationYield = 0.;
766 G4double StepEnergyDeposit = aStep.GetTotalEnergyDeposit();
767 G4double PreStepKineticEnergy = aStep.GetPreStepPoint()->GetKineticEnergy();
768
769 if(PreStepKineticEnergy <= yieldVector->GetMaxEnergy())
770 {
771 // G4double Yield1 = yieldVector->Value(PreStepKineticEnergy);
772 // G4double Yield2 = yieldVector->Value(PreStepKineticEnergy -
773 // StepEnergyDeposit); ScintillationYield = Yield1 - Yield2;
774 ScintillationYield =
775 yieldVector->Value(PreStepKineticEnergy) -
776 yieldVector->Value(PreStepKineticEnergy - StepEnergyDeposit);
777 }
778 else
779 {
781 ed << "\nG4Scintillation::GetScintillationYieldByParticleType(): Request\n"
782 << "for scintillation light yield above the available energy range\n"
783 << "specified in G4MaterialPropertiesTable. A linear interpolation\n"
784 << "will be performed to compute the scintillation light yield using\n"
785 << "(L_max / E_max) as the photon yield per unit energy." << G4endl;
786 G4String cmt = "\nScintillation yield may be unphysical!\n";
787 G4Exception("G4Scintillation::GetScintillationYieldByParticleType()",
788 "Scint03", JustWarning, ed, cmt);
789
790 // Units: [# scintillation photons]
791 ScintillationYield = yieldVector->GetMaxValue() /
792 yieldVector->GetMaxEnergy() * StepEnergyDeposit;
793 }
794
795#ifdef G4DEBUG_SCINTILLATION
796 // Increment track aggregators
797 ScintTrackYield += ScintillationYield;
798 ScintTrackEDep += StepEnergyDeposit;
799
800 G4cout << "\n--- G4Scintillation::GetScintillationYieldByParticleType() ---\n"
801 << "--\n"
802 << "-- Name = "
803 << aTrack.GetParticleDefinition()->GetParticleName() << "\n"
804 << "-- TrackID = " << aTrack.GetTrackID() << "\n"
805 << "-- ParentID = " << aTrack.GetParentID() << "\n"
806 << "-- Current KE = " << aTrack.GetKineticEnergy() / MeV
807 << " MeV\n"
808 << "-- Step EDep = " << aStep.GetTotalEnergyDeposit() / MeV
809 << " MeV\n"
810 << "-- Track EDep = " << ScintTrackEDep / MeV << " MeV\n"
811 << "-- Vertex KE = " << aTrack.GetVertexKineticEnergy() / MeV
812 << " MeV\n"
813 << "-- Step yield = " << ScintillationYield << " photons\n"
814 << "-- Track yield = " << ScintTrackYield << " photons\n"
815 << G4endl;
816
817 // The track has terminated within or has left the scintillator volume
818 if((aTrack.GetTrackStatus() == fStopButAlive) or
820 {
821 // Reset aggregators for the next track
822 ScintTrackEDep = 0.;
823 ScintTrackYield = 0.;
824 }
825#endif
826
827 return ScintillationYield;
828}
@ JustWarning
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
@ kELECTRONSCINTILLATIONYIELD
@ kALPHASCINTILLATIONYIELD
@ kPROTONSCINTILLATIONYIELD
@ kDEUTERONSCINTILLATIONYIELD
@ kIONSCINTILLATIONYIELD
@ kTRITONSCINTILLATIONYIELD
@ kTRITONSCINTILLATIONYIELD1
@ kDEUTERONSCINTILLATIONYIELD3
@ kIONSCINTILLATIONYIELD1
@ kDEUTERONSCINTILLATIONYIELD2
@ kTRITONSCINTILLATIONYIELD2
@ kALPHASCINTILLATIONYIELD2
@ kELECTRONSCINTILLATIONYIELD3
@ kALPHASCINTILLATIONYIELD1
@ kELECTRONSCINTILLATIONYIELD2
@ kIONSCINTILLATIONYIELD2
@ kIONSCINTILLATIONYIELD3
@ kPROTONSCINTILLATIONYIELD2
@ kDEUTERONSCINTILLATIONYIELD1
@ kTRITONSCINTILLATIONYIELD3
@ kPROTONSCINTILLATIONYIELD3
@ kELECTRONSCINTILLATIONYIELD1
@ kALPHASCINTILLATIONYIELD3
@ kPROTONSCINTILLATIONYIELD1
@ fGeomBoundary
Definition: G4StepStatus.hh:43
@ fStopButAlive
static G4Alpha * AlphaDefinition()
Definition: G4Alpha.cc:83
static G4Deuteron * DeuteronDefinition()
Definition: G4Deuteron.cc:88
G4ParticleDefinition * GetDefinition() const
G4bool ConstPropertyExists(const G4String &key) const
G4double GetConstProperty(const G4String &key) const
G4MaterialPropertiesTable * GetMaterialPropertiesTable() const
Definition: G4Material.hh:251
static G4Neutron * NeutronDefinition()
Definition: G4Neutron.cc:98
const G4String & GetParticleType() const
const G4String & GetParticleName() const
G4double GetMaxEnergy() const
G4double GetMaxValue() const
G4double Value(const G4double energy, std::size_t &lastidx) 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:88

Referenced by PostStepDoIt().

◆ GetStackPhotons()

G4bool G4Scintillation::GetStackPhotons ( ) const
inline

Definition at line 268 of file G4Scintillation.hh.

268{ return fStackingFlag; }

◆ GetTrackSecondariesFirst()

G4bool G4Scintillation::GetTrackSecondariesFirst ( ) const
inline

Definition at line 221 of file G4Scintillation.hh.

222{
223 return fTrackSecondariesFirst;
224}

◆ Initialise()

void G4Scintillation::Initialise ( )

Definition at line 166 of file G4Scintillation.cc.

167{
175}
G4int GetScintVerboseLevel() const
G4bool GetScintStackPhotons() const
static G4OpticalParameters * Instance()
G4bool GetScintByParticleType() const
G4bool GetScintFiniteRiseTime() const
G4bool GetScintTrackInfo() const
G4bool GetScintTrackSecondariesFirst() const
void SetTrackSecondariesFirst(const G4bool state)
void SetStackPhotons(const G4bool)
void SetVerboseLevel(G4int)
void SetScintillationTrackInfo(const G4bool trackType)
void SetFiniteRiseTime(const G4bool state)
void SetScintillationByParticleType(const G4bool)

Referenced by G4Scintillation(), and PreparePhysicsTable().

◆ IsApplicable()

G4bool G4Scintillation::IsApplicable ( const G4ParticleDefinition aParticleType)
overridevirtual

Reimplemented from G4VProcess.

Definition at line 150 of file G4Scintillation.cc.

151{
152 if(aParticleType.GetParticleName() == "opticalphoton")
153 return false;
154 if(aParticleType.IsShortLived())
155 return false;
156 return true;
157}

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

◆ operator=()

G4Scintillation & G4Scintillation::operator= ( const G4Scintillation right)
delete

◆ PostStepDoIt()

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

Reimplemented from G4VRestDiscreteProcess.

Definition at line 347 of file G4Scintillation.cc.

353{
355 fNumPhotons = 0;
356
357 const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
358 const G4Material* aMaterial = aTrack.GetMaterial();
359
360 G4StepPoint* pPreStepPoint = aStep.GetPreStepPoint();
361 G4StepPoint* pPostStepPoint = aStep.GetPostStepPoint();
362
363 G4ThreeVector x0 = pPreStepPoint->GetPosition();
364 G4ThreeVector p0 = aStep.GetDeltaPosition().unit();
365 G4double t0 = pPreStepPoint->GetGlobalTime();
366
367 G4double TotalEnergyDeposit = aStep.GetTotalEnergyDeposit();
368
370 if(!MPT)
371 return G4VRestDiscreteProcess::PostStepDoIt(aTrack, aStep);
372
373 G4int N_timeconstants = 1;
374
376 N_timeconstants = 3;
378 N_timeconstants = 2;
379 else if(!(MPT->GetProperty(kSCINTILLATIONCOMPONENT1)))
380 {
381 // no components were specified
382 return G4VRestDiscreteProcess::PostStepDoIt(aTrack, aStep);
383 }
384
385 G4double ResolutionScale = MPT->GetConstProperty(kRESOLUTIONSCALE);
386 G4double MeanNumberOfPhotons;
387
388 G4double yield1 = 0.;
389 G4double yield2 = 0.;
390 G4double yield3 = 0.;
391 G4double sum_yields = 0.;
392
393 if(fScintillationByParticleType)
394 {
395 MeanNumberOfPhotons = GetScintillationYieldByParticleType(
396 aTrack, aStep, yield1, yield2, yield3);
397 }
398 else
399 {
402 : 1.;
405 : 0.;
408 : 0.;
409 // The default linear scintillation process
410 // Units: [# scintillation photons / MeV]
411 MeanNumberOfPhotons = MPT->GetConstProperty(kSCINTILLATIONYIELD);
412 // Birk's correction via fEmSaturation and specifying scintillation by
413 // by particle type are physically mutually exclusive
414 if(fEmSaturation)
415 MeanNumberOfPhotons *=
416 (fEmSaturation->VisibleEnergyDepositionAtAStep(&aStep));
417 else
418 MeanNumberOfPhotons *= TotalEnergyDeposit;
419 }
420 sum_yields = yield1 + yield2 + yield3;
421
422 if(MeanNumberOfPhotons > 10.)
423 {
424 G4double sigma = ResolutionScale * std::sqrt(MeanNumberOfPhotons);
425 fNumPhotons = G4int(G4RandGauss::shoot(MeanNumberOfPhotons, sigma) + 0.5);
426 }
427 else
428 {
429 fNumPhotons = G4int(G4Poisson(MeanNumberOfPhotons));
430 }
431
432 if(fNumPhotons <= 0 || !fStackingFlag)
433 {
434 // return unchanged particle and no secondaries
436 return G4VRestDiscreteProcess::PostStepDoIt(aTrack, aStep);
437 }
438
440
441 if(fTrackSecondariesFirst)
442 {
443 if(aTrack.GetTrackStatus() == fAlive)
445 }
446
447 G4int materialIndex = (G4int)aMaterial->GetIndex();
448
449 // Retrieve the Scintillation Integral for this material
450 // new G4PhysicsFreeVector allocated to hold CII's
451 std::size_t numPhot = fNumPhotons;
452 G4double scintTime = 0.;
453 G4double riseTime = 0.;
454 G4PhysicsFreeVector* scintIntegral = nullptr;
455 G4ScintillationType scintType = Slow;
456
457 for(G4int scnt = 0; scnt < N_timeconstants; ++scnt)
458 {
459 // if there is 1 time constant it is #1, etc.
460 if(scnt == 0)
461 {
462 if(N_timeconstants == 1)
463 {
464 numPhot = fNumPhotons;
465 }
466 else
467 {
468 numPhot = yield1 / sum_yields * fNumPhotons;
469 }
471 if(fFiniteRiseTime)
472 {
474 }
475 scintType = Fast;
476 scintIntegral =
477 (G4PhysicsFreeVector*) ((*fIntegralTable1)(materialIndex));
478 }
479 else if(scnt == 1)
480 {
481 // to be consistent with old version (due to double->int conversion)
482 if(N_timeconstants == 2)
483 {
484 numPhot = fNumPhotons - numPhot;
485 }
486 else
487 {
488 numPhot = yield2 / sum_yields * fNumPhotons;
489 }
491 if(fFiniteRiseTime)
492 {
494 }
495 scintType = Medium;
496 scintIntegral =
497 (G4PhysicsFreeVector*) ((*fIntegralTable2)(materialIndex));
498 }
499 else if(scnt == 2)
500 {
501 numPhot = yield3 / sum_yields * fNumPhotons;
503 if(fFiniteRiseTime)
504 {
506 }
507 scintType = Slow;
508 scintIntegral =
509 (G4PhysicsFreeVector*) ((*fIntegralTable3)(materialIndex));
510 }
511
512 if(!scintIntegral)
513 continue;
514
515 G4double CIImax = scintIntegral->GetMaxValue();
516 for(std::size_t i = 0; i < numPhot; ++i)
517 {
518 // Determine photon energy
519 G4double CIIvalue = G4UniformRand() * CIImax;
520 G4double sampledEnergy = scintIntegral->GetEnergy(CIIvalue);
521
522 if(verboseLevel > 1)
523 {
524 G4cout << "sampledEnergy = " << sampledEnergy << G4endl;
525 G4cout << "CIIvalue = " << CIIvalue << G4endl;
526 }
527
528 // Generate random photon direction
529 G4double cost = 1. - 2. * G4UniformRand();
530 G4double sint = std::sqrt((1. - cost) * (1. + cost));
531 G4double phi = twopi * G4UniformRand();
532 G4double sinp = std::sin(phi);
533 G4double cosp = std::cos(phi);
534 G4ParticleMomentum photonMomentum(sint * cosp, sint * sinp, cost);
535
536 // Determine polarization of new photon
537 G4ThreeVector photonPolarization(cost * cosp, cost * sinp, -sint);
538 G4ThreeVector perp = photonMomentum.cross(photonPolarization);
539 phi = twopi * G4UniformRand();
540 sinp = std::sin(phi);
541 cosp = std::cos(phi);
542 photonPolarization = (cosp * photonPolarization + sinp * perp).unit();
543
544 // Generate a new photon:
545 auto scintPhoton = new G4DynamicParticle(opticalphoton, photonMomentum);
546 scintPhoton->SetPolarization(photonPolarization);
547 scintPhoton->SetKineticEnergy(sampledEnergy);
548
549 // Generate new G4Track object:
550 G4double rand = G4UniformRand();
551 if(aParticle->GetDefinition()->GetPDGCharge() == 0)
552 {
553 rand = 1.0;
554 }
555
556 // emission time distribution
557 G4double delta = rand * aStep.GetStepLength();
558 G4double deltaTime =
559 delta /
560 (pPreStepPoint->GetVelocity() +
561 rand * (pPostStepPoint->GetVelocity() - pPreStepPoint->GetVelocity()) /
562 2.);
563 if(riseTime == 0.0)
564 {
565 deltaTime -= scintTime * std::log(G4UniformRand());
566 }
567 else
568 {
569 deltaTime += sample_time(riseTime, scintTime);
570 }
571
572 G4double secTime = t0 + deltaTime;
573 G4ThreeVector secPosition = x0 + rand * aStep.GetDeltaPosition();
574
575 G4Track* secTrack = new G4Track(scintPhoton, secTime, secPosition);
576 secTrack->SetTouchableHandle(
578 secTrack->SetParentID(aTrack.GetTrackID());
579 secTrack->SetCreatorModelID(secID);
580 if(fScintillationTrackInfo)
581 secTrack->SetUserInformation(
582 new G4ScintillationTrackInformation(scintType));
584 }
585 }
586
587 if(verboseLevel > 1)
588 {
589 G4cout << "\n Exiting from G4Scintillation::DoIt -- NumberOfSecondaries = "
591 }
592
593 return G4VRestDiscreteProcess::PostStepDoIt(aTrack, aStep);
594}
@ kSCINTILLATIONTIMECONSTANT1
@ kSCINTILLATIONRISETIME2
@ kSCINTILLATIONRISETIME1
@ 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
G4double VisibleEnergyDepositionAtAStep(const G4Step *) const
size_t GetIndex() const
Definition: G4Material.hh:255
void AddSecondary(G4Track *aSecondary)
void Initialize(const G4Track &) override
G4double GetPDGCharge() const
G4double GetEnergy(const G4double value) const
G4double GetScintillationYieldByParticleType(const G4Track &aTrack, const G4Step &aStep, G4double &yield1, G4double &yield2, G4double &yield3)
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 SetCreatorModelID(const G4int id)
void SetParentID(const G4int aValue)
void ProposeTrackStatus(G4TrackStatus status)
G4int GetNumberOfSecondaries() const
void SetNumberOfSecondaries(G4int totSecondaries)
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:331
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)

Referenced by AtRestDoIt().

◆ PreparePhysicsTable()

void G4Scintillation::PreparePhysicsTable ( const G4ParticleDefinition part)
overridevirtual

Reimplemented from G4VProcess.

Definition at line 160 of file G4Scintillation.cc.

161{
162 Initialise();
163}

◆ ProcessDescription()

void G4Scintillation::ProcessDescription ( std::ostream &  out) const
overridevirtual

Reimplemented from G4VProcess.

Definition at line 133 of file G4Scintillation.cc.

134{
135 out << "Scintillation simulates production of optical photons produced\n"
136 "by a high energy particle traversing matter.\n"
137 "Various material properties need to be defined.\n";
139
141 out << "Track secondaries first: " << params->GetScintTrackSecondariesFirst();
142 out << "Finite rise time: " << params->GetScintFiniteRiseTime();
143 out << "Scintillation by particle type: " << params->GetScintByParticleType();
144 out << "Save track information: " << params->GetScintTrackInfo();
145 out << "Stack photons: " << params->GetScintStackPhotons();
146 out << "Verbose level: " << params->GetScintVerboseLevel();
147}
virtual void DumpInfo() const
Definition: G4VProcess.cc:173

Referenced by DumpInfo().

◆ RemoveSaturation()

void G4Scintillation::RemoveSaturation ( )
inline

Definition at line 251 of file G4Scintillation.hh.

251{ fEmSaturation = nullptr; }

Referenced by SetScintillationByParticleType().

◆ SetFiniteRiseTime()

void G4Scintillation::SetFiniteRiseTime ( const G4bool  state)

Definition at line 865 of file G4Scintillation.cc.

866{
867 fFiniteRiseTime = state;
869}
void SetScintFiniteRiseTime(G4bool)

Referenced by Initialise().

◆ SetScintillationByParticleType()

void G4Scintillation::SetScintillationByParticleType ( const G4bool  scintType)

Definition at line 872 of file G4Scintillation.cc.

873{
874 if(fEmSaturation && scintType)
875 {
876 G4Exception("G4Scintillation::SetScintillationByParticleType", "Scint02",
878 "Redefinition: Birks Saturation is replaced by "
879 "ScintillationByParticleType!");
881 }
882 fScintillationByParticleType = scintType;
884 fScintillationByParticleType);
885}
void SetScintByParticleType(G4bool)

Referenced by Initialise().

◆ SetScintillationTrackInfo()

void G4Scintillation::SetScintillationTrackInfo ( const G4bool  trackType)

Definition at line 888 of file G4Scintillation.cc.

889{
890 fScintillationTrackInfo = trackType;
891 G4OpticalParameters::Instance()->SetScintTrackInfo(fScintillationTrackInfo);
892}

Referenced by Initialise().

◆ SetStackPhotons()

void G4Scintillation::SetStackPhotons ( const G4bool  stackingFlag)

Definition at line 895 of file G4Scintillation.cc.

896{
897 fStackingFlag = stackingFlag;
899}
void SetScintStackPhotons(G4bool)

Referenced by Initialise().

◆ SetTrackSecondariesFirst()

void G4Scintillation::SetTrackSecondariesFirst ( const G4bool  state)

Definition at line 857 of file G4Scintillation.cc.

858{
859 fTrackSecondariesFirst = state;
861 fTrackSecondariesFirst);
862}
void SetScintTrackSecondariesFirst(G4bool)

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

◆ SetVerboseLevel()

void G4Scintillation::SetVerboseLevel ( G4int  verbose)

Definition at line 902 of file G4Scintillation.cc.

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


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