Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions | Static Public Member Functions | Protected Attributes | Static Protected Attributes
G4OpWLS Class Reference

#include <G4OpWLS.hh>

Inheritance diagram for G4OpWLS:
G4VDiscreteProcess G4VProcess

Public Member Functions

 G4OpWLS (const G4String &processName="OpWLS", G4ProcessType type=fOptical)
 
 ~G4OpWLS ()
 
G4bool IsApplicable (const G4ParticleDefinition &aParticleType)
 
void BuildPhysicsTable (const G4ParticleDefinition &aParticleType)
 
G4double GetMeanFreePath (const G4Track &aTrack, G4double, G4ForceCondition *)
 
G4VParticleChangePostStepDoIt (const G4Track &aTrack, const G4Step &aStep)
 
G4PhysicsTableGetIntegralTable () const
 
void DumpPhysicsTable () const
 
- 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 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
 
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 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
 
virtual void SetMasterProcess (G4VProcess *masterP)
 
const G4VProcessGetMasterProcess () const
 
virtual void BuildWorkerPhysicsTable (const G4ParticleDefinition &part)
 
virtual void PrepareWorkerPhysicsTable (const G4ParticleDefinition &)
 

Static Public Member Functions

static void UseTimeProfile (const G4String name)
 
- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 

Protected Attributes

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
 

Static Protected Attributes

static G4VWLSTimeGeneratorProfileWLSTimeGeneratorProfile = 0
 

Additional Inherited Members

- 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 ( const G4String processName = "OpWLS",
G4ProcessType  type = fOptical 
)

Definition at line 70 of file G4OpWLS.cc.

References fOpWLS, G4cout, G4endl, G4VProcess::GetProcessName(), G4VProcess::SetProcessSubType(), theIntegralTable, G4VProcess::verboseLevel, and WLSTimeGeneratorProfile.

71  : G4VDiscreteProcess(processName, type)
72 {
74 
75  theIntegralTable = NULL;
77  { WLSTimeGeneratorProfile = new G4WLSTimeGeneratorProfileDelta("WLSTimeGeneratorProfileDelta"); }
78 
79  if (verboseLevel>0) {
80  G4cout << GetProcessName() << " is created " << G4endl;
81  }
82 }
G4int verboseLevel
Definition: G4VProcess.hh:368
static G4VWLSTimeGeneratorProfile * WLSTimeGeneratorProfile
Definition: G4OpWLS.hh:142
G4GLOB_DLL std::ostream G4cout
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:432
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
#define G4endl
Definition: G4ios.hh:61
G4PhysicsTable * theIntegralTable
Definition: G4OpWLS.hh:143
G4OpWLS::~G4OpWLS ( )

Definition at line 88 of file G4OpWLS.cc.

References G4PhysicsTable::clearAndDestroy(), and theIntegralTable.

89 {
90  if (theIntegralTable != 0) {
92  delete theIntegralTable;
93  }
94 }
G4PhysicsTable * theIntegralTable
Definition: G4OpWLS.hh:143
void clearAndDestroy()

Member Function Documentation

void G4OpWLS::BuildPhysicsTable ( const G4ParticleDefinition aParticleType)
virtual

Reimplemented from G4VProcess.

Definition at line 100 of file G4OpWLS.cc.

References theIntegralTable.

101 {
102  if (!theIntegralTable) BuildThePhysicsTable();
103 }
G4PhysicsTable * theIntegralTable
Definition: G4OpWLS.hh:143
void G4OpWLS::DumpPhysicsTable ( ) const
inline

Definition at line 164 of file G4OpWLS.hh.

References G4PhysicsOrderedFreeVector::DumpValues(), G4PhysicsTable::entries(), theIntegralTable, and test::v.

165 {
166  G4int PhysicsTableSize = theIntegralTable->entries();
168 
169  for (G4int i = 0 ; i < PhysicsTableSize ; i++ )
170  {
172  v->DumpValues();
173  }
174 }
int G4int
Definition: G4Types.hh:78
G4PhysicsTable * theIntegralTable
Definition: G4OpWLS.hh:143
size_t entries() const
G4PhysicsTable * G4OpWLS::GetIntegralTable ( ) const
inline

Definition at line 158 of file G4OpWLS.hh.

References theIntegralTable.

159 {
160  return theIntegralTable;
161 }
G4PhysicsTable * theIntegralTable
Definition: G4OpWLS.hh:143
G4double G4OpWLS::GetMeanFreePath ( const G4Track aTrack,
G4double  ,
G4ForceCondition  
)
virtual

Implements G4VDiscreteProcess.

Definition at line 387 of file G4OpWLS.cc.

References DBL_MAX, G4Track::GetDynamicParticle(), G4Track::GetMaterial(), G4Material::GetMaterialPropertiesTable(), and G4DynamicParticle::GetTotalEnergy().

390 {
391  const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
392  const G4Material* aMaterial = aTrack.GetMaterial();
393 
394  G4double thePhotonEnergy = aParticle->GetTotalEnergy();
395 
396  G4MaterialPropertiesTable* aMaterialPropertyTable;
397  G4MaterialPropertyVector* AttenuationLengthVector;
398 
399  G4double AttenuationLength = DBL_MAX;
400 
401  aMaterialPropertyTable = aMaterial->GetMaterialPropertiesTable();
402 
403  if ( aMaterialPropertyTable ) {
404  AttenuationLengthVector = aMaterialPropertyTable->
405  GetProperty("WLSABSLENGTH");
406  if ( AttenuationLengthVector ){
407  AttenuationLength = AttenuationLengthVector->
408  Value(thePhotonEnergy);
409  }
410  else {
411  // G4cout << "No WLS absorption length specified" << G4endl;
412  }
413  }
414  else {
415  // G4cout << "No WLS absortion length specified" << G4endl;
416  }
417 
418  return AttenuationLength;
419 }
G4double GetTotalEnergy() const
const G4DynamicParticle * GetDynamicParticle() const
G4Material * GetMaterial() const
G4MaterialPropertiesTable * GetMaterialPropertiesTable() const
Definition: G4Material.hh:250
double G4double
Definition: G4Types.hh:76
#define DBL_MAX
Definition: templates.hh:83
G4bool G4OpWLS::IsApplicable ( const G4ParticleDefinition aParticleType)
inlinevirtual

Reimplemented from G4VProcess.

Definition at line 152 of file G4OpWLS.hh.

References G4OpticalPhoton::OpticalPhoton().

153 {
154  return ( &aParticleType == G4OpticalPhoton::OpticalPhoton() );
155 }
static G4OpticalPhoton * OpticalPhoton()
G4VParticleChange * G4OpWLS::PostStepDoIt ( const G4Track aTrack,
const G4Step aStep 
)
virtual

Reimplemented from G4VDiscreteProcess.

Definition at line 109 of file G4OpWLS.cc.

References G4ParticleChange::AddSecondary(), G4VProcess::aParticleChange, G4MaterialPropertiesTable::ConstPropertyExists(), CLHEP::Hep3Vector::cross(), fStopAndKill, G4cout, G4endl, G4Poisson(), G4UniformRand, G4VWLSTimeGeneratorProfile::GenerateTime(), G4Track::GetDynamicParticle(), G4PhysicsOrderedFreeVector::GetEnergy(), G4StepPoint::GetGlobalTime(), G4Material::GetIndex(), G4DynamicParticle::GetKineticEnergy(), G4Track::GetMaterial(), G4Material::GetMaterialPropertiesTable(), G4PhysicsOrderedFreeVector::GetMaxValue(), G4VParticleChange::GetNumberOfSecondaries(), G4StepPoint::GetPosition(), G4Step::GetPostStepPoint(), G4MaterialPropertiesTable::GetProperty(), G4Track::GetTouchableHandle(), G4Track::GetTrackID(), G4ParticleChange::Initialize(), ns, G4OpticalPhoton::OpticalPhoton(), G4VDiscreteProcess::PostStepDoIt(), G4VParticleChange::ProposeTrackStatus(), G4DynamicParticle::SetKineticEnergy(), G4VParticleChange::SetNumberOfSecondaries(), G4Track::SetParentID(), G4DynamicParticle::SetPolarization(), G4Track::SetTouchableHandle(), theIntegralTable, python.hepunit::twopi, CLHEP::Hep3Vector::unit(), G4VProcess::verboseLevel, WLSTimeGeneratorProfile, CLHEP::Hep3Vector::x(), CLHEP::Hep3Vector::y(), and CLHEP::Hep3Vector::z().

110 {
111  aParticleChange.Initialize(aTrack);
112 
114 
115  if (verboseLevel>0) {
116  G4cout << "\n** Photon absorbed! **" << G4endl;
117  }
118 
119  const G4Material* aMaterial = aTrack.GetMaterial();
120 
121  G4StepPoint* pPostStepPoint = aStep.GetPostStepPoint();
122 
123  G4MaterialPropertiesTable* aMaterialPropertiesTable =
124  aMaterial->GetMaterialPropertiesTable();
125  if (!aMaterialPropertiesTable)
126  return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
127 
128  const G4MaterialPropertyVector* WLS_Intensity =
129  aMaterialPropertiesTable->GetProperty("WLSCOMPONENT");
130 
131  if (!WLS_Intensity)
132  return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
133 
134  G4int NumPhotons = 1;
135 
136  if (aMaterialPropertiesTable->ConstPropertyExists("WLSMEANNUMBERPHOTONS")) {
137 
138  G4double MeanNumberOfPhotons = aMaterialPropertiesTable->
139  GetConstProperty("WLSMEANNUMBERPHOTONS");
140 
141  NumPhotons = G4int(G4Poisson(MeanNumberOfPhotons));
142 
143  if (NumPhotons <= 0) {
144 
145  // return unchanged particle and no secondaries
146 
148 
149  return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
150 
151  }
152 
153  }
154 
156 
157  G4double primaryEnergy = aTrack.GetDynamicParticle()->GetKineticEnergy();
158 
159  G4int materialIndex = aMaterial->GetIndex();
160 
161  // Retrieve the WLS Integral for this material
162  // new G4PhysicsOrderedFreeVector allocated to hold CII's
163 
164  G4double WLSTime = 0.*ns;
165  G4PhysicsOrderedFreeVector* WLSIntegral = 0;
166 
167  WLSTime = aMaterialPropertiesTable->
168  GetConstProperty("WLSTIMECONSTANT");
169  WLSIntegral =
170  (G4PhysicsOrderedFreeVector*)((*theIntegralTable)(materialIndex));
171 
172  // Max WLS Integral
173 
174  G4double CIImax = WLSIntegral->GetMaxValue();
175 
176  G4int NumberOfPhotons = NumPhotons;
177 
178  for (G4int i = 0; i < NumPhotons; i++) {
179 
180  G4double sampledEnergy;
181 
182  // Make sure the energy of the secondary is less than that of the primary
183 
184  for (G4int j = 1; j <= 100; j++) {
185 
186  // Determine photon energy
187 
188  G4double CIIvalue = G4UniformRand()*CIImax;
189  sampledEnergy = WLSIntegral->GetEnergy(CIIvalue);
190 
191  if (verboseLevel>1) {
192  G4cout << "sampledEnergy = " << sampledEnergy << G4endl;
193  G4cout << "CIIvalue = " << CIIvalue << G4endl;
194  }
195 
196  if (sampledEnergy <= primaryEnergy) break;
197  }
198 
199  // If no such energy can be sampled, return one less secondary, or none
200 
201  if (sampledEnergy > primaryEnergy) {
202  if (verboseLevel>1)
203  G4cout << " *** One less WLS photon will be returned ***" << G4endl;
204  NumberOfPhotons--;
205  aParticleChange.SetNumberOfSecondaries(NumberOfPhotons);
206  if (NumberOfPhotons == 0) {
207  if (verboseLevel>1)
208  G4cout << " *** No WLS photon can be sampled for this primary ***"
209  << G4endl;
210  // return unchanged particle and no secondaries
211  return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
212  }
213  continue;
214  }
215 
216  // Generate random photon direction
217 
218  G4double cost = 1. - 2.*G4UniformRand();
219  G4double sint = std::sqrt((1.-cost)*(1.+cost));
220 
221  G4double phi = twopi*G4UniformRand();
222  G4double sinp = std::sin(phi);
223  G4double cosp = std::cos(phi);
224 
225  G4double px = sint*cosp;
226  G4double py = sint*sinp;
227  G4double pz = cost;
228 
229  // Create photon momentum direction vector
230 
231  G4ParticleMomentum photonMomentum(px, py, pz);
232 
233  // Determine polarization of new photon
234 
235  G4double sx = cost*cosp;
236  G4double sy = cost*sinp;
237  G4double sz = -sint;
238 
239  G4ThreeVector photonPolarization(sx, sy, sz);
240 
241  G4ThreeVector perp = photonMomentum.cross(photonPolarization);
242 
243  phi = twopi*G4UniformRand();
244  sinp = std::sin(phi);
245  cosp = std::cos(phi);
246 
247  photonPolarization = cosp * photonPolarization + sinp * perp;
248 
249  photonPolarization = photonPolarization.unit();
250 
251  // Generate a new photon:
252 
253  G4DynamicParticle* aWLSPhoton =
255  photonMomentum);
256  aWLSPhoton->SetPolarization
257  (photonPolarization.x(),
258  photonPolarization.y(),
259  photonPolarization.z());
260 
261  aWLSPhoton->SetKineticEnergy(sampledEnergy);
262 
263  // Generate new G4Track object:
264 
265  // Must give position of WLS optical photon
266 
267  G4double TimeDelay = WLSTimeGeneratorProfile->GenerateTime(WLSTime);
268  G4double aSecondaryTime = (pPostStepPoint->GetGlobalTime()) + TimeDelay;
269 
270  G4ThreeVector aSecondaryPosition = pPostStepPoint->GetPosition();
271 
272  G4Track* aSecondaryTrack =
273  new G4Track(aWLSPhoton,aSecondaryTime,aSecondaryPosition);
274 
275  aSecondaryTrack->SetTouchableHandle(aTrack.GetTouchableHandle());
276  // aSecondaryTrack->SetTouchableHandle((G4VTouchable*)0);
277 
278  aSecondaryTrack->SetParentID(aTrack.GetTrackID());
279 
280  aParticleChange.AddSecondary(aSecondaryTrack);
281  }
282 
283  if (verboseLevel>0) {
284  G4cout << "\n Exiting from G4OpWLS::DoIt -- NumberOfSecondaries = "
286  }
287 
288  return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
289 }
G4long G4Poisson(G4double mean)
Definition: G4Poisson.hh:51
G4int GetNumberOfSecondaries() const
G4MaterialPropertyVector * GetProperty(const char *key)
G4int verboseLevel
Definition: G4VProcess.hh:368
G4double GetKineticEnergy() const
const G4DynamicParticle * GetDynamicParticle() const
size_t GetIndex() const
Definition: G4Material.hh:260
static G4VWLSTimeGeneratorProfile * WLSTimeGeneratorProfile
Definition: G4OpWLS.hh:142
void SetTouchableHandle(const G4TouchableHandle &apValue)
int G4int
Definition: G4Types.hh:78
#define G4UniformRand()
Definition: Randomize.hh:87
G4GLOB_DLL std::ostream G4cout
virtual G4double GenerateTime(const G4double time_constant)=0
const G4ThreeVector & GetPosition() const
G4int GetTrackID() const
void SetPolarization(G4double polX, G4double polY, G4double polZ)
void SetKineticEnergy(G4double aEnergy)
const G4TouchableHandle & GetTouchableHandle() const
G4double GetEnergy(G4double aValue)
G4Material * GetMaterial() const
static G4OpticalPhoton * OpticalPhoton()
virtual void Initialize(const G4Track &)
G4MaterialPropertiesTable * GetMaterialPropertiesTable() const
Definition: G4Material.hh:250
void SetNumberOfSecondaries(G4int totSecondaries)
Hep3Vector unit() const
void SetParentID(const G4int aValue)
G4StepPoint * GetPostStepPoint() const
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:289
G4bool ConstPropertyExists(const char *key)
void AddSecondary(G4Track *aSecondary)
#define G4endl
Definition: G4ios.hh:61
G4double GetGlobalTime() const
Hep3Vector cross(const Hep3Vector &) const
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
G4PhysicsTable * theIntegralTable
Definition: G4OpWLS.hh:143
#define ns
Definition: xmlparse.cc:597
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
void G4OpWLS::UseTimeProfile ( const G4String  name)
static

Definition at line 426 of file G4OpWLS.cc.

References FatalException, G4Exception(), and WLSTimeGeneratorProfile.

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

427 {
428  G4AutoLock l(&wlsCmdHandlingMutex);
429  if (name == "delta")
430  {
433  new G4WLSTimeGeneratorProfileDelta("delta");
434  }
435  else if (name == "exponential")
436  {
439  new G4WLSTimeGeneratorProfileExponential("exponential");
440  }
441  else
442  {
443  G4Exception("G4OpWLS::UseTimeProfile", "em0202",
445  "generator does not exist");
446  }
447 }
static G4VWLSTimeGeneratorProfile * WLSTimeGeneratorProfile
Definition: G4OpWLS.hh:142
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41

Field Documentation

G4PhysicsTable* G4OpWLS::theIntegralTable
protected
G4VWLSTimeGeneratorProfile * G4OpWLS::WLSTimeGeneratorProfile = 0
staticprotected

Definition at line 142 of file G4OpWLS.hh.

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


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