Geant4 11.2.2
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 G4VDiscreteProcess
- 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
G4int verboseLevel
void SetProcessSubType(G4int)
const G4String & GetProcessName() const

◆ ~G4OpWLS()

G4OpWLS::~G4OpWLS ( )
virtual

Definition at line 67 of file G4OpWLS.cc.

68{
70 {
72 delete theIntegralTable;
73 }
75}
void clearAndDestroy()

Member Function Documentation

◆ BuildPhysicsTable()

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

Reimplemented from G4VProcess.

Definition at line 229 of file G4OpWLS.cc.

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

◆ DumpPhysicsTable()

void G4OpWLS::DumpPhysicsTable ( ) const
inlinevirtual

Definition at line 114 of file G4OpWLS.hh.

115{
116 std::size_t PhysicsTableSize = theIntegralTable->entries();
118
119 for(std::size_t i = 0; i < PhysicsTableSize; ++i)
120 {
122 v->DumpValues();
123 }
124}
std::size_t entries() const
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 294 of file G4OpWLS.cc.

296{
297 G4double thePhotonEnergy = aTrack.GetDynamicParticle()->GetTotalEnergy();
298 G4double attLength = DBL_MAX;
301
302 if(MPT)
303 {
305 if(attVector)
306 {
307 attLength = attVector->Value(thePhotonEnergy, idx_wls);
308 }
309 }
310 return attLength;
311}
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{
86}
void SetVerboseLevel(G4int)
Definition G4OpWLS.cc:339
virtual void UseTimeProfile(const G4String name)
Definition G4OpWLS.cc:314
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 std::vector<G4Track*> proposedSecondaries;
95
96 if(verboseLevel > 1)
97 {
98 G4cout << "\n** G4OpWLS: Photon absorbed! **" << G4endl;
99 }
100
101 G4StepPoint* pPostStepPoint = aStep.GetPostStepPoint();
104 if(!MPT)
105 {
106 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
107 }
108 if(!MPT->GetProperty(kWLSCOMPONENT))
109 {
110 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
111 }
112
113 G4int NumPhotons = 1;
115 {
116 G4double MeanNumberOfPhotons = MPT->GetConstProperty(kWLSMEANNUMBERPHOTONS);
117 NumPhotons = G4int(G4Poisson(MeanNumberOfPhotons));
118 if(NumPhotons <= 0)
119 {
120 // return unchanged particle and no secondaries
122 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
123 }
124 }
125
126 // Retrieve the WLS Integral for this material
127 // new G4PhysicsFreeVector allocated to hold CII's
128 G4double primaryEnergy = aTrack.GetDynamicParticle()->GetKineticEnergy();
129 G4double WLSTime = 0.;
130 G4PhysicsFreeVector* WLSIntegral = nullptr;
131
132 WLSTime = MPT->GetConstProperty(kWLSTIMECONSTANT);
133 WLSIntegral = (G4PhysicsFreeVector*) ((*theIntegralTable)(
134 aTrack.GetMaterial()->GetIndex()));
135
136 // Max WLS Integral
137 G4double CIImax = WLSIntegral->GetMaxValue();
138 G4int NumberOfPhotons = NumPhotons;
139
140 for(G4int i = 0; i < NumPhotons; ++i)
141 {
142 G4double sampledEnergy;
143 // Make sure the energy of the secondary is less than that of the primary
144 for(G4int j = 1; j <= 100; ++j)
145 {
146 // Determine photon energy
147 G4double CIIvalue = G4UniformRand() * CIImax;
148 sampledEnergy = WLSIntegral->GetEnergy(CIIvalue);
149 if(sampledEnergy <= primaryEnergy)
150 break;
151 }
152 // If no such energy can be sampled, return one less secondary, or none
153 if(sampledEnergy > primaryEnergy)
154 {
155 if(verboseLevel > 1)
156 {
157 G4cout << " *** G4OpWLS: One less WLS photon will be returned ***"
158 << G4endl;
159 }
160 NumberOfPhotons--;
161 if(NumberOfPhotons == 0)
162 {
163 if(verboseLevel > 1)
164 {
165 G4cout
166 << " *** G4OpWLS: No WLS photon can be sampled for this primary ***"
167 << G4endl;
168 }
169 // return unchanged particle and no secondaries
171 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
172 }
173 continue;
174 }
175 else if(verboseLevel > 1)
176 {
177 G4cout << "G4OpWLS: Created photon with energy: " << sampledEnergy
178 << G4endl;
179 }
180
181 // Generate random photon direction
182 G4double cost = 1. - 2. * G4UniformRand();
183 G4double sint = std::sqrt((1. - cost) * (1. + cost));
184 G4double phi = twopi * G4UniformRand();
185 G4double sinp = std::sin(phi);
186 G4double cosp = std::cos(phi);
187 G4ParticleMomentum photonMomentum(sint * cosp, sint * sinp, cost);
188
189 G4ThreeVector photonPolarization(cost * cosp, cost * sinp, -sint);
190 G4ThreeVector perp = photonMomentum.cross(photonPolarization);
191
192 phi = twopi * G4UniformRand();
193 sinp = std::sin(phi);
194 cosp = std::cos(phi);
195 photonPolarization = (cosp * photonPolarization + sinp * perp).unit();
196
197 // Generate a new photon:
198 auto sec_dp =
200 sec_dp->SetPolarization(photonPolarization);
201 sec_dp->SetKineticEnergy(sampledEnergy);
202
203 G4double secTime = pPostStepPoint->GetGlobalTime() +
205 G4ThreeVector secPos = pPostStepPoint->GetPosition();
206 G4Track* secTrack = new G4Track(sec_dp, secTime, secPos);
207
208 secTrack->SetTouchableHandle(aTrack.GetTouchableHandle());
209 secTrack->SetParentID(aTrack.GetTrackID());
210
211 proposedSecondaries.push_back(secTrack);
212 }
213
214 aParticleChange.SetNumberOfSecondaries((G4int)proposedSecondaries.size());
215 for(auto sec : proposedSecondaries)
216 {
218 }
219 if(verboseLevel > 1)
220 {
221 G4cout << "\n Exiting from G4OpWLS::DoIt -- NumberOfSecondaries = "
223 }
224
225 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
226}
G4long G4Poisson(G4double mean)
Definition G4Poisson.hh:50
@ 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
void AddSecondary(G4Track *aSecondary)
void Initialize(const G4Track &) override
G4double GetEnergy(const G4double value) 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 &)
void ProposeTrackStatus(G4TrackStatus status)
G4int GetNumberOfSecondaries() const
void SetNumberOfSecondaries(G4int totSecondaries)
G4ParticleChange aParticleChange
virtual G4double GenerateTime(const G4double time_constant)=0

◆ 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 339 of file G4OpWLS.cc.

Referenced by Initialise().

◆ UseTimeProfile()

void G4OpWLS::UseTimeProfile ( const G4String name)
virtual

Definition at line 314 of file G4OpWLS.cc.

315{
317 {
319 WLSTimeGeneratorProfile = nullptr;
320 }
321 if(name == "delta")
322 {
324 }
325 else if(name == "exponential")
326 {
328 new G4WLSTimeGeneratorProfileExponential("exponential");
329 }
330 else
331 {
332 G4Exception("G4OpWLS::UseTimeProfile", "em0202", FatalException,
333 "generator does not exist");
334 }
336}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
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: