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

#include <G4OpWLS.hh>

Inheritance diagram for G4OpWLS:

Public Member Functions

 G4OpWLS (const G4String &processName="OpWLS", G4ProcessType type=fOptical)
virtual ~G4OpWLS ()
virtual G4bool IsApplicable (const G4ParticleDefinition &aParticleType) override
virtual void BuildPhysicsTable (const G4ParticleDefinition &aParticleType) override
virtual G4double GetMeanFreePath (const G4Track &aTrack, G4double, G4ForceCondition *) override
virtual G4VParticleChangePostStepDoIt (const G4Track &aTrack, const G4Step &aStep) override
virtual G4PhysicsTableGetIntegralTable () const
virtual void DumpPhysicsTable () const
virtual void UseTimeProfile (const G4String name)
virtual void PreparePhysicsTable (const G4ParticleDefinition &) override
virtual void Initialise ()
void SetVerboseLevel (G4int)
Public Member Functions inherited from G4VDiscreteProcess
 G4VDiscreteProcess (const G4String &aName, G4ProcessType aType=fNotDefined)
 G4VDiscreteProcess (G4VDiscreteProcess &)
virtual ~G4VDiscreteProcess ()
G4VDiscreteProcessoperator= (const G4VDiscreteProcess &)=delete
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
virtual G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &)
virtual G4VParticleChangeAlongStepDoIt (const G4Track &, const G4Step &)
virtual G4double GetCrossSection (const G4double, const G4MaterialCutsCouple *)
virtual G4double MinPrimaryEnergy (const G4ParticleDefinition *, const G4Material *)
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
virtual void DumpInfo () const
virtual void ProcessDescription (std::ostream &outfile) const
void SetVerboseLevel (G4int value)
G4int GetVerboseLevel () const
virtual void SetMasterProcess (G4VProcess *masterP)
const G4VProcessGetMasterProcess () const
virtual void BuildWorkerPhysicsTable (const G4ParticleDefinition &part)
virtual void PrepareWorkerPhysicsTable (const G4ParticleDefinition &)

Protected Attributes

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

Additional Inherited Members

Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double prevStepSize)
void ClearNumberOfInteractionLengthLeft ()

Detailed Description

Definition at line 51 of file G4OpWLS.hh.

Constructor & Destructor Documentation

◆ G4OpWLS()

G4OpWLS::G4OpWLS ( const G4String & processName = "OpWLS",
G4ProcessType type = fOptical )
explicit

Definition at line 54 of file G4OpWLS.cc.

55 : G4VDiscreteProcess(processName, type)
56{
58 Initialise();
60 theIntegralTable = nullptr;
61
62 if(verboseLevel > 0)
63 G4cout << GetProcessName() << " is created " << G4endl;
64}
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
virtual void Initialise()
Definition G4OpWLS.cc:81
G4PhysicsTable * theIntegralTable
Definition G4OpWLS.hh:91
G4VWLSTimeGeneratorProfile * WLSTimeGeneratorProfile
Definition G4OpWLS.hh:90
G4VDiscreteProcess(const G4String &aName, G4ProcessType aType=fNotDefined)
G4int verboseLevel
void SetProcessSubType(G4int)
const G4String & GetProcessName() const

◆ ~G4OpWLS()

G4OpWLS::~G4OpWLS ( )
virtual

Definition at line 67 of file G4OpWLS.cc.

68{
70 {
71 theIntegralTable->clearAndDestroy();
72 delete theIntegralTable;
73 }
75}

Member Function Documentation

◆ BuildPhysicsTable()

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

Reimplemented from G4VProcess.

Definition at line 216 of file G4OpWLS.cc.

217{
219 {
220 theIntegralTable->clearAndDestroy();
221 delete theIntegralTable;
222 theIntegralTable = nullptr;
223 }
224
225 const G4MaterialTable* materialTable = G4Material::GetMaterialTable();
226 std::size_t numOfMaterials = G4Material::GetNumberOfMaterials();
227 theIntegralTable = new G4PhysicsTable(numOfMaterials);
228
229 // loop for materials
230 for(std::size_t i = 0; i < numOfMaterials; ++i)
231 {
232 auto physVector = new G4PhysicsFreeVector();
233
234 // Retrieve vector of WLS wavelength intensity for
235 // the material from the material's optical properties table.
236 G4MaterialPropertiesTable* MPT =
237 (*materialTable)[i]->GetMaterialPropertiesTable();
238 if(MPT)
239 {
241 if(wlsVector)
242 {
243 // Retrieve the first intensity point in vector
244 // of (photon energy, intensity) pairs
245 G4double currentIN = (*wlsVector)[0];
246 if(currentIN >= 0.0)
247 {
248 // Create first (photon energy)
249 G4double currentPM = wlsVector->Energy(0);
250 G4double currentCII = 0.0;
251 physVector->InsertValues(currentPM, currentCII);
252
253 // Set previous values to current ones prior to loop
254 G4double prevPM = currentPM;
255 G4double prevCII = currentCII;
256 G4double prevIN = currentIN;
257
258 // loop over all (photon energy, intensity)
259 // pairs stored for this material
260 for(std::size_t j = 1; j < wlsVector->GetVectorLength(); ++j)
261 {
262 currentPM = wlsVector->Energy(j);
263 currentIN = (*wlsVector)[j];
264 currentCII =
265 prevCII + 0.5 * (currentPM - prevPM) * (prevIN + currentIN);
266
267 physVector->InsertValues(currentPM, currentCII);
268
269 prevPM = currentPM;
270 prevCII = currentCII;
271 prevIN = currentIN;
272 }
273 }
274 }
275 }
276 theIntegralTable->insertAt(i, physVector);
277 }
278}
G4PhysicsFreeVector G4MaterialPropertyVector
std::vector< G4Material * > G4MaterialTable
double G4double
Definition G4Types.hh:83
G4MaterialPropertyVector * GetProperty(const char *key) const
static std::size_t GetNumberOfMaterials()
static G4MaterialTable * GetMaterialTable()
G4double Energy(const std::size_t index) const
std::size_t GetVectorLength() const

◆ DumpPhysicsTable()

void G4OpWLS::DumpPhysicsTable ( ) const
inlinevirtual

Definition at line 114 of file G4OpWLS.hh.

115{
116 std::size_t PhysicsTableSize = theIntegralTable->entries();
117 G4PhysicsFreeVector* v;
118
119 for(std::size_t i = 0; i < PhysicsTableSize; ++i)
120 {
121 v = (G4PhysicsFreeVector*) (*theIntegralTable)[i];
122 v->DumpValues();
123 }
124}
void DumpValues(G4double unitE=1.0, G4double unitV=1.0) const

◆ GetIntegralTable()

G4PhysicsTable * G4OpWLS::GetIntegralTable ( ) const
inlinevirtual

Definition at line 109 of file G4OpWLS.hh.

110{
111 return theIntegralTable;
112}

◆ GetMeanFreePath()

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

Implements G4VDiscreteProcess.

Definition at line 281 of file G4OpWLS.cc.

283{
284 G4double thePhotonEnergy = aTrack.GetDynamicParticle()->GetTotalEnergy();
285 G4double attLength = DBL_MAX;
286 G4MaterialPropertiesTable* MPT =
288
289 if(MPT)
290 {
292 if(attVector)
293 {
294 attLength = attVector->Value(thePhotonEnergy, idx_wls);
295 }
296 }
297 return attLength;
298}
G4double GetTotalEnergy() const
G4MaterialPropertiesTable * GetMaterialPropertiesTable() const
G4double Value(const G4double energy, std::size_t &lastidx) const
G4Material * GetMaterial() const
const G4DynamicParticle * GetDynamicParticle() const
#define DBL_MAX
Definition templates.hh:62

◆ Initialise()

void G4OpWLS::Initialise ( )
virtual

Definition at line 81 of file G4OpWLS.cc.

82{
83 G4OpticalParameters* params = G4OpticalParameters::Instance();
86}
void SetVerboseLevel(G4int)
Definition G4OpWLS.cc:326
virtual void UseTimeProfile(const G4String name)
Definition G4OpWLS.cc:301
static G4OpticalParameters * Instance()
G4String GetWLSTimeProfile() const

Referenced by G4OpWLS(), and PreparePhysicsTable().

◆ IsApplicable()

G4bool G4OpWLS::IsApplicable ( const G4ParticleDefinition & aParticleType)
inlineoverridevirtual

Reimplemented from G4VProcess.

Definition at line 104 of file G4OpWLS.hh.

105{
106 return (&aParticleType == G4OpticalPhoton::OpticalPhoton());
107}
static G4OpticalPhoton * OpticalPhoton()

◆ PostStepDoIt()

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

Reimplemented from G4VDiscreteProcess.

Definition at line 89 of file G4OpWLS.cc.

91{
92 aParticleChange.Initialize(aTrack);
93 aParticleChange.ProposeTrackStatus(fStopAndKill);
94
95 if(verboseLevel > 1)
96 {
97 G4cout << "\n** G4OpWLS: Photon absorbed! **" << G4endl;
98 }
99
100 G4MaterialPropertiesTable* MPT =
102 if(!MPT)
103 {
104 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
105 }
106
107 G4PhysicsFreeVector* WLSIntegral = nullptr;
108 if(!MPT->GetProperty(kWLSCOMPONENT))
109 {
110 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
111 }
112 else
113 {
114 // Retrieve the WLS Integral for this material
115 WLSIntegral = (G4PhysicsFreeVector*) ((*theIntegralTable)
116 (aTrack.GetMaterial()->GetIndex()));
117 }
118
119 G4double primaryEnergy = aTrack.GetDynamicParticle()->GetKineticEnergy();
120 // No WLS photons are produced if the primary photon's energy is below
121 // the lower bound of the WLS integral range
122 if(primaryEnergy < WLSIntegral->GetMinValue())
123 {
124 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
125 }
126
127 G4int NumPhotons = 1;
129 {
130 G4double MeanNumberOfPhotons = MPT->GetConstProperty(kWLSMEANNUMBERPHOTONS);
131 NumPhotons = G4int(G4Poisson(MeanNumberOfPhotons));
132 if(NumPhotons <= 0)
133 {
134 // return unchanged particle and no secondaries
135 aParticleChange.SetNumberOfSecondaries(0);
136 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
137 }
138 }
139
141 G4StepPoint* pPostStepPoint = aStep.GetPostStepPoint();
142 std::vector<G4Track*> proposedSecondaries;
143
144 // Max WLS Integral - force sampling of a WLS photon with energy below
145 // the primary photon's energy
146 G4double CIImax = (primaryEnergy > WLSIntegral->GetMaxEnergy()) ?
147 WLSIntegral->GetMaxValue() : WLSIntegral->Value(primaryEnergy);
148
149 for(G4int i = 0; i < NumPhotons; ++i)
150 {
151 // Determine photon energy
152 G4double sampledEnergy = WLSIntegral->GetEnergy(G4UniformRand() * CIImax);
153
154 // Make sure the energy of the secondary is less than that of the primary
155 if(sampledEnergy > primaryEnergy)
156 {
158 ed << "Sampled photon energy " << sampledEnergy << " is greater than "
159 << "the primary photon energy " << primaryEnergy << G4endl;
160 G4Exception("G4OpWLS::PostStepDoIt", "WSL01", FatalException, ed);
161 }
162 else if(verboseLevel > 1)
163 {
164 G4cout << "G4OpWLS: Created photon with energy: " << sampledEnergy
165 << G4endl;
166 }
167
168 // Generate random photon direction
169 G4double cost = 1. - 2. * G4UniformRand();
170 G4double sint = std::sqrt((1. - cost) * (1. + cost));
171 G4double phi = twopi * G4UniformRand();
172 G4double sinp = std::sin(phi);
173 G4double cosp = std::cos(phi);
174 G4ParticleMomentum photonMomentum(sint * cosp, sint * sinp, cost);
175
176 G4ThreeVector photonPolarization(cost * cosp, cost * sinp, -sint);
177 G4ThreeVector perp = photonMomentum.cross(photonPolarization);
178
179 phi = twopi * G4UniformRand();
180 sinp = std::sin(phi);
181 cosp = std::cos(phi);
182 photonPolarization = (cosp * photonPolarization + sinp * perp).unit();
183
184 // Generate a new photon:
185 auto sec_dp =
186 new G4DynamicParticle(G4OpticalPhoton::OpticalPhoton(), photonMomentum);
187 sec_dp->SetPolarization(photonPolarization);
188 sec_dp->SetKineticEnergy(sampledEnergy);
189
190 G4double secTime = pPostStepPoint->GetGlobalTime() +
191 WLSTimeGeneratorProfile->GenerateTime(WLSTime);
192 G4ThreeVector secPos = pPostStepPoint->GetPosition();
193 G4Track* secTrack = new G4Track(sec_dp, secTime, secPos);
194
195 secTrack->SetTouchableHandle(aTrack.GetTouchableHandle());
196 secTrack->SetParentID(aTrack.GetTrackID());
197
198 proposedSecondaries.push_back(secTrack);
199 }
200
201 aParticleChange.SetNumberOfSecondaries((G4int)proposedSecondaries.size());
202 for(auto sec : proposedSecondaries)
203 {
204 aParticleChange.AddSecondary(sec);
205 }
206 if(verboseLevel > 1)
207 {
208 G4cout << "\n Exiting from G4OpWLS::DoIt -- NumberOfSecondaries = "
209 << aParticleChange.GetNumberOfSecondaries() << G4endl;
210 }
211
212 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
213}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
G4ThreeVector G4ParticleMomentum
G4long G4Poisson(G4double mean)
Definition G4Poisson.hh:50
CLHEP::Hep3Vector G4ThreeVector
@ fStopAndKill
int G4int
Definition G4Types.hh:85
#define G4UniformRand()
Definition Randomize.hh:52
Hep3Vector cross(const Hep3Vector &) const
G4double GetKineticEnergy() const
G4bool ConstPropertyExists(const G4String &key) const
G4double GetConstProperty(const G4String &key) const
std::size_t GetIndex() const
G4double GetEnergy(const G4double value) const
G4double GetMaxEnergy() const
G4double GetMaxValue() const
G4double GetGlobalTime() const
const G4ThreeVector & GetPosition() const
G4StepPoint * GetPostStepPoint() const
G4int GetTrackID() const
void SetTouchableHandle(const G4TouchableHandle &apValue)
const G4TouchableHandle & GetTouchableHandle() const
void SetParentID(const G4int aValue)
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
G4ParticleChange aParticleChange

◆ PreparePhysicsTable()

void G4OpWLS::PreparePhysicsTable ( const G4ParticleDefinition & )
overridevirtual

Reimplemented from G4VProcess.

Definition at line 78 of file G4OpWLS.cc.

78{ Initialise(); }

◆ SetVerboseLevel()

void G4OpWLS::SetVerboseLevel ( G4int verbose)

Definition at line 326 of file G4OpWLS.cc.

Referenced by Initialise().

◆ UseTimeProfile()

void G4OpWLS::UseTimeProfile ( const G4String name)
virtual

Definition at line 301 of file G4OpWLS.cc.

302{
304 {
306 WLSTimeGeneratorProfile = nullptr;
307 }
308 if(name == "delta")
309 {
310 WLSTimeGeneratorProfile = new G4WLSTimeGeneratorProfileDelta("delta");
311 }
312 else if(name == "exponential")
313 {
315 new G4WLSTimeGeneratorProfileExponential("exponential");
316 }
317 else
318 {
319 G4Exception("G4OpWLS::UseTimeProfile", "em0202", FatalException,
320 "generator does not exist");
321 }
323}
void SetWLSTimeProfile(const G4String &)

Referenced by Initialise().

Member Data Documentation

◆ theIntegralTable

G4PhysicsTable* G4OpWLS::theIntegralTable
protected

◆ WLSTimeGeneratorProfile

G4VWLSTimeGeneratorProfile* G4OpWLS::WLSTimeGeneratorProfile
protected

Definition at line 90 of file G4OpWLS.hh.

Referenced by G4OpWLS(), PostStepDoIt(), UseTimeProfile(), and ~G4OpWLS().


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