Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4OpMieHG.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26//
27////////////////////////////////////////////////////////////////////////
28//
29// File G4OpMieHG.hh
30// Description: Discrete Process -- Mie Scattering of Optical Photons
31// Created: 2010-07-03
32// Author: Xin Qian
33// Based on work from Vlasios Vasileiou
34//
35// This subroutine will mimic the Mie scattering based on
36// Henyey-Greenstein phase function
37// Forward and backward angles are treated separately.
38//
39////////////////////////////////////////////////////////////////////////
40
41#include "G4OpMieHG.hh"
44#include "G4OpProcessSubType.hh"
45
46//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
48 : G4VDiscreteProcess(processName, type)
49{
50 Initialise();
51 if(verboseLevel > 0)
52 {
53 G4cout << GetProcessName() << " is created " << G4endl;
54 }
56}
57
58//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
59G4OpMieHG::~G4OpMieHG() = default;
60
61//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
66
67//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
69{
70 SetVerboseLevel(G4OpticalParameters::Instance()->GetMieVerboseLevel());
71}
72
73//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
75 const G4Step& aStep)
76{
77 aParticleChange.Initialize(aTrack);
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 G4bool isForward = (G4UniformRand() <= forwardRatio);
95 G4double gg =
97
98 // Sample the direction
100 G4double costh;
101
102 if(gg != 0.)
103 {
104 G4double fact = (1. + gg)/(1. - gg + 2. * gg * r);
105 costh = 2. * r * fact * fact * (1. - gg + gg * r) - 1.;
106 }
107 else
108 {
109 costh = 2. * r - 1.;
110 }
111 G4double sinth = std::sqrt(std::max(0., (1. - costh) * ( 1. + costh)));
112
113 if(!isForward)
114 costh *= -1.0; // backward scattering
115
116 G4double phi = G4UniformRand() * twopi;
117 G4ThreeVector newMomDir(sinth * std::cos(phi), sinth * std::sin(phi), costh);
118 newMomDir.rotateUz(aParticle->GetMomentumDirection());
119 newMomDir = newMomDir.unit();
120
121 G4ThreeVector oldPol = aParticle->GetPolarization();
122 G4ThreeVector newPol = (newMomDir - oldPol / newMomDir.dot(oldPol)).unit();
123
124 if(newPol.mag() == 0.)
125 {
126 r = G4UniformRand() * twopi;
127 newPol.set(std::cos(r), std::sin(r), 0.);
128 newPol.rotateUz(newMomDir);
129 }
130 else
131 {
132 // There are two directions perpendicular to new momentum direction
133 if(G4UniformRand() < 0.5)
134 newPol = -newPol;
135 }
136
137 aParticleChange.ProposePolarization(newPol);
138 aParticleChange.ProposeMomentumDirection(newMomDir);
139
140 if(verboseLevel > 1)
141 {
142 G4cout << "OpMie New Polarization: " << newPol << G4endl
143 << " Polarization Change: " << *(aParticleChange.GetPolarization())
144 << G4endl << " New Momentum Direction: " << newMomDir << G4endl
145 << " Momentum Change: " << *(aParticleChange.GetMomentumDirection())
146 << G4endl;
147 }
148
149 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
150}
151
152//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
155{
156 G4double attLength = DBL_MAX;
159 if(MPT)
160 {
162 if(attVector)
163 {
164 attLength = attVector->Value(
165 aTrack.GetDynamicParticle()->GetTotalEnergy(), idx_mie);
166 }
167 }
168 return attLength;
169}
170
171//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
G4ForceCondition
G4PhysicsFreeVector G4MaterialPropertyVector
G4ProcessType
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
#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
G4double GetTotalEnergy() const
const G4ThreeVector & GetPolarization() const
G4double GetConstProperty(const G4String &key) const
G4MaterialPropertyVector * GetProperty(const char *key) const
G4MaterialPropertiesTable * GetMaterialPropertiesTable() const
virtual G4double GetMeanFreePath(const G4Track &aTrack, G4double, G4ForceCondition *) override
Definition G4OpMieHG.cc:153
G4OpMieHG(const G4String &processName="OpMieHG", G4ProcessType type=fOptical)
Definition G4OpMieHG.cc:47
virtual void PreparePhysicsTable(const G4ParticleDefinition &) override
Definition G4OpMieHG.cc:62
void SetVerboseLevel(G4int)
Definition G4OpMieHG.cc:172
virtual void Initialise()
Definition G4OpMieHG.cc:68
virtual G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep) override
Definition G4OpMieHG.cc:74
virtual ~G4OpMieHG()
static G4OpticalParameters * Instance()
G4double Value(const G4double energy, std::size_t &lastidx) const
G4Material * GetMaterial() const
const G4DynamicParticle * GetDynamicParticle() const
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
G4VDiscreteProcess(const G4String &aName, G4ProcessType aType=fNotDefined)
G4ParticleChange aParticleChange
G4int verboseLevel
void SetProcessSubType(G4int)
const G4String & GetProcessName() const
#define DBL_MAX
Definition templates.hh:62