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

#include <G4OpRayleigh.hh>

+ Inheritance diagram for G4OpRayleigh:

Public Member Functions

 G4OpRayleigh (const G4String &processName="OpRayleigh", G4ProcessType type=fOptical)
 
 ~G4OpRayleigh ()
 
G4bool IsApplicable (const G4ParticleDefinition &aParticleType)
 
G4double GetMeanFreePath (const G4Track &aTrack, G4double, G4ForceCondition *)
 
G4VParticleChangePostStepDoIt (const G4Track &aTrack, const G4Step &aStep)
 
G4PhysicsTableGetPhysicsTable () 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 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

G4PhysicsTablethePhysicsTable
 
- 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 76 of file G4OpRayleigh.hh.

Constructor & Destructor Documentation

◆ G4OpRayleigh()

G4OpRayleigh::G4OpRayleigh ( const G4String processName = "OpRayleigh",
G4ProcessType  type = fOptical 
)

Definition at line 82 of file G4OpRayleigh.cc.

83 : G4VDiscreteProcess(processName, type)
84{
86
88
89 DefaultWater = false;
90
91 if (verboseLevel>0) {
92 G4cout << GetProcessName() << " is created " << G4endl;
93 }
94
95 BuildThePhysicsTable();
96}
@ fOpRayleigh
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
G4PhysicsTable * thePhysicsTable
G4int verboseLevel
Definition: G4VProcess.hh:368
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:403
const G4String & GetProcessName() const
Definition: G4VProcess.hh:379

◆ ~G4OpRayleigh()

G4OpRayleigh::~G4OpRayleigh ( )

Definition at line 106 of file G4OpRayleigh.cc.

107{
108 if (thePhysicsTable!= 0) {
110 delete thePhysicsTable;
111 }
112}
void clearAndDestroy()

Member Function Documentation

◆ DumpPhysicsTable()

void G4OpRayleigh::DumpPhysicsTable ( ) const
inline

Definition at line 163 of file G4OpRayleigh.hh.

165{
166 G4int PhysicsTableSize = thePhysicsTable->entries();
168
169 for (G4int i = 0 ; i < PhysicsTableSize ; i++ )
170 {
172 v->DumpValues();
173 }
174}
int G4int
Definition: G4Types.hh:66
size_t entries() const

◆ GetMeanFreePath()

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

Implements G4VDiscreteProcess.

Definition at line 262 of file G4OpRayleigh.cc.

265{
266 const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
267 const G4Material* aMaterial = aTrack.GetMaterial();
268
269 G4double thePhotonEnergy = aParticle->GetTotalEnergy();
270
271 G4double AttenuationLength = DBL_MAX;
272
273 if (aMaterial->GetName() == "Water" && DefaultWater){
274
275 G4bool isOutRange;
276
277 AttenuationLength =
278 (*thePhysicsTable)(aMaterial->GetIndex())->
279 GetValue(thePhotonEnergy, isOutRange);
280 }
281 else {
282
283 G4MaterialPropertiesTable* aMaterialPropertyTable =
284 aMaterial->GetMaterialPropertiesTable();
285
286 if(aMaterialPropertyTable){
287 G4MaterialPropertyVector* AttenuationLengthVector =
288 aMaterialPropertyTable->GetProperty("RAYLEIGH");
289 if(AttenuationLengthVector){
290 AttenuationLength = AttenuationLengthVector ->
291 Value(thePhotonEnergy);
292 }
293 else{
294// G4cout << "No Rayleigh scattering length specified" << G4endl;
295 }
296 }
297 else{
298// G4cout << "No Rayleigh scattering length specified" << G4endl;
299 }
300 }
301
302 return AttenuationLength;
303}
double G4double
Definition: G4Types.hh:64
bool G4bool
Definition: G4Types.hh:67
G4double GetTotalEnergy() const
G4MaterialPropertyVector * GetProperty(const char *key)
G4MaterialPropertiesTable * GetMaterialPropertiesTable() const
Definition: G4Material.hh:251
const G4String & GetName() const
Definition: G4Material.hh:177
size_t GetIndex() const
Definition: G4Material.hh:261
G4Material * GetMaterial() const
const G4DynamicParticle * GetDynamicParticle() const
#define DBL_MAX
Definition: templates.hh:83

◆ GetPhysicsTable()

G4PhysicsTable * G4OpRayleigh::GetPhysicsTable ( ) const
inline

Definition at line 176 of file G4OpRayleigh.hh.

177{
178 return thePhysicsTable;
179}

◆ IsApplicable()

G4bool G4OpRayleigh::IsApplicable ( const G4ParticleDefinition aParticleType)
inlinevirtual

Reimplemented from G4VProcess.

Definition at line 157 of file G4OpRayleigh.hh.

158{
159 return ( &aParticleType == G4OpticalPhoton::OpticalPhoton() );
160}
static G4OpticalPhoton * OpticalPhoton()

◆ PostStepDoIt()

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

Reimplemented from G4VDiscreteProcess.

Definition at line 122 of file G4OpRayleigh.cc.

123{
125
126 const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
127
128 if (verboseLevel>0) {
129 G4cout << "Scattering Photon!" << G4endl;
130 G4cout << "Old Momentum Direction: "
131 << aParticle->GetMomentumDirection() << G4endl;
132 G4cout << "Old Polarization: "
133 << aParticle->GetPolarization() << G4endl;
134 }
135
136 G4double cosTheta;
137 G4ThreeVector OldMomentumDirection, NewMomentumDirection;
138 G4ThreeVector OldPolarization, NewPolarization;
139
140 do {
141 // Try to simulate the scattered photon momentum direction
142 // w.r.t. the initial photon momentum direction
143
144 G4double CosTheta = G4UniformRand();
145 G4double SinTheta = std::sqrt(1.-CosTheta*CosTheta);
146 // consider for the angle 90-180 degrees
147 if (G4UniformRand() < 0.5) CosTheta = -CosTheta;
148
149 // simulate the phi angle
150 G4double rand = twopi*G4UniformRand();
151 G4double SinPhi = std::sin(rand);
152 G4double CosPhi = std::cos(rand);
153
154 // start constructing the new momentum direction
155 G4double unit_x = SinTheta * CosPhi;
156 G4double unit_y = SinTheta * SinPhi;
157 G4double unit_z = CosTheta;
158 NewMomentumDirection.set (unit_x,unit_y,unit_z);
159
160 // Rotate the new momentum direction into global reference system
161 OldMomentumDirection = aParticle->GetMomentumDirection();
162 OldMomentumDirection = OldMomentumDirection.unit();
163 NewMomentumDirection.rotateUz(OldMomentumDirection);
164 NewMomentumDirection = NewMomentumDirection.unit();
165
166 // calculate the new polarization direction
167 // The new polarization needs to be in the same plane as the new
168 // momentum direction and the old polarization direction
169 OldPolarization = aParticle->GetPolarization();
170 G4double constant = -1./NewMomentumDirection.dot(OldPolarization);
171
172 NewPolarization = NewMomentumDirection + constant*OldPolarization;
173 NewPolarization = NewPolarization.unit();
174
175 // There is a corner case, where the Newmomentum direction
176 // is the same as oldpolariztion direction:
177 // random generate the azimuthal angle w.r.t. Newmomentum direction
178 if (NewPolarization.mag() == 0.) {
179 rand = G4UniformRand()*twopi;
180 NewPolarization.set(std::cos(rand),std::sin(rand),0.);
181 NewPolarization.rotateUz(NewMomentumDirection);
182 } else {
183 // There are two directions which are perpendicular
184 // to the new momentum direction
185 if (G4UniformRand() < 0.5) NewPolarization = -NewPolarization;
186 }
187
188 // simulate according to the distribution cos^2(theta)
189 cosTheta = NewPolarization.dot(OldPolarization);
190 } while (std::pow(cosTheta,2) < G4UniformRand());
191
192 aParticleChange.ProposePolarization(NewPolarization);
193 aParticleChange.ProposeMomentumDirection(NewMomentumDirection);
194
195 if (verboseLevel>0) {
196 G4cout << "New Polarization: "
197 << NewPolarization << G4endl;
198 G4cout << "Polarization Change: "
200 G4cout << "New Momentum Direction: "
201 << NewMomentumDirection << G4endl;
202 G4cout << "Momentum Change: "
204 }
205
206 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
207}
#define G4UniformRand()
Definition: Randomize.hh:53
Hep3Vector unit() const
double dot(const Hep3Vector &) const
double mag() const
void set(double x, double y, double z)
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:72
const G4ThreeVector & GetMomentumDirection() const
const G4ThreeVector & GetPolarization() const
void ProposePolarization(G4double Px, G4double Py, G4double Pz)
const G4ThreeVector * GetPolarization() const
const G4ThreeVector * GetMomentumDirection() const
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
virtual void Initialize(const G4Track &)
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:289

Member Data Documentation

◆ thePhysicsTable

G4PhysicsTable* G4OpRayleigh::thePhysicsTable
protected

Definition at line 141 of file G4OpRayleigh.hh.

Referenced by DumpPhysicsTable(), G4OpRayleigh(), GetPhysicsTable(), and ~G4OpRayleigh().


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