Geant4 9.6.0
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)
 
 ~G4OpMieHG ()
 
G4bool IsApplicable (const G4ParticleDefinition &aParticleType)
 
G4double GetMeanFreePath (const G4Track &aTrack, G4double, G4ForceCondition *)
 
G4VParticleChangePostStepDoIt (const G4Track &aTrack, const G4Step &aStep)
 
- 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
 

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

Detailed Description

Definition at line 48 of file G4OpMieHG.hh.

Constructor & Destructor Documentation

◆ G4OpMieHG()

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

Definition at line 47 of file G4OpMieHG.cc.

48 : G4VDiscreteProcess(processName, type)
49{
50 if (verboseLevel>0) {
51 G4cout << GetProcessName() << " is created " << G4endl;
52 }
53
55}
@ fOpMieHG
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
G4int verboseLevel
Definition: G4VProcess.hh:368
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:403
const G4String & GetProcessName() const
Definition: G4VProcess.hh:379

◆ ~G4OpMieHG()

G4OpMieHG::~G4OpMieHG ( )

Definition at line 57 of file G4OpMieHG.cc.

57{}

Member Function Documentation

◆ GetMeanFreePath()

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

Implements G4VDiscreteProcess.

Definition at line 159 of file G4OpMieHG.cc.

162{
163 const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
164 const G4Material* aMaterial = aTrack.GetMaterial();
165
166 G4double thePhotonEnergy = aParticle->GetTotalEnergy();
167
168 G4double AttenuationLength = DBL_MAX;
169
170 G4MaterialPropertiesTable* aMaterialPropertyTable =
171 aMaterial->GetMaterialPropertiesTable();
172
173 if (aMaterialPropertyTable) {
174 G4MaterialPropertyVector* AttenuationLengthVector =
175 aMaterialPropertyTable->GetProperty("MIEHG");
176 if (AttenuationLengthVector) {
177 AttenuationLength = AttenuationLengthVector ->
178 Value(thePhotonEnergy);
179 } else {
180// G4cout << "No Mie scattering length specified" << G4endl;
181 }
182 } else {
183// G4cout << "No Mie scattering length specified" << G4endl;
184 }
185
186// G4cout << thePhotonEnergy/GeV << " \t" << AttenuationLength/m << G4endl;
187
188 return AttenuationLength;
189}
double G4double
Definition: G4Types.hh:64
G4double GetTotalEnergy() const
G4MaterialPropertyVector * GetProperty(const char *key)
G4MaterialPropertiesTable * GetMaterialPropertiesTable() const
Definition: G4Material.hh:251
G4Material * GetMaterial() const
const G4DynamicParticle * GetDynamicParticle() const
#define DBL_MAX
Definition: templates.hh:83

◆ IsApplicable()

G4bool G4OpMieHG::IsApplicable ( const G4ParticleDefinition aParticleType)
inlinevirtual

Reimplemented from G4VProcess.

Definition at line 92 of file G4OpMieHG.hh.

93{
94 return ( &aParticleType == G4OpticalPhoton::OpticalPhoton() );
95}
static G4OpticalPhoton * OpticalPhoton()

◆ PostStepDoIt()

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

Reimplemented from G4VDiscreteProcess.

Definition at line 67 of file G4OpMieHG.cc.

68{
70
71 const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
72 const G4Material* aMaterial = aTrack.GetMaterial();
73 G4MaterialPropertiesTable* aMaterialPropertyTable =
74 aMaterial->GetMaterialPropertiesTable();
75
76 G4double forward_g =
77 aMaterialPropertyTable->GetConstProperty("MIEHG_FORWARD");
78 G4double backward_g =
79 aMaterialPropertyTable->GetConstProperty("MIEHG_BACKWARD");
80 G4double ForwardRatio =
81 aMaterialPropertyTable->GetConstProperty("MIEHG_FORWARD_RATIO");
82
83 if (verboseLevel>0) {
84 G4cout << "MIE Scattering Photon!" << G4endl;
85 G4cout << "MIE Old Momentum Direction: "
86 << aParticle->GetMomentumDirection() << G4endl;
87 G4cout << "MIE Old Polarization: "
88 << aParticle->GetPolarization() << G4endl;
89 }
90
91 G4double gg;
92 G4int direction;
93 if (G4UniformRand()<=ForwardRatio){
94 gg = forward_g;
95 direction = 1;
96 } else {
97 gg = backward_g;
98 direction = -1;
99 }
100
102
103 G4double Theta;
104 //sample the direction
105 if (gg!=0) {
106 Theta = std::acos(2*r*(1+gg)*(1+gg)*(1-gg+gg*r)/((1-gg+2*gg*r)*(1-gg+2*gg*r)) -1);
107 } else {
108 Theta = std::acos(2*r-1.);
109 }
110 G4double Phi = G4UniformRand()*2*pi;
111
112 if (direction==-1) Theta = pi - Theta; //backward scattering
113
114 G4ThreeVector NewMomentumDirection, OldMomentumDirection;
115 G4ThreeVector OldPolarization, NewPolarization;
116
117 NewMomentumDirection.set
118 (std::sin(Theta)*std::cos(Phi), std::sin(Theta)*std::sin(Phi), std::cos(Theta));
119 OldMomentumDirection = aParticle->GetMomentumDirection();
120 NewMomentumDirection.rotateUz(OldMomentumDirection);
121 NewMomentumDirection = NewMomentumDirection.unit();
122
123 OldPolarization = aParticle->GetPolarization();
124 G4double constant = -1./NewMomentumDirection.dot(OldPolarization);
125
126 NewPolarization = NewMomentumDirection + constant*OldPolarization;
127 NewPolarization = NewPolarization.unit();
128
129 if (NewPolarization.mag()==0) {
130 r = G4UniformRand()*twopi;
131 NewPolarization.set(std::cos(r),std::sin(r),0.);
132 NewPolarization.rotateUz(NewMomentumDirection);
133 } else {
134 // There are two directions which perpendicular
135 // new momentum direction
136 if (G4UniformRand() < 0.5) NewPolarization = -NewPolarization;
137 }
138
139 aParticleChange.ProposePolarization(NewPolarization);
140 aParticleChange.ProposeMomentumDirection(NewMomentumDirection);
141
142 if (verboseLevel>0) {
143 G4cout << "MIE New Polarization: "
144 << NewPolarization << G4endl;
145 G4cout << "MIE Polarization Change: "
147 G4cout << "MIE New Momentum Direction: "
148 << NewMomentumDirection << G4endl;
149 G4cout << "MIE Momentum Change: "
151 }
152
153 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
154}
int G4int
Definition: G4Types.hh:66
#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
G4double GetConstProperty(const char *key)
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
const G4double pi

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