Geant4 11.2.2
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, G4double &timeconstant1, G4double &timeconstant2, G4double &timeconstant3)
 
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 G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
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
 
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 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
 
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)
 
- Protected Member Functions inherited from G4VRestDiscreteProcess
- 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:67
G4GLOB_DLL std::ostream G4cout
static G4int GetModelID(const G4int modelIndex)
G4int verboseLevel
void SetProcessSubType(G4int)
const G4String & GetProcessName() const

◆ ~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 250 of file G4Scintillation.hh.

251{
252 fEmSaturation = sat;
253}

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 std::size_t GetNumberOfMaterials()
static G4MaterialTable * GetMaterialTable()
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 953 of file G4Scintillation.cc.

954{
955 if(fIntegralTable1)
956 {
957 for(std::size_t i = 0; i < fIntegralTable1->entries(); ++i)
958 {
959 ((G4PhysicsFreeVector*) (*fIntegralTable1)[i])->DumpValues();
960 }
961 }
962 if(fIntegralTable2)
963 {
964 for(std::size_t i = 0; i < fIntegralTable2->entries(); ++i)
965 {
966 ((G4PhysicsFreeVector*) (*fIntegralTable2)[i])->DumpValues();
967 }
968 }
969 if(fIntegralTable3)
970 {
971 for(std::size_t i = 0; i < fIntegralTable3->entries(); ++i)
972 {
973 ((G4PhysicsFreeVector*) (*fIntegralTable3)[i])->DumpValues();
974 }
975 }
976}
std::size_t entries() const

◆ GetFiniteRiseTime()

G4bool G4Scintillation::GetFiniteRiseTime ( ) const
inline

Definition at line 230 of file G4Scintillation.hh.

231{
232 return fFiniteRiseTime;
233}

◆ GetIntegralTable1()

G4PhysicsTable * G4Scintillation::GetIntegralTable1 ( ) const
inline

Definition at line 235 of file G4Scintillation.hh.

236{
237 return fIntegralTable1;
238}

◆ GetIntegralTable2()

G4PhysicsTable * G4Scintillation::GetIntegralTable2 ( ) const
inline

Definition at line 240 of file G4Scintillation.hh.

241{
242 return fIntegralTable2;
243}

◆ GetIntegralTable3()

G4PhysicsTable * G4Scintillation::GetIntegralTable3 ( ) const
inline

Definition at line 245 of file G4Scintillation.hh.

246{
247 return fIntegralTable3;
248}

◆ GetMeanFreePath()

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

Implements G4VRestDiscreteProcess.

Definition at line 622 of file G4Scintillation.cc.

624{
626 return DBL_MAX;
627}
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 630 of file G4Scintillation.cc.

632{
633 *condition = Forced;
634 return DBL_MAX;
635}
@ Forced

◆ GetNumPhotons()

G4int G4Scintillation::GetNumPhotons ( ) const
inline

Definition at line 274 of file G4Scintillation.hh.

274{ return fNumPhotons; }

◆ GetSaturation()

G4EmSaturation * G4Scintillation::GetSaturation ( ) const
inline

Definition at line 257 of file G4Scintillation.hh.

258{
259 return fEmSaturation;
260}

◆ GetScintillationByParticleType()

G4bool G4Scintillation::GetScintillationByParticleType ( ) const
inline

Definition at line 262 of file G4Scintillation.hh.

263{
264 return fScintillationByParticleType;
265}

◆ GetScintillationTrackInfo()

G4bool G4Scintillation::GetScintillationTrackInfo ( ) const
inline

Definition at line 267 of file G4Scintillation.hh.

268{
269 return fScintillationTrackInfo;
270}

◆ GetScintillationYieldByParticleType()

G4double G4Scintillation::GetScintillationYieldByParticleType ( const G4Track & aTrack,
const G4Step & aStep,
G4double & yield1,
G4double & yield2,
G4double & yield3,
G4double & timeconstant1,
G4double & timeconstant2,
G4double & timeconstant3 )

Definition at line 655 of file G4Scintillation.cc.

659{
660 // new in 10.7, allow multiple time constants with ScintByParticleType
661 // Get the G4MaterialPropertyVector containing the scintillation
662 // yield as a function of the energy deposited and particle type
663 // In 11.2, allow different time constants for different particles
664
666 G4MaterialPropertyVector* yieldVector = nullptr;
669
670 // Protons
671 if(pDef == G4Proton::ProtonDefinition())
672 {
673 yieldVector = MPT->GetProperty(kPROTONSCINTILLATIONYIELD);
676 : 1.;
679 : 0.;
682 : 0.;
686 if(yield2 > 0.)
687 {
691 }
692 if(yield3 > 0.)
693 {
697 }
698 }
699
700 // Deuterons
701 else if(pDef == G4Deuteron::DeuteronDefinition())
702 {
703 yieldVector = MPT->GetProperty(kDEUTERONSCINTILLATIONYIELD);
706 : 1.;
709 : 0.;
712 : 0.;
716 if(yield2 > 0.)
717 {
721 }
722 if(yield3 > 0.)
723 {
727 }
728 }
729
730 // Tritons
731 else if(pDef == G4Triton::TritonDefinition())
732 {
733 yieldVector = MPT->GetProperty(kTRITONSCINTILLATIONYIELD);
736 : 1.;
739 : 0.;
742 : 0.;
746 if(yield2 > 0.)
747 {
751 }
752 if(yield3 > 0.)
753 {
757 }
758 }
759
760 // Alphas
761 else if(pDef == G4Alpha::AlphaDefinition())
762 {
763 yieldVector = MPT->GetProperty(kALPHASCINTILLATIONYIELD);
766 : 1.;
769 : 0.;
772 : 0.;
776 if(yield2 > 0.)
777 {
781 }
782 if(yield3 > 0.)
783 {
787 }
788 }
789
790 // Ions (particles derived from G4VIon and G4Ions) and recoil ions
791 // below the production cut from neutrons after hElastic
792 else if(pDef->GetParticleType() == "nucleus" ||
794 {
795 yieldVector = MPT->GetProperty(kIONSCINTILLATIONYIELD);
798 : 1.;
801 : 0.;
804 : 0.;
808 if(yield2 > 0.)
809 {
813 }
814 if(yield3 > 0.)
815 {
819 }
820 }
821
822 // Electrons (must also account for shell-binding energy
823 // attributed to gamma from standard photoelectric effect)
824 // and, default for particles not enumerated/listed above
825 else
826 {
827 yieldVector = MPT->GetProperty(kELECTRONSCINTILLATIONYIELD);
830 : 1.;
833 : 0.;
836 : 0.;
840 if(yield2 > 0.)
841 {
845 }
846 if(yield3 > 0.)
847 {
851 }
852 }
853
854 // Throw an exception if no scintillation yield vector is found
855 if(yieldVector == nullptr)
856 {
858 ed << "\nG4Scintillation::PostStepDoIt(): "
859 << "Request for scintillation yield for energy deposit and particle\n"
860 << "type without correct entry in MaterialPropertiesTable. A material\n"
861 << "property (vector) with name like PARTICLESCINTILLATIONYIELD is\n"
862 << "needed (hint: PARTICLE might not be the primary particle."
863 << G4endl;
864 G4String comments = "Missing MaterialPropertiesTable entry - No correct "
865 "entry in MaterialPropertiesTable";
866 G4Exception("G4Scintillation::PostStepDoIt", "Scint01", FatalException, ed,
867 comments);
868 return 0.; // NOLINT: required to help Coverity recognise this as exit point
869 }
870
871 ///////////////////////////////////////
872 // Calculate the scintillation light //
873 ///////////////////////////////////////
874 // To account for potential nonlinearity and scintillation photon
875 // density along the track, light (L) is produced according to:
876 // L_currentStep = L(PreStepKE) - L(PreStepKE - EDep)
877
878 G4double ScintillationYield = 0.;
879 G4double StepEnergyDeposit = aStep.GetTotalEnergyDeposit();
880 G4double PreStepKineticEnergy = aStep.GetPreStepPoint()->GetKineticEnergy();
881
882 if(PreStepKineticEnergy <= yieldVector->GetMaxEnergy())
883 {
884 // G4double Yield1 = yieldVector->Value(PreStepKineticEnergy);
885 // G4double Yield2 = yieldVector->Value(PreStepKineticEnergy -
886 // StepEnergyDeposit); ScintillationYield = Yield1 - Yield2;
887 ScintillationYield =
888 yieldVector->Value(PreStepKineticEnergy) -
889 yieldVector->Value(PreStepKineticEnergy - StepEnergyDeposit);
890 }
891 else
892 {
893 ++fNumEnergyWarnings;
894 if(verboseLevel > 0 && fNumEnergyWarnings <= 10)
895 {
897 ed << "\nG4Scintillation::GetScintillationYieldByParticleType(): Request\n"
898 << "for scintillation light yield above the available energy range\n"
899 << "specified in G4MaterialPropertiesTable. A linear interpolation\n"
900 << "will be performed to compute the scintillation light yield using\n"
901 << "(L_max / E_max) as the photon yield per unit energy." << G4endl;
902 G4String cmt = "\nScintillation yield may be unphysical!\n";
903
904 if(fNumEnergyWarnings == 10)
905 {
906 ed << G4endl << "*** Scintillation energy warnings stopped.";
907 }
908 G4Exception("G4Scintillation::GetScintillationYieldByParticleType()",
909 "Scint03", JustWarning, ed, cmt);
910 }
911
912 // Units: [# scintillation photons]
913 ScintillationYield = yieldVector->GetMaxValue() /
914 yieldVector->GetMaxEnergy() * StepEnergyDeposit;
915 }
916
917#ifdef G4DEBUG_SCINTILLATION
918 // Increment track aggregators
919 ScintTrackYield += ScintillationYield;
920 ScintTrackEDep += StepEnergyDeposit;
921
922 G4cout << "\n--- G4Scintillation::GetScintillationYieldByParticleType() ---\n"
923 << "--\n"
924 << "-- Name = "
925 << aTrack.GetParticleDefinition()->GetParticleName() << "\n"
926 << "-- TrackID = " << aTrack.GetTrackID() << "\n"
927 << "-- ParentID = " << aTrack.GetParentID() << "\n"
928 << "-- Current KE = " << aTrack.GetKineticEnergy() / MeV
929 << " MeV\n"
930 << "-- Step EDep = " << aStep.GetTotalEnergyDeposit() / MeV
931 << " MeV\n"
932 << "-- Track EDep = " << ScintTrackEDep / MeV << " MeV\n"
933 << "-- Vertex KE = " << aTrack.GetVertexKineticEnergy() / MeV
934 << " MeV\n"
935 << "-- Step yield = " << ScintillationYield << " photons\n"
936 << "-- Track yield = " << ScintTrackYield << " photons\n"
937 << G4endl;
938
939 // The track has terminated within or has left the scintillator volume
940 if((aTrack.GetTrackStatus() == fStopButAlive) or
942 {
943 // Reset aggregators for the next track
944 ScintTrackEDep = 0.;
945 ScintTrackYield = 0.;
946 }
947#endif
948
949 return ScintillationYield;
950}
@ JustWarning
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
@ kELECTRONSCINTILLATIONYIELD
@ kALPHASCINTILLATIONYIELD
@ kPROTONSCINTILLATIONYIELD
@ kDEUTERONSCINTILLATIONYIELD
@ kTRITONSCINTILLATIONYIELD
@ kSCINTILLATIONTIMECONSTANT1
@ kTRITONSCINTILLATIONYIELD1
@ kDEUTERONSCINTILLATIONYIELD3
@ kPROTONSCINTILLATIONTIMECONSTANT2
@ kALPHASCINTILLATIONTIMECONSTANT1
@ kELECTRONSCINTILLATIONTIMECONSTANT2
@ kDEUTERONSCINTILLATIONYIELD2
@ kELECTRONSCINTILLATIONTIMECONSTANT3
@ kDEUTERONSCINTILLATIONTIMECONSTANT1
@ kTRITONSCINTILLATIONTIMECONSTANT2
@ kTRITONSCINTILLATIONYIELD2
@ kALPHASCINTILLATIONYIELD2
@ kELECTRONSCINTILLATIONYIELD3
@ kTRITONSCINTILLATIONTIMECONSTANT1
@ kALPHASCINTILLATIONYIELD1
@ kALPHASCINTILLATIONTIMECONSTANT2
@ kPROTONSCINTILLATIONTIMECONSTANT3
@ kDEUTERONSCINTILLATIONTIMECONSTANT3
@ kELECTRONSCINTILLATIONYIELD2
@ kPROTONSCINTILLATIONYIELD2
@ kDEUTERONSCINTILLATIONYIELD1
@ kTRITONSCINTILLATIONTIMECONSTANT3
@ kTRITONSCINTILLATIONYIELD3
@ kSCINTILLATIONTIMECONSTANT3
@ kIONSCINTILLATIONTIMECONSTANT1
@ kIONSCINTILLATIONTIMECONSTANT3
@ kPROTONSCINTILLATIONYIELD3
@ kIONSCINTILLATIONTIMECONSTANT2
@ kALPHASCINTILLATIONTIMECONSTANT3
@ kELECTRONSCINTILLATIONYIELD1
@ kALPHASCINTILLATIONYIELD3
@ kSCINTILLATIONTIMECONSTANT2
@ kPROTONSCINTILLATIONYIELD1
@ kDEUTERONSCINTILLATIONTIMECONSTANT2
@ kPROTONSCINTILLATIONTIMECONSTANT1
@ kELECTRONSCINTILLATIONTIMECONSTANT1
@ fGeomBoundary
@ fStopButAlive
static G4Alpha * AlphaDefinition()
Definition G4Alpha.cc:78
static G4Deuteron * DeuteronDefinition()
Definition G4Deuteron.cc:85
G4ParticleDefinition * GetDefinition() const
G4bool ConstPropertyExists(const G4String &key) const
G4double GetConstProperty(const G4String &key) const
G4MaterialPropertiesTable * GetMaterialPropertiesTable() const
static G4Neutron * NeutronDefinition()
Definition G4Neutron.cc:96
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:85
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:85

Referenced by PostStepDoIt().

◆ GetStackPhotons()

G4bool G4Scintillation::GetStackPhotons ( ) const
inline

Definition at line 272 of file G4Scintillation.hh.

272{ return fStackingFlag; }

◆ GetTrackSecondariesFirst()

G4bool G4Scintillation::GetTrackSecondariesFirst ( ) const
inline

Definition at line 225 of file G4Scintillation.hh.

226{
227 return fTrackSecondariesFirst;
228}

◆ 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 timeconstant1 = 0.;
392 G4double timeconstant2 = 0.;
393 G4double timeconstant3 = 0.;
394 G4double sum_yields = 0.;
395
396 if(fScintillationByParticleType)
397 {
398 MeanNumberOfPhotons = GetScintillationYieldByParticleType(
399 aTrack, aStep, yield1, yield2, yield3, timeconstant1, timeconstant2,
400 timeconstant3);
401 }
402 else
403 {
406 : 1.;
409 : 0.;
412 : 0.;
413 // The default linear scintillation process
414 // Units: [# scintillation photons / MeV]
415 MeanNumberOfPhotons = MPT->GetConstProperty(kSCINTILLATIONYIELD);
416 // Birk's correction via fEmSaturation and specifying scintillation by
417 // by particle type are physically mutually exclusive
418 if(fEmSaturation)
419 MeanNumberOfPhotons *=
420 (fEmSaturation->VisibleEnergyDepositionAtAStep(&aStep));
421 else
422 MeanNumberOfPhotons *= TotalEnergyDeposit;
423 }
424 sum_yields = yield1 + yield2 + yield3;
425
426 if(MeanNumberOfPhotons > 10.)
427 {
428 G4double sigma = ResolutionScale * std::sqrt(MeanNumberOfPhotons);
429 fNumPhotons = G4int(G4RandGauss::shoot(MeanNumberOfPhotons, sigma) + 0.5);
430 }
431 else
432 {
433 fNumPhotons = G4int(G4Poisson(MeanNumberOfPhotons));
434 }
435
436 if(fNumPhotons <= 0 || !fStackingFlag)
437 {
438 // return unchanged particle and no secondaries
440 return G4VRestDiscreteProcess::PostStepDoIt(aTrack, aStep);
441 }
442
444
445 if(fTrackSecondariesFirst)
446 {
447 if(aTrack.GetTrackStatus() == fAlive)
449 }
450
451 G4int materialIndex = (G4int)aMaterial->GetIndex();
452
453 // Retrieve the Scintillation Integral for this material
454 // new G4PhysicsFreeVector allocated to hold CII's
455 std::size_t numPhot = fNumPhotons;
456 G4double scintTime = 0.;
457 G4double riseTime = 0.;
458 G4PhysicsFreeVector* scintIntegral = nullptr;
459 G4ScintillationType scintType = Slow;
460
461 for(G4int scnt = 0; scnt < N_timeconstants; ++scnt)
462 {
463 // if there is 1 time constant it is #1, etc.
464 if(scnt == 0)
465 {
466 if(N_timeconstants == 1)
467 {
468 numPhot = fNumPhotons;
469 }
470 else
471 {
472 numPhot = yield1 / sum_yields * fNumPhotons;
473 }
474 if(fScintillationByParticleType)
475 {
476 scintTime = timeconstant1;
477 }
478 else
479 {
481 }
482 if(fFiniteRiseTime)
483 {
485 }
486 scintType = Fast;
487 scintIntegral =
488 (G4PhysicsFreeVector*) ((*fIntegralTable1)(materialIndex));
489 }
490 else if(scnt == 1)
491 {
492 // to be consistent with old version (due to double->int conversion)
493 if(N_timeconstants == 2)
494 {
495 numPhot = fNumPhotons - numPhot;
496 }
497 else
498 {
499 numPhot = yield2 / sum_yields * fNumPhotons;
500 }
501 if(fScintillationByParticleType)
502 {
503 scintTime = timeconstant2;
504 }
505 else
506 {
508 }
509 if(fFiniteRiseTime)
510 {
512 }
513 scintType = Medium;
514 scintIntegral =
515 (G4PhysicsFreeVector*) ((*fIntegralTable2)(materialIndex));
516 }
517 else if(scnt == 2)
518 {
519 numPhot = yield3 / sum_yields * fNumPhotons;
520 if(fScintillationByParticleType)
521 {
522 scintTime = timeconstant3;
523 }
524 else
525 {
527 }
528 if(fFiniteRiseTime)
529 {
531 }
532 scintType = Slow;
533 scintIntegral =
534 (G4PhysicsFreeVector*) ((*fIntegralTable3)(materialIndex));
535 }
536
537 if(!scintIntegral)
538 continue;
539
540 G4double CIImax = scintIntegral->GetMaxValue();
541 for(std::size_t i = 0; i < numPhot; ++i)
542 {
543 // Determine photon energy
544 G4double CIIvalue = G4UniformRand() * CIImax;
545 G4double sampledEnergy = scintIntegral->GetEnergy(CIIvalue);
546
547 if(verboseLevel > 1)
548 {
549 G4cout << "sampledEnergy = " << sampledEnergy << G4endl;
550 G4cout << "CIIvalue = " << CIIvalue << G4endl;
551 }
552
553 // Generate random photon direction
554 G4double cost = 1. - 2. * G4UniformRand();
555 G4double sint = std::sqrt((1. - cost) * (1. + cost));
556 G4double phi = twopi * G4UniformRand();
557 G4double sinp = std::sin(phi);
558 G4double cosp = std::cos(phi);
559 G4ParticleMomentum photonMomentum(sint * cosp, sint * sinp, cost);
560
561 // Determine polarization of new photon
562 G4ThreeVector photonPolarization(cost * cosp, cost * sinp, -sint);
563 G4ThreeVector perp = photonMomentum.cross(photonPolarization);
564 phi = twopi * G4UniformRand();
565 sinp = std::sin(phi);
566 cosp = std::cos(phi);
567 photonPolarization = (cosp * photonPolarization + sinp * perp).unit();
568
569 // Generate a new photon:
570 auto scintPhoton = new G4DynamicParticle(opticalphoton, photonMomentum);
571 scintPhoton->SetPolarization(photonPolarization);
572 scintPhoton->SetKineticEnergy(sampledEnergy);
573
574 // Generate new G4Track object:
575 G4double rand = G4UniformRand();
576 if(aParticle->GetDefinition()->GetPDGCharge() == 0)
577 {
578 rand = 1.0;
579 }
580
581 // emission time distribution
582 G4double delta = rand * aStep.GetStepLength();
583 G4double deltaTime =
584 delta /
585 (pPreStepPoint->GetVelocity() +
586 rand * (pPostStepPoint->GetVelocity() - pPreStepPoint->GetVelocity()) /
587 2.);
588 if(riseTime == 0.0)
589 {
590 deltaTime -= scintTime * std::log(G4UniformRand());
591 }
592 else
593 {
594 deltaTime += sample_time(riseTime, scintTime);
595 }
596
597 G4double secTime = t0 + deltaTime;
598 G4ThreeVector secPosition = x0 + rand * aStep.GetDeltaPosition();
599
600 G4Track* secTrack = new G4Track(scintPhoton, secTime, secPosition);
601 secTrack->SetTouchableHandle(
603 secTrack->SetParentID(aTrack.GetTrackID());
604 secTrack->SetCreatorModelID(secID);
605 if(fScintillationTrackInfo)
606 secTrack->SetUserInformation(
607 new G4ScintillationTrackInformation(scintType));
609 }
610 }
611
612 if(verboseLevel > 1)
613 {
614 G4cout << "\n Exiting from G4Scintillation::DoIt -- NumberOfSecondaries = "
616 }
617
618 return G4VRestDiscreteProcess::PostStepDoIt(aTrack, aStep);
619}
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
std::size_t GetIndex() const
void AddSecondary(G4Track *aSecondary)
void Initialize(const G4Track &) override
G4double GetEnergy(const G4double value) const
G4double GetScintillationYieldByParticleType(const G4Track &aTrack, const G4Step &aStep, G4double &yield1, G4double &yield2, G4double &yield3, G4double &timeconstant1, G4double &timeconstant2, G4double &timeconstant3)
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
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

Referenced by DumpInfo().

◆ RemoveSaturation()

void G4Scintillation::RemoveSaturation ( )
inline

Definition at line 255 of file G4Scintillation.hh.

255{ fEmSaturation = nullptr; }

Referenced by SetScintillationByParticleType().

◆ SetFiniteRiseTime()

void G4Scintillation::SetFiniteRiseTime ( const G4bool state)

Definition at line 987 of file G4Scintillation.cc.

988{
989 fFiniteRiseTime = state;
991}
void SetScintFiniteRiseTime(G4bool)

Referenced by Initialise().

◆ SetScintillationByParticleType()

void G4Scintillation::SetScintillationByParticleType ( const G4bool scintType)

Definition at line 994 of file G4Scintillation.cc.

995{
996 if(fEmSaturation && scintType)
997 {
998 G4Exception("G4Scintillation::SetScintillationByParticleType", "Scint02",
1000 "Redefinition: Birks Saturation is replaced by "
1001 "ScintillationByParticleType!");
1003 }
1004 fScintillationByParticleType = scintType;
1006 fScintillationByParticleType);
1007}
void SetScintByParticleType(G4bool)

Referenced by Initialise().

◆ SetScintillationTrackInfo()

void G4Scintillation::SetScintillationTrackInfo ( const G4bool trackType)

Definition at line 1010 of file G4Scintillation.cc.

1011{
1012 fScintillationTrackInfo = trackType;
1013 G4OpticalParameters::Instance()->SetScintTrackInfo(fScintillationTrackInfo);
1014}

Referenced by Initialise().

◆ SetStackPhotons()

void G4Scintillation::SetStackPhotons ( const G4bool stackingFlag)

Definition at line 1017 of file G4Scintillation.cc.

1018{
1019 fStackingFlag = stackingFlag;
1021}

Referenced by Initialise().

◆ SetTrackSecondariesFirst()

void G4Scintillation::SetTrackSecondariesFirst ( const G4bool state)

Definition at line 979 of file G4Scintillation.cc.

980{
981 fTrackSecondariesFirst = state;
983 fTrackSecondariesFirst);
984}
void SetScintTrackSecondariesFirst(G4bool)

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

◆ SetVerboseLevel()

void G4Scintillation::SetVerboseLevel ( G4int verbose)

Definition at line 1024 of file G4Scintillation.cc.

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


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