Geant4 9.6.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)
 
 ~G4OpWLS ()
 
G4bool IsApplicable (const G4ParticleDefinition &aParticleType)
 
G4double GetMeanFreePath (const G4Track &aTrack, G4double, G4ForceCondition *)
 
G4VParticleChangePostStepDoIt (const G4Track &aTrack, const G4Step &aStep)
 
G4PhysicsTableGetIntegralTable () const
 
void DumpPhysicsTable () const
 
void UseTimeProfile (const G4String name)
 
- Public Member Functions inherited from G4VDiscreteProcess
 G4VDiscreteProcess (const G4String &, G4ProcessType aType=fNotDefined)
 
 G4VDiscreteProcess (G4VDiscreteProcess &)
 
virtual ~G4VDiscreteProcess ()
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &)
 
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 &)
 
- Public Member Functions inherited from G4VProcess
 G4VProcess (const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
 
 G4VProcess (const G4VProcess &right)
 
virtual ~G4VProcess ()
 
G4int operator== (const G4VProcess &right) const
 
G4int operator!= (const G4VProcess &right) const
 
virtual G4VParticleChangePostStepDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)=0
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &track, G4ForceCondition *condition)=0
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)=0
 
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 IsApplicable (const G4ParticleDefinition &)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void PreparePhysicsTable (const G4ParticleDefinition &)
 
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 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
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 

Protected Attributes

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

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 
virtual G4double GetMeanFreePath (const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)=0
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double previousStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

Detailed Description

Definition at line 80 of file G4OpWLS.hh.

Constructor & Destructor Documentation

◆ G4OpWLS()

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

Definition at line 64 of file G4OpWLS.cc.

65 : G4VDiscreteProcess(processName, type)
66{
68
70
71 if (verboseLevel>0) {
72 G4cout << GetProcessName() << " is created " << G4endl;
73 }
74
76 new G4WLSTimeGeneratorProfileDelta("WLSTimeGeneratorProfileDelta");
77
78 BuildThePhysicsTable();
79}
@ fOpWLS
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
G4PhysicsTable * theIntegralTable
Definition: G4OpWLS.hh:140
G4VWLSTimeGeneratorProfile * WLSTimeGeneratorProfile
Definition: G4OpWLS.hh:139
G4int verboseLevel
Definition: G4VProcess.hh:368
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:403
const G4String & GetProcessName() const
Definition: G4VProcess.hh:379

◆ ~G4OpWLS()

G4OpWLS::~G4OpWLS ( )

Definition at line 85 of file G4OpWLS.cc.

86{
87 if (theIntegralTable != 0) {
89 delete theIntegralTable;
90 }
92}
void clearAndDestroy()

Member Function Documentation

◆ DumpPhysicsTable()

void G4OpWLS::DumpPhysicsTable ( ) const
inline

Definition at line 161 of file G4OpWLS.hh.

162{
163 G4int PhysicsTableSize = theIntegralTable->entries();
165
166 for (G4int i = 0 ; i < PhysicsTableSize ; i++ )
167 {
169 v->DumpValues();
170 }
171}
int G4int
Definition: G4Types.hh:66
size_t entries() const

◆ GetIntegralTable()

G4PhysicsTable * G4OpWLS::GetIntegralTable ( ) const
inline

Definition at line 155 of file G4OpWLS.hh.

156{
157 return theIntegralTable;
158}

◆ GetMeanFreePath()

G4double G4OpWLS::GetMeanFreePath ( const G4Track aTrack,
G4double  ,
G4ForceCondition  
)
virtual

Implements G4VDiscreteProcess.

Definition at line 351 of file G4OpWLS.cc.

354{
355 const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
356 const G4Material* aMaterial = aTrack.GetMaterial();
357
358 G4double thePhotonEnergy = aParticle->GetTotalEnergy();
359
360 G4MaterialPropertiesTable* aMaterialPropertyTable;
361 G4MaterialPropertyVector* AttenuationLengthVector;
362
363 G4double AttenuationLength = DBL_MAX;
364
365 aMaterialPropertyTable = aMaterial->GetMaterialPropertiesTable();
366
367 if ( aMaterialPropertyTable ) {
368 AttenuationLengthVector = aMaterialPropertyTable->
369 GetProperty("WLSABSLENGTH");
370 if ( AttenuationLengthVector ){
371 AttenuationLength = AttenuationLengthVector->
372 Value(thePhotonEnergy);
373 }
374 else {
375 // G4cout << "No WLS absorption length specified" << G4endl;
376 }
377 }
378 else {
379 // G4cout << "No WLS absortion length specified" << G4endl;
380 }
381
382 return AttenuationLength;
383}
double G4double
Definition: G4Types.hh:64
G4double GetTotalEnergy() const
G4MaterialPropertiesTable * GetMaterialPropertiesTable() const
Definition: G4Material.hh:251
G4Material * GetMaterial() const
const G4DynamicParticle * GetDynamicParticle() const
#define DBL_MAX
Definition: templates.hh:83

◆ IsApplicable()

G4bool G4OpWLS::IsApplicable ( const G4ParticleDefinition aParticleType)
inlinevirtual

Reimplemented from G4VProcess.

Definition at line 149 of file G4OpWLS.hh.

150{
151 return ( &aParticleType == G4OpticalPhoton::OpticalPhoton() );
152}
static G4OpticalPhoton * OpticalPhoton()

◆ PostStepDoIt()

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

Reimplemented from G4VDiscreteProcess.

Definition at line 102 of file G4OpWLS.cc.

103{
105
107
108 if (verboseLevel>0) {
109 G4cout << "\n** Photon absorbed! **" << G4endl;
110 }
111
112 const G4Material* aMaterial = aTrack.GetMaterial();
113
114 G4StepPoint* pPostStepPoint = aStep.GetPostStepPoint();
115
116 G4MaterialPropertiesTable* aMaterialPropertiesTable =
117 aMaterial->GetMaterialPropertiesTable();
118 if (!aMaterialPropertiesTable)
119 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
120
121 const G4MaterialPropertyVector* WLS_Intensity =
122 aMaterialPropertiesTable->GetProperty("WLSCOMPONENT");
123
124 if (!WLS_Intensity)
125 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
126
127 G4int NumPhotons = 1;
128
129 if (aMaterialPropertiesTable->ConstPropertyExists("WLSMEANNUMBERPHOTONS")) {
130
131 G4double MeanNumberOfPhotons = aMaterialPropertiesTable->
132 GetConstProperty("WLSMEANNUMBERPHOTONS");
133
134 NumPhotons = G4int(G4Poisson(MeanNumberOfPhotons));
135
136 if (NumPhotons <= 0) {
137
138 // return unchanged particle and no secondaries
139
141
142 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
143
144 }
145
146 }
147
149
150 G4int materialIndex = aMaterial->GetIndex();
151
152 // Retrieve the WLS Integral for this material
153 // new G4PhysicsOrderedFreeVector allocated to hold CII's
154
155 G4double WLSTime = 0.*ns;
156 G4PhysicsOrderedFreeVector* WLSIntegral = 0;
157
158 WLSTime = aMaterialPropertiesTable->
159 GetConstProperty("WLSTIMECONSTANT");
160 WLSIntegral =
161 (G4PhysicsOrderedFreeVector*)((*theIntegralTable)(materialIndex));
162
163 // Max WLS Integral
164
165 G4double CIImax = WLSIntegral->GetMaxValue();
166
167 for (G4int i = 0; i < NumPhotons; i++) {
168
169 // Determine photon energy
170
171 G4double CIIvalue = G4UniformRand()*CIImax;
172 G4double sampledEnergy =
173 WLSIntegral->GetEnergy(CIIvalue);
174
175 if (verboseLevel>1) {
176 G4cout << "sampledEnergy = " << sampledEnergy << G4endl;
177 G4cout << "CIIvalue = " << CIIvalue << G4endl;
178 }
179
180 // Generate random photon direction
181
182 G4double cost = 1. - 2.*G4UniformRand();
183 G4double sint = std::sqrt((1.-cost)*(1.+cost));
184
185 G4double phi = twopi*G4UniformRand();
186 G4double sinp = std::sin(phi);
187 G4double cosp = std::cos(phi);
188
189 G4double px = sint*cosp;
190 G4double py = sint*sinp;
191 G4double pz = cost;
192
193 // Create photon momentum direction vector
194
195 G4ParticleMomentum photonMomentum(px, py, pz);
196
197 // Determine polarization of new photon
198
199 G4double sx = cost*cosp;
200 G4double sy = cost*sinp;
201 G4double sz = -sint;
202
203 G4ThreeVector photonPolarization(sx, sy, sz);
204
205 G4ThreeVector perp = photonMomentum.cross(photonPolarization);
206
207 phi = twopi*G4UniformRand();
208 sinp = std::sin(phi);
209 cosp = std::cos(phi);
210
211 photonPolarization = cosp * photonPolarization + sinp * perp;
212
213 photonPolarization = photonPolarization.unit();
214
215 // Generate a new photon:
216
217 G4DynamicParticle* aWLSPhoton =
219 photonMomentum);
220 aWLSPhoton->SetPolarization
221 (photonPolarization.x(),
222 photonPolarization.y(),
223 photonPolarization.z());
224
225 aWLSPhoton->SetKineticEnergy(sampledEnergy);
226
227 // Generate new G4Track object:
228
229 // Must give position of WLS optical photon
230
231 G4double TimeDelay = WLSTimeGeneratorProfile->GenerateTime(WLSTime);
232 G4double aSecondaryTime = (pPostStepPoint->GetGlobalTime()) + TimeDelay;
233
234 G4ThreeVector aSecondaryPosition = pPostStepPoint->GetPosition();
235
236 G4Track* aSecondaryTrack =
237 new G4Track(aWLSPhoton,aSecondaryTime,aSecondaryPosition);
238
239 aSecondaryTrack->SetTouchableHandle(aTrack.GetTouchableHandle());
240 // aSecondaryTrack->SetTouchableHandle((G4VTouchable*)0);
241
242 aSecondaryTrack->SetParentID(aTrack.GetTrackID());
243
244 aParticleChange.AddSecondary(aSecondaryTrack);
245 }
246
247 if (verboseLevel>0) {
248 G4cout << "\n Exiting from G4OpWLS::DoIt -- NumberOfSecondaries = "
250 }
251
252 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
253}
G4long G4Poisson(G4double mean)
Definition: G4Poisson.hh:50
@ fStopAndKill
#define G4UniformRand()
Definition: Randomize.hh:53
Hep3Vector unit() const
Hep3Vector cross(const Hep3Vector &) const
void SetPolarization(G4double polX, G4double polY, G4double polZ)
void SetKineticEnergy(G4double aEnergy)
G4MaterialPropertyVector * GetProperty(const char *key)
G4bool ConstPropertyExists(const char *key)
size_t GetIndex() const
Definition: G4Material.hh:261
void AddSecondary(G4Track *aSecondary)
virtual void Initialize(const G4Track &)
G4double GetEnergy(G4double aValue)
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
Definition: G4VProcess.hh:289
virtual G4double GenerateTime(const G4double time_constant)=0
#define ns
Definition: xmlparse.cc:597

◆ UseTimeProfile()

void G4OpWLS::UseTimeProfile ( const G4String  name)

Definition at line 385 of file G4OpWLS.cc.

386{
387 if (name == "delta")
388 {
392 }
393 else if (name == "exponential")
394 {
397 new G4WLSTimeGeneratorProfileExponential("exponential");
398 }
399 else
400 {
401 G4Exception("G4OpWLS::UseTimeProfile", "em0202",
403 "generator does not exist");
404 }
405}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41

Referenced by G4OpticalPhysics::ConstructProcess(), and G4OpticalPhysics::SetWLSTimeProfile().

Member Data Documentation

◆ theIntegralTable

G4PhysicsTable* G4OpWLS::theIntegralTable
protected

Definition at line 140 of file G4OpWLS.hh.

Referenced by DumpPhysicsTable(), G4OpWLS(), GetIntegralTable(), PostStepDoIt(), and ~G4OpWLS().

◆ WLSTimeGeneratorProfile

G4VWLSTimeGeneratorProfile* G4OpWLS::WLSTimeGeneratorProfile
protected

Definition at line 139 of file G4OpWLS.hh.

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


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