Geant4 11.2.2
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 ()
 
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

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)
 
- Protected Member Functions inherited from G4VDiscreteProcess
- 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:67
G4GLOB_DLL std::ostream G4cout
G4PhysicsTable * thePhysicsTable
virtual void Initialise()
G4int verboseLevel
void SetProcessSubType(G4int)
const G4String & GetProcessName() const

◆ ~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
89 delete thePhysicsTable;
90}

Member Function Documentation

◆ BuildPhysicsTable()

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

Reimplemented from G4VProcess.

Definition at line 188 of file G4OpRayleigh.cc.

189{
191 {
192 // thePhysicsTable->clearAndDestroy();
193 delete thePhysicsTable;
194 thePhysicsTable = nullptr;
195 }
196
197 const G4MaterialTable* theMaterialTable = G4Material::GetMaterialTable();
198 const size_t numOfMaterials = G4Material::GetNumberOfMaterials();
199 thePhysicsTable = new G4PhysicsTable(numOfMaterials);
200
201 for(size_t i = 0; i < numOfMaterials; ++i)
202 {
203 G4Material* material = (*theMaterialTable)[i];
205 G4PhysicsFreeVector* rayleigh = nullptr;
206 if(matProp)
207 {
208 rayleigh = matProp->GetProperty(kRAYLEIGH);
209 if(rayleigh == nullptr)
210 rayleigh = CalculateRayleighMeanFreePaths(material);
211 }
212 thePhysicsTable->insertAt(i, rayleigh);
213 }
214}
std::vector< G4Material * > G4MaterialTable
G4MaterialPropertyVector * GetProperty(const char *key) const
G4MaterialPropertiesTable * GetMaterialPropertiesTable() const
static std::size_t GetNumberOfMaterials()
static G4MaterialTable * GetMaterialTable()
void insertAt(std::size_t, G4PhysicsVector *)

◆ DumpPhysicsTable()

void G4OpRayleigh::DumpPhysicsTable ( ) const
inlinevirtual

Definition at line 112 of file G4OpRayleigh.hh.

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

◆ GetMeanFreePath()

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

Implements G4VDiscreteProcess.

Definition at line 217 of file G4OpRayleigh.cc.

219{
220 auto rayleigh = static_cast<G4PhysicsFreeVector*>(
221 (*thePhysicsTable)(aTrack.GetMaterial()->GetIndex()));
222
223 G4double rsLength = DBL_MAX;
224 if(rayleigh)
225 {
226 rsLength = rayleigh->Value(aTrack.GetDynamicParticle()->GetTotalMomentum(),
227 idx_rslength);
228 }
229 return rsLength;
230}
double G4double
Definition G4Types.hh:83
G4double GetTotalMomentum() const
G4Material * GetMaterial() const
const G4DynamicParticle * GetDynamicParticle() const
#define DBL_MAX
Definition templates.hh:62

◆ GetPhysicsTable()

G4PhysicsTable * G4OpRayleigh::GetPhysicsTable ( ) const
inlinevirtual

Definition at line 120 of file G4OpRayleigh.hh.

121{
122 return thePhysicsTable;
123}

◆ Initialise()

void G4OpRayleigh::Initialise ( )
virtual

Definition at line 99 of file G4OpRayleigh.cc.

100{
101 SetVerboseLevel(G4OpticalParameters::Instance()->GetRayleighVerboseLevel());
102}
void SetVerboseLevel(G4int)
static G4OpticalParameters * Instance()

Referenced by G4OpRayleigh(), and PreparePhysicsTable().

◆ IsApplicable()

G4bool G4OpRayleigh::IsApplicable ( const G4ParticleDefinition & aParticleType)
inlineoverridevirtual

Reimplemented from G4VProcess.

Definition at line 106 of file G4OpRayleigh.hh.

108{
109 return (&aParticleType == G4OpticalPhoton::OpticalPhoton());
110}
static G4OpticalPhoton * OpticalPhoton()

◆ PostStepDoIt()

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

Reimplemented from G4VDiscreteProcess.

Definition at line 105 of file G4OpRayleigh.cc.

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

◆ PreparePhysicsTable()

void G4OpRayleigh::PreparePhysicsTable ( const G4ParticleDefinition & )
overridevirtual

Reimplemented from G4VProcess.

Definition at line 93 of file G4OpRayleigh.cc.

94{
95 Initialise();
96}

◆ SetVerboseLevel()

void G4OpRayleigh::SetVerboseLevel ( G4int verbose)

Definition at line 308 of file G4OpRayleigh.cc.

Referenced by Initialise().

Member Data Documentation

◆ thePhysicsTable

G4PhysicsTable* G4OpRayleigh::thePhysicsTable
protected

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