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

#include <G4QuasiScintillation.hh>

Inheritance diagram for G4QuasiScintillation:

Public Member Functions

 G4QuasiScintillation (const G4String &procName="QausiScintillation", G4ProcessType type=fElectromagnetic)
 ~G4QuasiScintillation ()
 G4QuasiScintillation (const G4QuasiScintillation &right)=delete
G4QuasiScintillationoperator= (const G4QuasiScintillation &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
void SetOffloadPhotons (const G4bool)
G4bool GetOffloadPhotons () 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 62 of file G4QuasiScintillation.hh.

Constructor & Destructor Documentation

◆ G4QuasiScintillation() [1/2]

G4QuasiScintillation::G4QuasiScintillation ( const G4String & procName = "QausiScintillation",
G4ProcessType type = fElectromagnetic )
explicit

Definition at line 51 of file G4QuasiScintillation.cc.

53 : G4VRestDiscreteProcess(processName, type)
54 , fIntegralTable1(nullptr)
55 , fIntegralTable2(nullptr)
56 , fIntegralTable3(nullptr)
57 , fEmSaturation(nullptr)
58 , fNumPhotons(0)
59{
60 secID = G4PhysicsModelCatalog::GetModelID("model_QuasiScintillation");
62
63#ifdef G4DEBUG_SCINTILLATION
64 ScintTrackEDep = 0.;
65 ScintTrackYield = 0.;
66#endif
67
68 if(verboseLevel > 0)
69 {
70 G4cout << GetProcessName() << " is created " << G4endl;
71 }
72 Initialise();
73}
@ 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 G4QuasiScintillation(), and operator=().

◆ ~G4QuasiScintillation()

G4QuasiScintillation::~G4QuasiScintillation ( )

Definition at line 76 of file G4QuasiScintillation.cc.

77{
78 if(fIntegralTable1 != nullptr)
79 {
80 fIntegralTable1->clearAndDestroy();
81 delete fIntegralTable1;
82 }
83 if(fIntegralTable2 != nullptr)
84 {
85 fIntegralTable2->clearAndDestroy();
86 delete fIntegralTable2;
87 }
88 if(fIntegralTable3 != nullptr)
89 {
90 fIntegralTable3->clearAndDestroy();
91 delete fIntegralTable3;
92 }
93}

◆ G4QuasiScintillation() [2/2]

G4QuasiScintillation::G4QuasiScintillation ( const G4QuasiScintillation & right)
delete

Member Function Documentation

◆ AddSaturation()

void G4QuasiScintillation::AddSaturation ( G4EmSaturation * sat)
inline

Definition at line 255 of file G4QuasiScintillation.hh.

256{
257 fEmSaturation = sat;
258}

◆ AtRestDoIt()

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

Reimplemented from G4VRestDiscreteProcess.

Definition at line 219 of file G4QuasiScintillation.cc.

223{
224 return G4QuasiScintillation::PostStepDoIt(aTrack, aStep);
225}
G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep) override

◆ BuildPhysicsTable()

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

Reimplemented from G4VProcess.

Definition at line 143 of file G4QuasiScintillation.cc.

144{
145 const G4MaterialTable* materialTable = G4Material::GetMaterialTable();
146 std::size_t numOfMaterials = G4Material::GetNumberOfMaterials();
147
148 // Find the number of materials that have non-empty material property tables
149 std::size_t numOfMaterialsWithMPT = 0;
150 for(std::size_t i = 0; i < numOfMaterials; ++i)
151 {
152 if(((*materialTable)[i])->GetMaterialPropertiesTable())
153 {
154 ++numOfMaterialsWithMPT;
155 }
156 }
157
158 // create new physics table
159 fIntegralTable1 = new G4PhysicsTable(numOfMaterialsWithMPT);
160 fIntegralTable2 = new G4PhysicsTable(numOfMaterialsWithMPT);
161 fIntegralTable3 = new G4PhysicsTable(numOfMaterialsWithMPT);
162
163 std::size_t indexMPT = 0;
164 for(std::size_t i = 0; i < numOfMaterials; ++i)
165 {
166 // Retrieve vector of scintillation wavelength intensity for
167 // the material from the material's optical properties table.
168 G4MaterialPropertiesTable* MPT =
169 ((*materialTable)[i])->GetMaterialPropertiesTable();
170
171 if(MPT)
172 {
173 auto vector1 = new G4PhysicsFreeVector();
174 auto vector2 = new G4PhysicsFreeVector();
175 auto vector3 = new G4PhysicsFreeVector();
176
177 BuildInverseCdfTable(MPT->GetProperty(kSCINTILLATIONCOMPONENT1), vector1);
178 BuildInverseCdfTable(MPT->GetProperty(kSCINTILLATIONCOMPONENT2), vector2);
179 BuildInverseCdfTable(MPT->GetProperty(kSCINTILLATIONCOMPONENT3), vector3);
180
181 fIntegralTable1->insertAt(indexMPT, vector1);
182 fIntegralTable2->insertAt(indexMPT, vector2);
183 fIntegralTable3->insertAt(indexMPT, vector3);
184
185 fIndexMPT.insert(std::make_pair(i, indexMPT));
186 ++indexMPT;
187 }
188 }
189}
@ kSCINTILLATIONCOMPONENT1
@ kSCINTILLATIONCOMPONENT2
@ kSCINTILLATIONCOMPONENT3
std::vector< G4Material * > G4MaterialTable
G4MaterialPropertyVector * GetProperty(const char *key) const
static std::size_t GetNumberOfMaterials()
static G4MaterialTable * GetMaterialTable()

◆ DumpInfo()

void G4QuasiScintillation::DumpInfo ( ) const
inlineoverridevirtual

Reimplemented from G4VProcess.

Definition at line 81 of file G4QuasiScintillation.hh.

void ProcessDescription(std::ostream &) const override

◆ DumpPhysicsTable()

void G4QuasiScintillation::DumpPhysicsTable ( ) const

Definition at line 865 of file G4QuasiScintillation.cc.

866{
867 if(fIntegralTable1)
868 {
869 for(std::size_t i = 0; i < fIntegralTable1->entries(); ++i)
870 {
871 ((G4PhysicsFreeVector*) (*fIntegralTable1)[i])->DumpValues();
872 }
873 }
874 if(fIntegralTable2)
875 {
876 for(std::size_t i = 0; i < fIntegralTable2->entries(); ++i)
877 {
878 ((G4PhysicsFreeVector*) (*fIntegralTable2)[i])->DumpValues();
879 }
880 }
881 if(fIntegralTable3)
882 {
883 for(std::size_t i = 0; i < fIntegralTable3->entries(); ++i)
884 {
885 ((G4PhysicsFreeVector*) (*fIntegralTable3)[i])->DumpValues();
886 }
887 }
888}

◆ GetFiniteRiseTime()

G4bool G4QuasiScintillation::GetFiniteRiseTime ( ) const
inline

Definition at line 235 of file G4QuasiScintillation.hh.

236{
237 return fFiniteRiseTime;
238}

◆ GetIntegralTable1()

G4PhysicsTable * G4QuasiScintillation::GetIntegralTable1 ( ) const
inline

Definition at line 240 of file G4QuasiScintillation.hh.

241{
242 return fIntegralTable1;
243}

◆ GetIntegralTable2()

G4PhysicsTable * G4QuasiScintillation::GetIntegralTable2 ( ) const
inline

Definition at line 245 of file G4QuasiScintillation.hh.

246{
247 return fIntegralTable2;
248}

◆ GetIntegralTable3()

G4PhysicsTable * G4QuasiScintillation::GetIntegralTable3 ( ) const
inline

Definition at line 250 of file G4QuasiScintillation.hh.

251{
252 return fIntegralTable3;
253}

◆ GetMeanFreePath()

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

Implements G4VRestDiscreteProcess.

Definition at line 532 of file G4QuasiScintillation.cc.

534{
536 return DBL_MAX;
537}
G4double condition(const G4ErrorSymMatrix &m)
@ StronglyForced
#define DBL_MAX
Definition templates.hh:62

◆ GetMeanLifeTime()

G4double G4QuasiScintillation::GetMeanLifeTime ( const G4Track & aTrack,
G4ForceCondition * condition )
overridevirtual

Implements G4VRestDiscreteProcess.

Definition at line 540 of file G4QuasiScintillation.cc.

542{
543 *condition = Forced;
544 return DBL_MAX;
545}
@ Forced

◆ GetNumPhotons()

G4int G4QuasiScintillation::GetNumPhotons ( ) const
inline

Definition at line 290 of file G4QuasiScintillation.hh.

291{
292 return fNumPhotons;
293}

◆ GetOffloadPhotons()

G4bool G4QuasiScintillation::GetOffloadPhotons ( ) const
inline

Definition at line 285 of file G4QuasiScintillation.hh.

286{
287 return fOffloadingFlag;
288}

◆ GetSaturation()

G4EmSaturation * G4QuasiScintillation::GetSaturation ( ) const
inline

Definition at line 265 of file G4QuasiScintillation.hh.

266{
267 return fEmSaturation;
268}

◆ GetScintillationByParticleType()

G4bool G4QuasiScintillation::GetScintillationByParticleType ( ) const
inline

Definition at line 270 of file G4QuasiScintillation.hh.

271{
272 return fScintillationByParticleType;
273}

◆ GetScintillationTrackInfo()

G4bool G4QuasiScintillation::GetScintillationTrackInfo ( ) const
inline

Definition at line 275 of file G4QuasiScintillation.hh.

276{
277 return fScintillationTrackInfo;
278}

◆ GetScintillationYieldByParticleType()

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

Definition at line 565 of file G4QuasiScintillation.cc.

569{
570 // new in 10.7, allow multiple time constants with ScintByParticleType
571 // Get the G4MaterialPropertyVector containing the scintillation
572 // yield as a function of the energy deposited and particle type
573 // In 11.2, allow different time constants for different particles
574
575 G4ParticleDefinition* pDef = aTrack.GetDynamicParticle()->GetDefinition();
576 G4MaterialPropertyVector* yieldVector = nullptr;
577 G4MaterialPropertiesTable* MPT =
579
580 // Protons
581 if(pDef == G4Proton::ProtonDefinition())
582 {
583 yieldVector = MPT->GetProperty(kPROTONSCINTILLATIONYIELD);
586 : 1.;
589 : 0.;
592 : 0.;
596 if(yield2 > 0.)
597 {
601 }
602 if(yield3 > 0.)
603 {
607 }
608 }
609
610 // Deuterons
611 else if(pDef == G4Deuteron::DeuteronDefinition())
612 {
613 yieldVector = MPT->GetProperty(kDEUTERONSCINTILLATIONYIELD);
616 : 1.;
619 : 0.;
622 : 0.;
626 if(yield2 > 0.)
627 {
631 }
632 if(yield3 > 0.)
633 {
637 }
638 }
639
640 // Tritons
641 else if(pDef == G4Triton::TritonDefinition())
642 {
643 yieldVector = MPT->GetProperty(kTRITONSCINTILLATIONYIELD);
646 : 1.;
649 : 0.;
652 : 0.;
656 if(yield2 > 0.)
657 {
661 }
662 if(yield3 > 0.)
663 {
667 }
668 }
669
670 // Alphas
671 else if(pDef == G4Alpha::AlphaDefinition())
672 {
673 yieldVector = MPT->GetProperty(kALPHASCINTILLATIONYIELD);
676 : 1.;
679 : 0.;
682 : 0.;
686 if(yield2 > 0.)
687 {
691 }
692 if(yield3 > 0.)
693 {
697 }
698 }
699
700 // Ions (particles derived from G4VIon and G4Ions) and recoil ions
701 // below the production cut from neutrons after hElastic
702 else if(pDef->GetParticleType() == "nucleus" ||
704 {
705 yieldVector = MPT->GetProperty(kIONSCINTILLATIONYIELD);
708 : 1.;
711 : 0.;
714 : 0.;
718 if(yield2 > 0.)
719 {
723 }
724 if(yield3 > 0.)
725 {
729 }
730 }
731
732 // Electrons (must also account for shell-binding energy
733 // attributed to gamma from standard photoelectric effect)
734 // and, default for particles not enumerated/listed above
735 else
736 {
737 yieldVector = MPT->GetProperty(kELECTRONSCINTILLATIONYIELD);
740 : 1.;
743 : 0.;
746 : 0.;
750 if(yield2 > 0.)
751 {
755 }
756 if(yield3 > 0.)
757 {
761 }
762 }
763
764 // Throw an exception if no scintillation yield vector is found
765 if(yieldVector == nullptr)
766 {
768 ed << "\nG4QuasiScintillation::PostStepDoIt(): "
769 << "Request for scintillation yield for energy deposit and particle\n"
770 << "type without correct entry in MaterialPropertiesTable. A material\n"
771 << "property (vector) with name like PARTICLESCINTILLATIONYIELD is\n"
772 << "needed (hint: PARTICLE might not be the primary particle."
773 << G4endl;
774 G4String comments = "Missing MaterialPropertiesTable entry - No correct "
775 "entry in MaterialPropertiesTable";
776 G4Exception("G4QuasiScintillation::PostStepDoIt", "Scint01", FatalException, ed,
777 comments);
778 return 0.; // NOLINT: required to help Coverity recognise this as exit point
779 }
780
781 ///////////////////////////////////////
782 // Calculate the scintillation light //
783 ///////////////////////////////////////
784 // To account for potential nonlinearity and scintillation photon
785 // density along the track, light (L) is produced according to:
786 // L_currentStep = L(PreStepKE) - L(PreStepKE - EDep)
787
788 G4double ScintillationYield = 0.;
789 G4double StepEnergyDeposit = aStep.GetTotalEnergyDeposit();
790 G4double PreStepKineticEnergy = aStep.GetPreStepPoint()->GetKineticEnergy();
791
792 if(PreStepKineticEnergy <= yieldVector->GetMaxEnergy())
793 {
794 // G4double Yield1 = yieldVector->Value(PreStepKineticEnergy);
795 // G4double Yield2 = yieldVector->Value(PreStepKineticEnergy -
796 // StepEnergyDeposit); ScintillationYield = Yield1 - Yield2;
797 ScintillationYield =
798 yieldVector->Value(PreStepKineticEnergy) -
799 yieldVector->Value(PreStepKineticEnergy - StepEnergyDeposit);
800 }
801 else
802 {
803 ++fNumEnergyWarnings;
804 if(verboseLevel > 0 && fNumEnergyWarnings <= 10)
805 {
807 ed << "\nG4QuasiScintillation::GetScintillationYieldByParticleType(): "
808 << "Request\n"
809 << "for scintillation light yield above the available energy range\n"
810 << "specified in G4MaterialPropertiesTable. A linear interpolation\n"
811 << "will be performed to compute the scintillation light yield using\n"
812 << "(L_max / E_max) as the photon yield per unit energy." << G4endl;
813 G4String cmt = "\nScintillation yield may be unphysical!\n";
814
815 if(fNumEnergyWarnings == 10)
816 {
817 ed << G4endl << "*** Scintillation energy warnings stopped.";
818 }
819 G4Exception("G4QuasiScintillation::GetScintillationYieldByParticleType()",
820 "Scint03", JustWarning, ed, cmt);
821 }
822
823 // Units: [# scintillation photons]
824 ScintillationYield = yieldVector->GetMaxValue() /
825 yieldVector->GetMaxEnergy() * StepEnergyDeposit;
826 }
827
828#ifdef G4DEBUG_SCINTILLATION
829 // Increment track aggregators
830 ScintTrackYield += ScintillationYield;
831 ScintTrackEDep += StepEnergyDeposit;
832
833 G4cout << "\n-- "
834 << "G4QuasiScintillation::GetScintillationYieldByParticleType() --\n"
835 << "--\n"
836 << "-- Name = "
837 << aTrack.GetParticleDefinition()->GetParticleName() << "\n"
838 << "-- TrackID = " << aTrack.GetTrackID() << "\n"
839 << "-- ParentID = " << aTrack.GetParentID() << "\n"
840 << "-- Current KE = " << aTrack.GetKineticEnergy() / MeV
841 << " MeV\n"
842 << "-- Step EDep = " << aStep.GetTotalEnergyDeposit() / MeV
843 << " MeV\n"
844 << "-- Track EDep = " << ScintTrackEDep / MeV << " MeV\n"
845 << "-- Vertex KE = " << aTrack.GetVertexKineticEnergy() / MeV
846 << " MeV\n"
847 << "-- Step yield = " << ScintillationYield << " photons\n"
848 << "-- Track yield = " << ScintTrackYield << " photons\n"
849 << G4endl;
850
851 // The track has terminated within or has left the scintillator volume
852 if((aTrack.GetTrackStatus() == fStopButAlive) or
854 {
855 // Reset aggregators for the next track
856 ScintTrackEDep = 0.;
857 ScintTrackYield = 0.;
858 }
859#endif
860
861 return ScintillationYield;
862}
@ 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 G4QuasiScintillation::GetStackPhotons ( ) const
inline

Definition at line 280 of file G4QuasiScintillation.hh.

281{
282 return fStackingFlag;
283}

◆ GetTrackSecondariesFirst()

G4bool G4QuasiScintillation::GetTrackSecondariesFirst ( ) const
inline

Definition at line 230 of file G4QuasiScintillation.hh.

231{
232 return fTrackSecondariesFirst;
233}

◆ Initialise()

void G4QuasiScintillation::Initialise ( )

Definition at line 130 of file G4QuasiScintillation.cc.

131{
132 G4OpticalParameters* params = G4OpticalParameters::Instance();
140}
G4bool GetScintOffloadPhotons() const
G4int GetScintVerboseLevel() const
G4bool GetScintStackPhotons() const
static G4OpticalParameters * Instance()
G4bool GetScintByParticleType() const
G4bool GetScintFiniteRiseTime() const
G4bool GetScintTrackInfo() const
G4bool GetScintTrackSecondariesFirst() const
void SetStackPhotons(const G4bool)
void SetScintillationTrackInfo(const G4bool trackType)
void SetTrackSecondariesFirst(const G4bool state)
void SetScintillationByParticleType(const G4bool)
void SetOffloadPhotons(const G4bool)
void SetFiniteRiseTime(const G4bool state)

Referenced by G4QuasiScintillation(), and PreparePhysicsTable().

◆ IsApplicable()

G4bool G4QuasiScintillation::IsApplicable ( const G4ParticleDefinition & aParticleType)
overridevirtual

Reimplemented from G4VProcess.

Definition at line 114 of file G4QuasiScintillation.cc.

115{
116 if(aParticleType.GetParticleName() == "opticalphoton")
117 return false;
118 if(aParticleType.IsShortLived())
119 return false;
120 return true;
121}

◆ operator=()

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

◆ PostStepDoIt()

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

Reimplemented from G4VRestDiscreteProcess.

Definition at line 228 of file G4QuasiScintillation.cc.

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

Reimplemented from G4VProcess.

Definition at line 124 of file G4QuasiScintillation.cc.

125{
126 Initialise();
127}

◆ ProcessDescription()

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

Reimplemented from G4VProcess.

Definition at line 96 of file G4QuasiScintillation.cc.

97{
98 out << "Scintillation simulates production of optical photons produced\n"
99 "by a high energy particle traversing matter.\n"
100 "Various material properties need to be defined.\n";
102
103 G4OpticalParameters* params = G4OpticalParameters::Instance();
104 out << "Track secondaries first: " << params->GetScintTrackSecondariesFirst();
105 out << "Finite rise time: " << params->GetScintFiniteRiseTime();
106 out << "Scintillation by particle type: " << params->GetScintByParticleType();
107 out << "Save track information: " << params->GetScintTrackInfo();
108 out << "Stack photons: " << params->GetScintStackPhotons();
109 out << "Verbose level: " << params->GetScintVerboseLevel();
110}
virtual void DumpInfo() const

Referenced by DumpInfo().

◆ RemoveSaturation()

void G4QuasiScintillation::RemoveSaturation ( )
inline

Definition at line 260 of file G4QuasiScintillation.hh.

261{
262 fEmSaturation = nullptr;
263}

Referenced by SetScintillationByParticleType().

◆ SetFiniteRiseTime()

void G4QuasiScintillation::SetFiniteRiseTime ( const G4bool state)

Definition at line 899 of file G4QuasiScintillation.cc.

900{
901 fFiniteRiseTime = state;
903}
void SetScintFiniteRiseTime(G4bool)

Referenced by Initialise().

◆ SetOffloadPhotons()

void G4QuasiScintillation::SetOffloadPhotons ( const G4bool offloadingFlag)

Definition at line 936 of file G4QuasiScintillation.cc.

937{
938 fOffloadingFlag = offloadingFlag;
940}
void SetScintOffloadPhotons(G4bool)

Referenced by Initialise().

◆ SetScintillationByParticleType()

void G4QuasiScintillation::SetScintillationByParticleType ( const G4bool scintType)

Definition at line 906 of file G4QuasiScintillation.cc.

907{
908 if(fEmSaturation && scintType)
909 {
910 G4Exception("G4QuasiScintillation::SetScintillationByParticleType",
911 "Scint02", JustWarning,
912 "Redefinition: Birks Saturation is replaced by "
913 "ScintillationByParticleType!");
915 }
916 fScintillationByParticleType = scintType;
918 fScintillationByParticleType);
919}
void SetScintByParticleType(G4bool)

Referenced by Initialise().

◆ SetScintillationTrackInfo()

void G4QuasiScintillation::SetScintillationTrackInfo ( const G4bool trackType)

Definition at line 922 of file G4QuasiScintillation.cc.

923{
924 fScintillationTrackInfo = trackType;
925 G4OpticalParameters::Instance()->SetScintTrackInfo(fScintillationTrackInfo);
926}

Referenced by Initialise().

◆ SetStackPhotons()

void G4QuasiScintillation::SetStackPhotons ( const G4bool stackingFlag)

Definition at line 929 of file G4QuasiScintillation.cc.

930{
931 fStackingFlag = stackingFlag;
933}

Referenced by Initialise().

◆ SetTrackSecondariesFirst()

void G4QuasiScintillation::SetTrackSecondariesFirst ( const G4bool state)

Definition at line 891 of file G4QuasiScintillation.cc.

892{
893 fTrackSecondariesFirst = state;
895 fTrackSecondariesFirst);
896}
void SetScintTrackSecondariesFirst(G4bool)

Referenced by Initialise().

◆ SetVerboseLevel()

void G4QuasiScintillation::SetVerboseLevel ( G4int verbose)

Definition at line 943 of file G4QuasiScintillation.cc.

Referenced by Initialise().


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