Geant4 10.7.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)
 
virtual ~G4OpRayleigh ()
 
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 G4PhysicsTableGetPhysicsTable () const
 
virtual void DumpPhysicsTable () const
 
virtual void PreparePhysicsTable (const G4ParticleDefinition &) override
 
virtual void Initialise ()
 
- 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 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 ()
 
G4VProcessoperator= (const G4VProcess &)=delete
 
G4bool operator== (const G4VProcess &right) const
 
G4bool 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
 
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

G4PhysicsTablethePhysicsTable
 
- 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)
 
virtual G4double GetMeanFreePath (const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)=0
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double prevStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

Detailed Description

Definition at line 53 of file G4OpRayleigh.hh.

Constructor & Destructor Documentation

◆ G4OpRayleigh()

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

Definition at line 71 of file G4OpRayleigh.cc.

72 : G4VDiscreteProcess(processName, type)
73{
74 Initialise();
76 thePhysicsTable = nullptr;
77
78 if(verboseLevel > 0)
79 {
80 G4cout << GetProcessName() << " is created " << G4endl;
81 }
82}
@ fOpRayleigh
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4PhysicsTable * thePhysicsTable
Definition: G4OpRayleigh.hh:86
virtual void Initialise()
G4int verboseLevel
Definition: G4VProcess.hh:356
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:406
const G4String & GetProcessName() const
Definition: G4VProcess.hh:382

◆ ~G4OpRayleigh()

G4OpRayleigh::~G4OpRayleigh ( )
virtual

Definition at line 85 of file G4OpRayleigh.cc.

86{
87 // VI: inside this PhysicsTable all properties are unique
88 // it is not possible to destroy
90 {
91 delete thePhysicsTable;
92 }
93}

Member Function Documentation

◆ BuildPhysicsTable()

void G4OpRayleigh::BuildPhysicsTable ( const G4ParticleDefinition aParticleType)
overridevirtual

Reimplemented from G4VProcess.

Definition at line 191 of file G4OpRayleigh.cc.

192{
194 {
195 // thePhysicsTable->clearAndDestroy();
196 delete thePhysicsTable;
197 thePhysicsTable = nullptr;
198 }
199
200 const G4MaterialTable* theMaterialTable = G4Material::GetMaterialTable();
201 const size_t numOfMaterials = G4Material::GetNumberOfMaterials();
202 thePhysicsTable = new G4PhysicsTable(numOfMaterials);
203
204 for(size_t i = 0; i < numOfMaterials; ++i)
205 {
206 G4Material* material = (*theMaterialTable)[i];
208 G4PhysicsOrderedFreeVector* rayleigh = nullptr;
209 if(matProp)
210 {
211 rayleigh = matProp->GetProperty(kRAYLEIGH);
212 if(rayleigh == nullptr)
213 rayleigh = CalculateRayleighMeanFreePaths(material);
214 }
215 thePhysicsTable->insertAt(i, rayleigh);
216 }
217}
std::vector< G4Material * > G4MaterialTable
G4MaterialPropertyVector * GetProperty(const char *key, G4bool warning=false)
G4MaterialPropertiesTable * GetMaterialPropertiesTable() const
Definition: G4Material.hh:254
static size_t GetNumberOfMaterials()
Definition: G4Material.cc:644
static G4MaterialTable * GetMaterialTable()
Definition: G4Material.cc:637
void insertAt(std::size_t, G4PhysicsVector *)

◆ DumpPhysicsTable()

void G4OpRayleigh::DumpPhysicsTable ( ) const
inlinevirtual

Definition at line 110 of file G4OpRayleigh.hh.

111{
112 for(size_t i = 0; i < thePhysicsTable->entries(); ++i)
113 {
114 ((G4PhysicsOrderedFreeVector*) (*thePhysicsTable)[i])->DumpValues();
115 }
116}
std::size_t entries() const

◆ GetMeanFreePath()

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

Implements G4VDiscreteProcess.

Definition at line 220 of file G4OpRayleigh.cc.

222{
224 static_cast<G4PhysicsOrderedFreeVector*>(
225 (*thePhysicsTable)(aTrack.GetMaterial()->GetIndex()));
226
227 G4double rsLength = DBL_MAX;
228 if(rayleigh)
229 {
230 rsLength = rayleigh->Value(aTrack.GetDynamicParticle()->GetTotalMomentum(),
231 idx_rslength);
232 }
233 return rsLength;
234}
double G4double
Definition: G4Types.hh:83
G4double GetTotalMomentum() const
size_t GetIndex() const
Definition: G4Material.hh:258
G4double Value(G4double theEnergy, std::size_t &lastidx) const
G4Material * GetMaterial() const
const G4DynamicParticle * GetDynamicParticle() const
#define DBL_MAX
Definition: templates.hh:62

◆ GetPhysicsTable()

G4PhysicsTable * G4OpRayleigh::GetPhysicsTable ( ) const
inlinevirtual

Definition at line 118 of file G4OpRayleigh.hh.

119{
120 return thePhysicsTable;
121}

◆ Initialise()

void G4OpRayleigh::Initialise ( )
virtual

Definition at line 102 of file G4OpRayleigh.cc.

103{
104 SetVerboseLevel(G4OpticalParameters::Instance()->GetRayleighVerboseLevel());
105}
static G4OpticalParameters * Instance()
void SetVerboseLevel(G4int value)
Definition: G4VProcess.hh:412

Referenced by G4OpRayleigh(), and PreparePhysicsTable().

◆ IsApplicable()

G4bool G4OpRayleigh::IsApplicable ( const G4ParticleDefinition aParticleType)
inlineoverridevirtual

Reimplemented from G4VProcess.

Definition at line 104 of file G4OpRayleigh.hh.

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

◆ PostStepDoIt()

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

Reimplemented from G4VDiscreteProcess.

Definition at line 108 of file G4OpRayleigh.cc.

110{
112 const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
113
114 if(verboseLevel > 1)
115 {
116 G4cout << "OpRayleigh: Scattering Photon!" << G4endl
117 << "Old Momentum Direction: " << aParticle->GetMomentumDirection()
118 << G4endl << "Old Polarization: " << aParticle->GetPolarization()
119 << G4endl;
120 }
121
122 G4double cosTheta;
123 G4ThreeVector oldMomDir, newMomDir;
124 G4ThreeVector oldPol, newPol;
125 G4double rand;
126 G4double cost, sint, sinphi, cosphi;
127
128 do
129 {
130 // Try to simulate the scattered photon momentum direction
131 // w.r.t. the initial photon momentum direction
132 cost = G4UniformRand();
133 sint = std::sqrt(1. - cost * cost);
134 // consider for the angle 90-180 degrees
135 if(G4UniformRand() < 0.5)
136 cost = -cost;
137
138 // simulate the phi angle
139 rand = twopi * G4UniformRand();
140 sinphi = std::sin(rand);
141 cosphi = std::cos(rand);
142
143 // construct the new momentum direction
144 newMomDir.set(sint * cosphi, sint * sinphi, cost);
145 oldMomDir = aParticle->GetMomentumDirection();
146 newMomDir.rotateUz(oldMomDir);
147
148 // calculate the new polarization direction
149 // The new polarization needs to be in the same plane as the new
150 // momentum direction and the old polarization direction
151 oldPol = aParticle->GetPolarization();
152 newPol = (oldPol - newMomDir.dot(oldPol) * newMomDir).unit();
153
154 // There is a corner case, where the new momentum direction
155 // is the same as old polarization direction:
156 // random generate the azimuthal angle w.r.t. new momentum direction
157 if(newPol.mag() == 0.)
158 {
159 rand = G4UniformRand() * twopi;
160 newPol.set(std::cos(rand), std::sin(rand), 0.);
161 newPol.rotateUz(newMomDir);
162 }
163 else
164 {
165 // There are two directions perpendicular to the new momentum direction
166 if(G4UniformRand() < 0.5)
167 newPol = -newPol;
168 }
169
170 // simulate according to the distribution cos^2(theta)
171 cosTheta = newPol.dot(oldPol);
172 // Loop checking, 13-Aug-2015, Peter Gumplinger
173 } while(std::pow(cosTheta, 2) < G4UniformRand());
174
177
178 if(verboseLevel > 1)
179 {
180 G4cout << "New Polarization: " << newPol << G4endl
181 << "Polarization Change: " << *(aParticleChange.GetPolarization())
182 << G4endl << "New Momentum Direction: " << newMomDir << G4endl
183 << "Momentum Change: " << *(aParticleChange.GetMomentumDirection())
184 << G4endl;
185 }
186
187 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
188}
#define G4UniformRand()
Definition: Randomize.hh:52
double dot(const Hep3Vector &) const
double mag() const
void set(double x, double y, double z)
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:33
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:327

◆ PreparePhysicsTable()

void G4OpRayleigh::PreparePhysicsTable ( const G4ParticleDefinition )
overridevirtual

Reimplemented from G4VProcess.

Definition at line 96 of file G4OpRayleigh.cc.

97{
98 Initialise();
99}

Member Data Documentation

◆ thePhysicsTable

G4PhysicsTable* G4OpRayleigh::thePhysicsTable
protected

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