Geant4 11.4.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 ()
 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 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 69 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
G4VRestDiscreteProcess(const G4String &, G4ProcessType aType=fNotDefined)

Referenced by G4Scintillation(), and operator=().

◆ ~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}

◆ G4Scintillation() [2/2]

G4Scintillation::G4Scintillation ( const G4Scintillation & right)
delete

Member Function Documentation

◆ AddSaturation()

void G4Scintillation::AddSaturation ( G4EmSaturation * sat)
inline

Definition at line 259 of file G4Scintillation.hh.

260{
261 fEmSaturation = sat;
262}

◆ AtRestDoIt()

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

Reimplemented from G4VRestDiscreteProcess.

Definition at line 249 of file G4Scintillation.cc.

253{
254 return G4Scintillation::PostStepDoIt(aTrack, aStep);
255}
G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep) override

◆ BuildPhysicsTable()

void G4Scintillation::BuildPhysicsTable ( const G4ParticleDefinition & aParticleType)
overridevirtual

Reimplemented from G4VProcess.

Definition at line 174 of file G4Scintillation.cc.

175{
176 const G4MaterialTable* materialTable = G4Material::GetMaterialTable();
177 std::size_t numOfMaterials = G4Material::GetNumberOfMaterials();
178
179 // Find the number of materials that have non-empty material property tables
180 std::size_t numOfMaterialsWithMPT = 0;
181 for(std::size_t i = 0; i < numOfMaterials; ++i)
182 {
183 if(((*materialTable)[i])->GetMaterialPropertiesTable())
184 {
185 ++numOfMaterialsWithMPT;
186 }
187 }
188
189 // create new physics table
190 fIntegralTable1 = new G4PhysicsTable(numOfMaterialsWithMPT);
191 fIntegralTable2 = new G4PhysicsTable(numOfMaterialsWithMPT);
192 fIntegralTable3 = new G4PhysicsTable(numOfMaterialsWithMPT);
193
194 std::size_t indexMPT = 0;
195 for(std::size_t i = 0; i < numOfMaterials; ++i)
196 {
197 // Retrieve vector of scintillation wavelength intensity for
198 // the material from the material's optical properties table.
199 G4MaterialPropertiesTable* MPT =
200 ((*materialTable)[i])->GetMaterialPropertiesTable();
201
202 if(MPT)
203 {
204 auto vector1 = new G4PhysicsFreeVector();
205 auto vector2 = new G4PhysicsFreeVector();
206 auto vector3 = new G4PhysicsFreeVector();
207
208 BuildInverseCdfTable(MPT->GetProperty(kSCINTILLATIONCOMPONENT1), vector1);
209 BuildInverseCdfTable(MPT->GetProperty(kSCINTILLATIONCOMPONENT2), vector2);
210 BuildInverseCdfTable(MPT->GetProperty(kSCINTILLATIONCOMPONENT3), vector3);
211
212 fIntegralTable1->insertAt(indexMPT, vector1);
213 fIntegralTable2->insertAt(indexMPT, vector2);
214 fIntegralTable3->insertAt(indexMPT, vector3);
215
216 fIndexMPT.insert(std::make_pair(i, indexMPT));
217 ++indexMPT;
218 }
219 }
220}
@ kSCINTILLATIONCOMPONENT1
@ kSCINTILLATIONCOMPONENT2
@ kSCINTILLATIONCOMPONENT3
std::vector< G4Material * > G4MaterialTable
G4MaterialPropertyVector * GetProperty(const char *key) const
static std::size_t GetNumberOfMaterials()
static G4MaterialTable * GetMaterialTable()

◆ DumpInfo()

void G4Scintillation::DumpInfo ( ) const
inlineoverridevirtual

Reimplemented from G4VProcess.

Definition at line 88 of file G4Scintillation.hh.

void ProcessDescription(std::ostream &) const override

◆ DumpPhysicsTable()

void G4Scintillation::DumpPhysicsTable ( ) const

Definition at line 872 of file G4Scintillation.cc.

873{
874 if(fIntegralTable1)
875 {
876 for(std::size_t i = 0; i < fIntegralTable1->entries(); ++i)
877 {
878 ((G4PhysicsFreeVector*) (*fIntegralTable1)[i])->DumpValues();
879 }
880 }
881 if(fIntegralTable2)
882 {
883 for(std::size_t i = 0; i < fIntegralTable2->entries(); ++i)
884 {
885 ((G4PhysicsFreeVector*) (*fIntegralTable2)[i])->DumpValues();
886 }
887 }
888 if(fIntegralTable3)
889 {
890 for(std::size_t i = 0; i < fIntegralTable3->entries(); ++i)
891 {
892 ((G4PhysicsFreeVector*) (*fIntegralTable3)[i])->DumpValues();
893 }
894 }
895}

◆ GetFiniteRiseTime()

G4bool G4Scintillation::GetFiniteRiseTime ( ) const
inline

Definition at line 239 of file G4Scintillation.hh.

240{
241 return fFiniteRiseTime;
242}

◆ GetIntegralTable1()

G4PhysicsTable * G4Scintillation::GetIntegralTable1 ( ) const
inline

Definition at line 244 of file G4Scintillation.hh.

245{
246 return fIntegralTable1;
247}

◆ GetIntegralTable2()

G4PhysicsTable * G4Scintillation::GetIntegralTable2 ( ) const
inline

Definition at line 249 of file G4Scintillation.hh.

250{
251 return fIntegralTable2;
252}

◆ GetIntegralTable3()

G4PhysicsTable * G4Scintillation::GetIntegralTable3 ( ) const
inline

Definition at line 254 of file G4Scintillation.hh.

255{
256 return fIntegralTable3;
257}

◆ GetMeanFreePath()

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

Implements G4VRestDiscreteProcess.

Definition at line 541 of file G4Scintillation.cc.

543{
545 return DBL_MAX;
546}
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 549 of file G4Scintillation.cc.

551{
552 *condition = Forced;
553 return DBL_MAX;
554}
@ Forced

◆ GetNumPhotons()

G4int G4Scintillation::GetNumPhotons ( ) const
inline

Definition at line 283 of file G4Scintillation.hh.

283{ return fNumPhotons; }

◆ GetSaturation()

G4EmSaturation * G4Scintillation::GetSaturation ( ) const
inline

Definition at line 266 of file G4Scintillation.hh.

267{
268 return fEmSaturation;
269}

◆ GetScintillationByParticleType()

G4bool G4Scintillation::GetScintillationByParticleType ( ) const
inline

Definition at line 271 of file G4Scintillation.hh.

272{
273 return fScintillationByParticleType;
274}

◆ GetScintillationTrackInfo()

G4bool G4Scintillation::GetScintillationTrackInfo ( ) const
inline

Definition at line 276 of file G4Scintillation.hh.

277{
278 return fScintillationTrackInfo;
279}

◆ 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 574 of file G4Scintillation.cc.

578{
579 // new in 10.7, allow multiple time constants with ScintByParticleType
580 // Get the G4MaterialPropertyVector containing the scintillation
581 // yield as a function of the energy deposited and particle type
582 // In 11.2, allow different time constants for different particles
583
584 G4ParticleDefinition* pDef = aTrack.GetDynamicParticle()->GetDefinition();
585 G4MaterialPropertyVector* yieldVector = nullptr;
586 G4MaterialPropertiesTable* MPT =
588
589 // Protons
590 if(pDef == G4Proton::ProtonDefinition())
591 {
592 yieldVector = MPT->GetProperty(kPROTONSCINTILLATIONYIELD);
595 : 1.;
598 : 0.;
601 : 0.;
605 if(yield2 > 0.)
606 {
610 }
611 if(yield3 > 0.)
612 {
616 }
617 }
618
619 // Deuterons
620 else if(pDef == G4Deuteron::DeuteronDefinition())
621 {
622 yieldVector = MPT->GetProperty(kDEUTERONSCINTILLATIONYIELD);
625 : 1.;
628 : 0.;
631 : 0.;
635 if(yield2 > 0.)
636 {
640 }
641 if(yield3 > 0.)
642 {
646 }
647 }
648
649 // Tritons
650 else if(pDef == G4Triton::TritonDefinition())
651 {
652 yieldVector = MPT->GetProperty(kTRITONSCINTILLATIONYIELD);
655 : 1.;
658 : 0.;
661 : 0.;
665 if(yield2 > 0.)
666 {
670 }
671 if(yield3 > 0.)
672 {
676 }
677 }
678
679 // Alphas
680 else if(pDef == G4Alpha::AlphaDefinition())
681 {
682 yieldVector = MPT->GetProperty(kALPHASCINTILLATIONYIELD);
685 : 1.;
688 : 0.;
691 : 0.;
695 if(yield2 > 0.)
696 {
700 }
701 if(yield3 > 0.)
702 {
706 }
707 }
708
709 // Ions (particles derived from G4VIon and G4Ions) and recoil ions
710 // below the production cut from neutrons after hElastic
711 else if(pDef->GetParticleType() == "nucleus" ||
713 {
714 yieldVector = MPT->GetProperty(kIONSCINTILLATIONYIELD);
717 : 1.;
720 : 0.;
723 : 0.;
727 if(yield2 > 0.)
728 {
732 }
733 if(yield3 > 0.)
734 {
738 }
739 }
740
741 // Electrons (must also account for shell-binding energy
742 // attributed to gamma from standard photoelectric effect)
743 // and, default for particles not enumerated/listed above
744 else
745 {
746 yieldVector = MPT->GetProperty(kELECTRONSCINTILLATIONYIELD);
749 : 1.;
752 : 0.;
755 : 0.;
759 if(yield2 > 0.)
760 {
764 }
765 if(yield3 > 0.)
766 {
770 }
771 }
772
773 // Throw an exception if no scintillation yield vector is found
774 if(yieldVector == nullptr)
775 {
777 ed << "\nG4Scintillation::PostStepDoIt(): "
778 << "Request for scintillation yield for energy deposit and particle\n"
779 << "type without correct entry in MaterialPropertiesTable. A material\n"
780 << "property (vector) with name like PARTICLESCINTILLATIONYIELD is\n"
781 << "needed (hint: PARTICLE might not be the primary particle."
782 << G4endl;
783 G4String comments = "Missing MaterialPropertiesTable entry - No correct "
784 "entry in MaterialPropertiesTable";
785 G4Exception("G4Scintillation::PostStepDoIt", "Scint01", FatalException, ed,
786 comments);
787 return 0.; // NOLINT: required to help Coverity recognise this as exit point
788 }
789
790 ///////////////////////////////////////
791 // Calculate the scintillation light //
792 ///////////////////////////////////////
793 // To account for potential nonlinearity and scintillation photon
794 // density along the track, light (L) is produced according to:
795 // L_currentStep = L(PreStepKE) - L(PreStepKE - EDep)
796
797 G4double ScintillationYield = 0.;
798 G4double StepEnergyDeposit = aStep.GetTotalEnergyDeposit();
799 G4double PreStepKineticEnergy = aStep.GetPreStepPoint()->GetKineticEnergy();
800
801 if(PreStepKineticEnergy <= yieldVector->GetMaxEnergy())
802 {
803 // G4double Yield1 = yieldVector->Value(PreStepKineticEnergy);
804 // G4double Yield2 = yieldVector->Value(PreStepKineticEnergy -
805 // StepEnergyDeposit); ScintillationYield = Yield1 - Yield2;
806 ScintillationYield =
807 yieldVector->Value(PreStepKineticEnergy) -
808 yieldVector->Value(PreStepKineticEnergy - StepEnergyDeposit);
809 }
810 else
811 {
812 ++fNumEnergyWarnings;
813 if(verboseLevel > 0 && fNumEnergyWarnings <= 10)
814 {
816 ed << "\nG4Scintillation::GetScintillationYieldByParticleType(): Request\n"
817 << "for scintillation light yield above the available energy range\n"
818 << "specified in G4MaterialPropertiesTable. A linear interpolation\n"
819 << "will be performed to compute the scintillation light yield using\n"
820 << "(L_max / E_max) as the photon yield per unit energy." << G4endl;
821 G4String cmt = "\nScintillation yield may be unphysical!\n";
822
823 if(fNumEnergyWarnings == 10)
824 {
825 ed << G4endl << "*** Scintillation energy warnings stopped.";
826 }
827 G4Exception("G4Scintillation::GetScintillationYieldByParticleType()",
828 "Scint03", JustWarning, ed, cmt);
829 }
830
831 // Units: [# scintillation photons]
832 ScintillationYield = yieldVector->GetMaxValue() /
833 yieldVector->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)
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
G4PhysicsFreeVector G4MaterialPropertyVector
@ fGeomBoundary
@ fStopButAlive
double G4double
Definition G4Types.hh:83
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 281 of file G4Scintillation.hh.

281{ return fStackingFlag; }

◆ GetTrackSecondariesFirst()

G4bool G4Scintillation::GetTrackSecondariesFirst ( ) const
inline

Definition at line 234 of file G4Scintillation.hh.

235{
236 return fTrackSecondariesFirst;
237}

◆ Initialise()

void G4Scintillation::Initialise ( )

Definition at line 162 of file G4Scintillation.cc.

163{
164 G4OpticalParameters* params = G4OpticalParameters::Instance();
171}
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 return (!aParticleType.IsShortLived());
153}

Referenced by LBE::ConstructOp().

◆ 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 258 of file G4Scintillation.cc.

264{
265 aParticleChange.Initialize(aTrack);
266 fNumPhotons = 0;
267
268 const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
269 const G4Material* aMaterial = aTrack.GetMaterial();
270
271 G4StepPoint* pPreStepPoint = aStep.GetPreStepPoint();
272 G4StepPoint* pPostStepPoint = aStep.GetPostStepPoint();
273
274 G4ThreeVector x0 = pPreStepPoint->GetPosition();
275 G4ThreeVector p0 = aStep.GetDeltaPosition().unit();
276 G4double t0 = pPreStepPoint->GetGlobalTime();
277
278 G4double TotalEnergyDeposit = aStep.GetTotalEnergyDeposit();
279 if (0.0 >= TotalEnergyDeposit) {
281 }
282
283 G4MaterialPropertiesTable* MPT = aMaterial->GetMaterialPropertiesTable();
284 if(!MPT)
285 return G4VRestDiscreteProcess::PostStepDoIt(aTrack, aStep);
286
287 G4int N_timeconstants = 1;
288
290 N_timeconstants = 3;
292 N_timeconstants = 2;
293 else if(!(MPT->GetProperty(kSCINTILLATIONCOMPONENT1)))
294 {
295 // no components were specified
296 return G4VRestDiscreteProcess::PostStepDoIt(aTrack, aStep);
297 }
298
299 G4double ResolutionScale = MPT->GetConstProperty(kRESOLUTIONSCALE);
300 G4double MeanNumberOfPhotons;
301
302 G4double yield1 = 0.;
303 G4double yield2 = 0.;
304 G4double yield3 = 0.;
305 G4double timeconstant1 = 0.;
306 G4double timeconstant2 = 0.;
307 G4double timeconstant3 = 0.;
308 G4double sum_yields = 0.;
309
310 if(fScintillationByParticleType)
311 {
312 MeanNumberOfPhotons = GetScintillationYieldByParticleType(
313 aTrack, aStep, yield1, yield2, yield3, timeconstant1, timeconstant2,
314 timeconstant3);
315 }
316 else
317 {
320 : 1.;
323 : 0.;
326 : 0.;
327 // The default linear scintillation process
328 // Units: [# scintillation photons / MeV]
329 MeanNumberOfPhotons = MPT->GetConstProperty(kSCINTILLATIONYIELD);
330 // Birk's correction via fEmSaturation and specifying scintillation by
331 // by particle type are physically mutually exclusive
332 if(fEmSaturation)
333 MeanNumberOfPhotons *=
334 (fEmSaturation->VisibleEnergyDepositionAtAStep(&aStep));
335 else
336 MeanNumberOfPhotons *= TotalEnergyDeposit;
337 }
338 sum_yields = yield1 + yield2 + yield3;
339
340 if(MeanNumberOfPhotons > 10.)
341 {
342 G4double sigma = ResolutionScale * std::sqrt(MeanNumberOfPhotons);
343 fNumPhotons = G4int(G4RandGauss::shoot(MeanNumberOfPhotons, sigma) + 0.5);
344 }
345 else
346 {
347 fNumPhotons = G4int(G4Poisson(MeanNumberOfPhotons));
348 }
349
350 if(fNumPhotons <= 0 || !fStackingFlag)
351 {
352 // return unchanged particle and no secondaries
353 aParticleChange.SetNumberOfSecondaries(0);
354 return G4VRestDiscreteProcess::PostStepDoIt(aTrack, aStep);
355 }
356
357 aParticleChange.SetNumberOfSecondaries(fNumPhotons);
358
359 if(fTrackSecondariesFirst)
360 {
361 if(aTrack.GetTrackStatus() == fAlive)
362 aParticleChange.ProposeTrackStatus(fSuspend);
363 }
364
365 std::size_t materialIndex = aMaterial->GetIndex();
366 auto it = fIndexMPT.find(materialIndex);
367
368 std::size_t indexMPT = 0;
369 if(it != fIndexMPT.end())
370 {
371 indexMPT = it->second;
372 }
373 else
374 {
376 ed << "G4MaterialPropertiesTable for " << aMaterial->GetName()
377 << " is not found!" << G4endl;
378 G4Exception("G4Scintillation::PostStepDoIt", "Scint04", FatalException, ed);
379 }
380
381 // Retrieve the Scintillation Integral for this material
382 // new G4PhysicsFreeVector allocated to hold CII's
383 G4int numPhot = fNumPhotons;
384 G4double scintTime = 0.;
385 G4double riseTime = 0.;
386 G4PhysicsFreeVector* scintIntegral = nullptr;
387 G4ScintillationType scintType = Slow;
388
389 G4bool isNeutral = (aParticle->GetDefinition()->GetPDGCharge() == 0);
390 G4double deltaVelocity = pPostStepPoint->GetVelocity() -
391 pPreStepPoint->GetVelocity();
392 auto touchableHandle = aStep.GetPreStepPoint()->GetTouchableHandle();
393
394 for(G4int scnt = 0; scnt < N_timeconstants; ++scnt)
395 {
396 // if there is 1 time constant it is #1, etc.
397 if(scnt == 0)
398 {
399 if(N_timeconstants == 1)
400 {
401 numPhot = fNumPhotons;
402 }
403 else
404 {
405 numPhot = yield1 / sum_yields * fNumPhotons;
406 }
407 if(fScintillationByParticleType)
408 {
409 scintTime = timeconstant1;
410 }
411 else
412 {
414 }
415 if(fFiniteRiseTime)
416 {
418 }
419 scintType = Fast;
420 scintIntegral = (G4PhysicsFreeVector*) ((*fIntegralTable1)(indexMPT));
421 }
422 else if(scnt == 1)
423 {
424 // to be consistent with old version (due to double->int conversion)
425 if(N_timeconstants == 2)
426 {
427 numPhot = fNumPhotons - numPhot;
428 }
429 else
430 {
431 numPhot = yield2 / sum_yields * fNumPhotons;
432 }
433 if(fScintillationByParticleType)
434 {
435 scintTime = timeconstant2;
436 }
437 else
438 {
440 }
441 if(fFiniteRiseTime)
442 {
444 }
445 scintType = Medium;
446 scintIntegral = (G4PhysicsFreeVector*) ((*fIntegralTable2)(indexMPT));
447 }
448 else if(scnt == 2)
449 {
450 numPhot = yield3 / sum_yields * fNumPhotons;
451 if(fScintillationByParticleType)
452 {
453 scintTime = timeconstant3;
454 }
455 else
456 {
458 }
459 if(fFiniteRiseTime)
460 {
462 }
463 scintType = Slow;
464 scintIntegral = (G4PhysicsFreeVector*) ((*fIntegralTable3)(indexMPT));
465 }
466
467 if(!scintIntegral)
468 continue;
469
470 for(G4int i = 0; i < numPhot; ++i)
471 {
472 // Determine photon energy
473 G4double sampledEnergy = scintIntegral->Value(G4UniformRand());
474
475 if(verboseLevel > 1)
476 {
477 G4cout << "sampledEnergy = " << sampledEnergy << G4endl;
478 }
479
480 // Generate random photon direction
481 G4double cost = 1. - 2. * G4UniformRand();
482 G4double sint = std::sqrt((1. - cost) * (1. + cost));
483 G4double phi = twopi * G4UniformRand();
484 G4double sinp = std::sin(phi);
485 G4double cosp = std::cos(phi);
486 G4ParticleMomentum photonMomentum(sint * cosp, sint * sinp, cost);
487
488 // Determine polarization of new photon
489 G4ThreeVector photonPolarization(cost * cosp, cost * sinp, -sint);
490 G4ThreeVector perp = photonMomentum.cross(photonPolarization);
491 phi = twopi * G4UniformRand();
492 sinp = std::sin(phi);
493 cosp = std::cos(phi);
494 photonPolarization = (cosp * photonPolarization + sinp * perp).unit();
495
496 // Generate a new photon:
497 auto scintPhoton = new G4DynamicParticle(opticalphoton, photonMomentum);
498 scintPhoton->SetPolarization(photonPolarization);
499 scintPhoton->SetKineticEnergy(sampledEnergy);
500
501 // Generate new G4Track object:
502 G4double rand = (isNeutral) ? 1.0 : G4UniformRand();
503
504 // emission time distribution
505 G4double delta = rand * aStep.GetStepLength();
506 G4double deltaTime =
507 delta / (pPreStepPoint->GetVelocity() + 0.5 * rand * deltaVelocity);
508 if(riseTime == 0.0)
509 {
510 deltaTime -= scintTime * std::log(G4UniformRand());
511 }
512 else
513 {
514 deltaTime += sample_time(riseTime, scintTime);
515 }
516
517 G4double secTime = t0 + deltaTime;
518 G4ThreeVector secPosition = x0 + rand * aStep.GetDeltaPosition();
519
520 G4Track* secTrack = new G4Track(scintPhoton, secTime, secPosition);
521 secTrack->SetTouchableHandle(touchableHandle);
522 secTrack->SetParentID(aTrack.GetTrackID());
523 secTrack->SetCreatorModelID(secID);
524 if(fScintillationTrackInfo)
525 secTrack->SetUserInformation(
526 new G4ScintillationTrackInformation(scintType));
527 aParticleChange.AddSecondary(secTrack);
528 }
529 }
530
531 if(verboseLevel > 1)
532 {
533 G4cout << "\n Exiting from G4Scintillation::DoIt -- NumberOfSecondaries = "
534 << aParticleChange.GetNumberOfSecondaries() << G4endl;
535 }
536
537 return G4VRestDiscreteProcess::PostStepDoIt(aTrack, aStep);
538}
G4ThreeVector G4ParticleMomentum
G4long G4Poisson(G4double mean)
Definition G4Poisson.hh:50
CLHEP::Hep3Vector G4ThreeVector
@ fSuspend
@ fAlive
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
#define G4UniformRand()
Definition Randomize.hh:52
Hep3Vector unit() const
Hep3Vector cross(const Hep3Vector &) const
std::size_t GetIndex() const
const G4String & GetName() 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)
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 156 of file G4Scintillation.cc.

157{
158 Initialise();
159}

◆ 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
140 G4OpticalParameters* params = G4OpticalParameters::Instance();
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 264 of file G4Scintillation.hh.

264{ fEmSaturation = nullptr; }

Referenced by SetScintillationByParticleType().

◆ SetFiniteRiseTime()

void G4Scintillation::SetFiniteRiseTime ( const G4bool state)

Definition at line 906 of file G4Scintillation.cc.

907{
908 fFiniteRiseTime = state;
910}
void SetScintFiniteRiseTime(G4bool)

Referenced by Initialise().

◆ SetScintillationByParticleType()

void G4Scintillation::SetScintillationByParticleType ( const G4bool scintType)

Definition at line 913 of file G4Scintillation.cc.

914{
915 if(fEmSaturation && scintType)
916 {
917 G4Exception("G4Scintillation::SetScintillationByParticleType", "Scint02",
919 "Redefinition: Birks Saturation is replaced by "
920 "ScintillationByParticleType!");
922 }
923 fScintillationByParticleType = scintType;
925 fScintillationByParticleType);
926}
void SetScintByParticleType(G4bool)

Referenced by Initialise().

◆ SetScintillationTrackInfo()

void G4Scintillation::SetScintillationTrackInfo ( const G4bool trackType)

Definition at line 929 of file G4Scintillation.cc.

930{
931 fScintillationTrackInfo = trackType;
932 G4OpticalParameters::Instance()->SetScintTrackInfo(fScintillationTrackInfo);
933}

Referenced by Initialise().

◆ SetStackPhotons()

void G4Scintillation::SetStackPhotons ( const G4bool stackingFlag)

Definition at line 936 of file G4Scintillation.cc.

937{
938 fStackingFlag = stackingFlag;
940}

Referenced by Initialise().

◆ SetTrackSecondariesFirst()

void G4Scintillation::SetTrackSecondariesFirst ( const G4bool state)

Definition at line 898 of file G4Scintillation.cc.

899{
900 fTrackSecondariesFirst = state;
902 fTrackSecondariesFirst);
903}
void SetScintTrackSecondariesFirst(G4bool)

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

◆ SetVerboseLevel()

void G4Scintillation::SetVerboseLevel ( G4int verbose)

Definition at line 943 of file G4Scintillation.cc.

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


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