Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4OpRayleigh.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////////////////////////////////////////////////////////////////////////
30// Optical Photon Rayleigh Scattering Class Implementation
31////////////////////////////////////////////////////////////////////////
32//
33// File: G4OpRayleigh.cc
34// Description: Discrete Process -- Rayleigh scattering of optical
35// photons
36// Version: 1.0
37// Created: 1996-05-31
38// Author: Juliet Armstrong
39// Updated: 2014-10-10 - This version calculates the Rayleigh scattering
40// length for more materials than just Water (although the Water
41// default is kept). To do this the user would need to specify the
42// ISOTHERMAL_COMPRESSIBILITY as a material property and
43// optionally an RS_SCALE_LENGTH (useful for testing). Code comes
44// from Philip Graham (Queen Mary University of London).
45// 2010-06-11 - Fix Bug 207; Thanks to Xin Qian
46// (Kellogg Radiation Lab of Caltech)
47// 2005-07-28 - add G4ProcessType to constructor
48// 2001-10-18 by Peter Gumplinger
49// eliminate unused variable warning on Linux (gcc-2.95.2)
50// 2001-09-18 by mma
51// >numOfMaterials=G4Material::GetNumberOfMaterials() in BuildPhy
52// 2001-01-30 by Peter Gumplinger
53// > allow for positiv and negative CosTheta and force the
54// > new momentum direction to be in the same plane as the
55// > new and old polarization vectors
56// 2001-01-29 by Peter Gumplinger
57// > fix calculation of SinTheta (from CosTheta)
58// 1997-04-09 by Peter Gumplinger
59// > new physics/tracking scheme
60//
61////////////////////////////////////////////////////////////////////////
62
63#include "G4OpRayleigh.hh"
64#include "G4ios.hh"
66#include "G4SystemOfUnits.hh"
68#include "G4OpProcessSubType.hh"
69
70//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
72 : G4VDiscreteProcess(processName, type)
73{
74 Initialise();
76 thePhysicsTable = nullptr;
77
78 if(verboseLevel > 0)
79 {
80 G4cout << GetProcessName() << " is created " << G4endl;
81 }
82}
83
84//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
86{
87 // VI: inside this PhysicsTable all properties are unique
88 // it is not possible to destroy
89 delete thePhysicsTable;
90}
91
92//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
94{
95 Initialise();
96}
97
98//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
100{
101 SetVerboseLevel(G4OpticalParameters::Instance()->GetRayleighVerboseLevel());
102}
103
104//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
106 const G4Step& aStep)
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}
186
187//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
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}
215
216//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
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}
231
232//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
233G4PhysicsFreeVector* G4OpRayleigh::CalculateRayleighMeanFreePaths(
234 const G4Material* material) const
235{
237
238 // Retrieve the beta_T or isothermal compressibility value. For backwards
239 // compatibility use a constant if the material is "Water". If the material
240 // doesn't have an ISOTHERMAL_COMPRESSIBILITY constant then return
241 G4double betat;
242 if(material->GetName() == "Water")
243 {
244 betat = 7.658e-23 * m3 / MeV;
245 }
247 {
249 }
250 else
251 {
252 return nullptr;
253 }
254
255 // If the material doesn't have a RINDEX property vector then return
257 if(rIndex == nullptr)
258 return nullptr;
259
260 // Retrieve the optional scale factor (scales the scattering length)
261 G4double scaleFactor = 1.0;
263 {
264 scaleFactor = MPT->GetConstProperty(kRS_SCALE_FACTOR);
265 }
266
267 // Retrieve the material temperature. For backwards compatibility use a
268 // constant if the material is "Water"
269 G4double temperature;
270 if(material->GetName() == "Water")
271 {
272 temperature =
273 283.15 * kelvin; // Temperature of water is 10 degrees celsius
274 }
275 else
276 {
277 temperature = material->GetTemperature();
278 }
279
280 auto rayleighMFPs = new G4PhysicsFreeVector();
281 // This calculates the meanFreePath via the Einstein-Smoluchowski formula
282 const G4double c1 =
283 scaleFactor * betat * temperature * k_Boltzmann / (6.0 * pi);
284
285 for(size_t uRIndex = 0; uRIndex < rIndex->GetVectorLength(); ++uRIndex)
286 {
287 const G4double energy = rIndex->Energy(uRIndex);
288 const G4double rIndexSquared = (*rIndex)[uRIndex] * (*rIndex)[uRIndex];
289 const G4double xlambda = h_Planck * c_light / energy;
290 const G4double c2 = std::pow(twopi / xlambda, 4);
291 const G4double c3 =
292 std::pow(((rIndexSquared - 1.0) * (rIndexSquared + 2.0) / 3.0), 2);
293
294 const G4double meanFreePath = 1.0 / (c1 * c2 * c3);
295
296 if(verboseLevel > 0)
297 {
298 G4cout << energy << "MeV\t" << meanFreePath << "mm" << G4endl;
299 }
300
301 rayleighMFPs->InsertValues(energy, meanFreePath);
302 }
303
304 return rayleighMFPs;
305}
306
307//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
309{
310 verboseLevel = verbose;
312}
G4ForceCondition
@ kISOTHERMAL_COMPRESSIBILITY
std::vector< G4Material * > G4MaterialTable
@ fOpRayleigh
G4ProcessType
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
#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
G4double GetTotalMomentum() const
const G4ThreeVector & GetPolarization() const
G4bool ConstPropertyExists(const G4String &key) const
G4double GetConstProperty(const G4String &key) const
G4MaterialPropertyVector * GetProperty(const char *key) const
G4MaterialPropertiesTable * GetMaterialPropertiesTable() const
Definition: G4Material.hh:251
static size_t GetNumberOfMaterials()
Definition: G4Material.cc:684
G4double GetTemperature() const
Definition: G4Material.hh:177
static G4MaterialTable * GetMaterialTable()
Definition: G4Material.cc:677
const G4String & GetName() const
Definition: G4Material.hh:172
size_t GetIndex() const
Definition: G4Material.hh:255
void SetVerboseLevel(G4int)
G4PhysicsTable * thePhysicsTable
Definition: G4OpRayleigh.hh:88
virtual void BuildPhysicsTable(const G4ParticleDefinition &aParticleType) override
virtual G4double GetMeanFreePath(const G4Track &aTrack, G4double, G4ForceCondition *) override
virtual ~G4OpRayleigh()
Definition: G4OpRayleigh.cc:85
virtual void PreparePhysicsTable(const G4ParticleDefinition &) override
Definition: G4OpRayleigh.cc:93
G4OpRayleigh(const G4String &processName="OpRayleigh", G4ProcessType type=fOptical)
Definition: G4OpRayleigh.cc:71
virtual G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep) override
virtual void Initialise()
Definition: G4OpRayleigh.cc:99
void SetRayleighVerboseLevel(G4int)
static G4OpticalParameters * Instance()
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)
void insertAt(std::size_t, G4PhysicsVector *)
G4double Energy(const std::size_t index) const
std::size_t GetVectorLength() const
Definition: G4Step.hh:62
G4Material * GetMaterial() const
const G4DynamicParticle * GetDynamicParticle() const
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:331
G4int verboseLevel
Definition: G4VProcess.hh:360
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:410
const G4String & GetProcessName() const
Definition: G4VProcess.hh:386
G4double energy(const ThreeVector &p, const G4double m)
const G4double pi
#define DBL_MAX
Definition: templates.hh:62