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

#include <G4OpWLS2.hh>

Inheritance diagram for G4OpWLS2:

Public Member Functions

 G4OpWLS2 (const G4String &processName="OpWLS2", G4ProcessType type=fOptical)
virtual ~G4OpWLS2 ()
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 G4OpWLS2.hh.

Constructor & Destructor Documentation

◆ G4OpWLS2()

G4OpWLS2::G4OpWLS2 ( const G4String & processName = "OpWLS2",
G4ProcessType type = fOptical )
explicit

Definition at line 54 of file G4OpWLS2.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 G4OpWLS2.cc:84
G4PhysicsTable * theIntegralTable
Definition G4OpWLS2.hh:91
G4VWLSTimeGeneratorProfile * WLSTimeGeneratorProfile
Definition G4OpWLS2.hh:90
G4VDiscreteProcess(const G4String &aName, G4ProcessType aType=fNotDefined)
G4int verboseLevel
void SetProcessSubType(G4int)
const G4String & GetProcessName() const

◆ ~G4OpWLS2()

G4OpWLS2::~G4OpWLS2 ( )
virtual

Definition at line 67 of file G4OpWLS2.cc.

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

Member Function Documentation

◆ BuildPhysicsTable()

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

Reimplemented from G4VProcess.

Definition at line 220 of file G4OpWLS2.cc.

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

Definition at line 114 of file G4OpWLS2.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 * G4OpWLS2::GetIntegralTable ( ) const
inlinevirtual

Definition at line 109 of file G4OpWLS2.hh.

110{
111 return theIntegralTable;
112}

◆ GetMeanFreePath()

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

Implements G4VDiscreteProcess.

Definition at line 285 of file G4OpWLS2.cc.

287{
288 G4double thePhotonEnergy = aTrack.GetDynamicParticle()->GetTotalEnergy();
289 G4double attLength = DBL_MAX;
290 G4MaterialPropertiesTable* MPT =
292
293 if(MPT)
294 {
296 if(attVector)
297 {
298 attLength = attVector->Value(thePhotonEnergy, idx_wls2);
299 }
300 }
301 return attLength;
302}
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 G4OpWLS2::Initialise ( )
virtual

Definition at line 84 of file G4OpWLS2.cc.

85{
86 G4OpticalParameters* params = G4OpticalParameters::Instance();
89}
void SetVerboseLevel(G4int)
Definition G4OpWLS2.cc:330
virtual void UseTimeProfile(const G4String name)
Definition G4OpWLS2.cc:305
G4int GetWLS2VerboseLevel() const
G4String GetWLS2TimeProfile() const
static G4OpticalParameters * Instance()

Referenced by G4OpWLS2(), and PreparePhysicsTable().

◆ IsApplicable()

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

Reimplemented from G4VProcess.

Definition at line 104 of file G4OpWLS2.hh.

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

◆ PostStepDoIt()

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

Reimplemented from G4VDiscreteProcess.

Definition at line 92 of file G4OpWLS2.cc.

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

Reimplemented from G4VProcess.

Definition at line 78 of file G4OpWLS2.cc.

79{
80 Initialise();
81}

◆ SetVerboseLevel()

void G4OpWLS2::SetVerboseLevel ( G4int verbose)

Definition at line 330 of file G4OpWLS2.cc.

Referenced by Initialise().

◆ UseTimeProfile()

void G4OpWLS2::UseTimeProfile ( const G4String name)
virtual

Definition at line 305 of file G4OpWLS2.cc.

306{
308 {
310 WLSTimeGeneratorProfile = nullptr;
311 }
312 if(name == "delta")
313 {
314 WLSTimeGeneratorProfile = new G4WLSTimeGeneratorProfileDelta("delta");
315 }
316 else if(name == "exponential")
317 {
319 new G4WLSTimeGeneratorProfileExponential("exponential");
320 }
321 else
322 {
323 G4Exception("G4OpWLS::UseTimeProfile", "em0202", FatalException,
324 "generator does not exist");
325 }
327}
void SetWLS2TimeProfile(const G4String &)

Referenced by Initialise().

Member Data Documentation

◆ theIntegralTable

G4PhysicsTable* G4OpWLS2::theIntegralTable
protected

◆ WLSTimeGeneratorProfile

G4VWLSTimeGeneratorProfile* G4OpWLS2::WLSTimeGeneratorProfile
protected

Definition at line 90 of file G4OpWLS2.hh.

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


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