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

#include <G4OpMieHG.hh>

+ Inheritance diagram for G4OpMieHG:

Public Member Functions

 G4OpMieHG (const G4String &processName="OpMieHG", G4ProcessType type=fOptical)
 
virtual ~G4OpMieHG ()
 
virtual G4bool IsApplicable (const G4ParticleDefinition &aParticleType) override
 
virtual G4double GetMeanFreePath (const G4Track &aTrack, G4double, G4ForceCondition *) override
 
virtual G4VParticleChangePostStepDoIt (const G4Track &aTrack, const G4Step &aStep) override
 
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 void BuildPhysicsTable (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 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 &)
 

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

Detailed Description

Definition at line 46 of file G4OpMieHG.hh.

Constructor & Destructor Documentation

◆ G4OpMieHG()

G4OpMieHG::G4OpMieHG ( const G4String & processName = "OpMieHG",
G4ProcessType type = fOptical )
explicit

Definition at line 47 of file G4OpMieHG.cc.

48 : G4VDiscreteProcess(processName, type)
49{
50 Initialise();
51 if(verboseLevel > 0)
52 {
53 G4cout << GetProcessName() << " is created " << G4endl;
54 }
56}
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
virtual void Initialise()
Definition G4OpMieHG.cc:68
G4int verboseLevel
void SetProcessSubType(G4int)
const G4String & GetProcessName() const

◆ ~G4OpMieHG()

G4OpMieHG::~G4OpMieHG ( )
virtualdefault

Member Function Documentation

◆ GetMeanFreePath()

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

Implements G4VDiscreteProcess.

Definition at line 168 of file G4OpMieHG.cc.

170{
171 G4double attLength = DBL_MAX;
174 if(MPT)
175 {
177 if(attVector)
178 {
179 attLength = attVector->Value(
180 aTrack.GetDynamicParticle()->GetTotalEnergy(), idx_mie);
181 }
182 }
183 return attLength;
184}
double G4double
Definition G4Types.hh:83
G4double GetTotalEnergy() const
G4MaterialPropertyVector * GetProperty(const char *key) const
G4MaterialPropertiesTable * GetMaterialPropertiesTable() const
G4double Value(const G4double energy, std::size_t &lastidx) const
G4Material * GetMaterial() const
const G4DynamicParticle * GetDynamicParticle() const
#define DBL_MAX
Definition templates.hh:62

◆ Initialise()

void G4OpMieHG::Initialise ( )
virtual

Definition at line 68 of file G4OpMieHG.cc.

69{
70 SetVerboseLevel(G4OpticalParameters::Instance()->GetMieVerboseLevel());
71}
void SetVerboseLevel(G4int)
Definition G4OpMieHG.cc:187
static G4OpticalParameters * Instance()

Referenced by G4OpMieHG(), and PreparePhysicsTable().

◆ IsApplicable()

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

Reimplemented from G4VProcess.

Definition at line 77 of file G4OpMieHG.hh.

78{
79 return (&aParticleType == G4OpticalPhoton::OpticalPhoton());
80}
static G4OpticalPhoton * OpticalPhoton()

◆ PostStepDoIt()

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

Reimplemented from G4VDiscreteProcess.

Definition at line 74 of file G4OpMieHG.cc.

76{
78
79 const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
80 const G4MaterialPropertiesTable* MPT =
82
84
85 if(verboseLevel > 1)
86 {
87 G4cout << "OpMie Scattering Photon!" << G4endl
88 << " Old Momentum Direction: " << aParticle->GetMomentumDirection()
89 << G4endl
90 << " MIE Old Polarization: " << aParticle->GetPolarization()
91 << G4endl;
92 }
93
94 G4double gg;
95 G4int direction;
96 if(G4UniformRand() <= forwardRatio)
97 {
99 direction = 1;
100 }
101 else
102 {
104 direction = -1;
105 }
106
108
109 // sample the direction
110 G4double theta;
111 if(gg != 0.)
112 {
113 theta = std::acos(2. * r * (1. + gg) * (1. + gg) * (1. - gg + gg * r) /
114 ((1. - gg + 2. * gg * r) * (1. - gg + 2. * gg * r)) -
115 1.);
116 }
117 else
118 {
119 theta = std::acos(2. * r - 1.);
120 }
121 G4double phi = G4UniformRand() * twopi;
122
123 if(direction == -1)
124 theta = pi - theta; // backward scattering
125
126 G4ThreeVector newMomDir, oldMomDir;
127 G4ThreeVector newPol, oldPol;
128
129 G4double sinth = std::sin(theta);
130 newMomDir.set(sinth * std::cos(phi), sinth * std::sin(phi), std::cos(theta));
131 oldMomDir = aParticle->GetMomentumDirection();
132 newMomDir.rotateUz(oldMomDir);
133 newMomDir = newMomDir.unit();
134
135 oldPol = aParticle->GetPolarization();
136 newPol = newMomDir - oldPol / newMomDir.dot(oldPol);
137 newPol = newPol.unit();
138
139 if(newPol.mag() == 0.)
140 {
141 r = G4UniformRand() * twopi;
142 newPol.set(std::cos(r), std::sin(r), 0.);
143 newPol.rotateUz(newMomDir);
144 }
145 else
146 {
147 // There are two directions perpendicular to new momentum direction
148 if(G4UniformRand() < 0.5)
149 newPol = -newPol;
150 }
151
154
155 if(verboseLevel > 1)
156 {
157 G4cout << "OpMie New Polarization: " << newPol << G4endl
158 << " Polarization Change: " << *(aParticleChange.GetPolarization())
159 << G4endl << " New Momentum Direction: " << newMomDir << G4endl
160 << " Momentum Change: " << *(aParticleChange.GetMomentumDirection())
161 << G4endl;
162 }
163
164 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
165}
int G4int
Definition G4Types.hh:85
#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
G4double GetConstProperty(const G4String &key) 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 G4OpMieHG::PreparePhysicsTable ( const G4ParticleDefinition & )
overridevirtual

Reimplemented from G4VProcess.

Definition at line 62 of file G4OpMieHG.cc.

63{
64 Initialise();
65}

◆ SetVerboseLevel()

void G4OpMieHG::SetVerboseLevel ( G4int verbose)

Definition at line 187 of file G4OpMieHG.cc.

Referenced by Initialise().


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