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
G4XrayRayleighModel.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// Author: Vladimir Grichine
29//
30// History:
31//
32// 14.10.12 V.Grichine, update of xsc and angular distribution
33// 25.05.2011 first implementation
34
37#include "G4SystemOfUnits.hh"
38
39//////////////////////////////////////////////////////////////////////////////////
40
41const G4double G4XrayRayleighModel::fCofA = 2.*pi2*Bohr_radius*Bohr_radius;
42
43const G4double G4XrayRayleighModel::fCofR = 8.*pi*classic_electr_radius*classic_electr_radius/3.;
44
45//////////////////////////////////////////////////////////////////////////////////
46
48 const G4String& nam)
49 :G4VEmModel(nam),isInitialised(false)
50{
51 fParticleChange = nullptr;
52 lowEnergyLimit = 250*eV;
53 highEnergyLimit = 10.*MeV;
54 fFormFactor = 0.0;
55
56 // SetLowEnergyLimit(lowEnergyLimit);
57 SetHighEnergyLimit(highEnergyLimit);
58 //
59 verboseLevel= 0;
60 // Verbosity scale:
61 // 0 = nothing
62 // 1 = warning for energy non-conservation
63 // 2 = details of energy budget
64 // 3 = calculation of cross sections, file openings, sampling of atoms
65 // 4 = entering in methods
66
67 if(verboseLevel > 0)
68 {
69 G4cout << "Xray Rayleigh is constructed " << G4endl
70 << "Energy range: "
71 << lowEnergyLimit / eV << " eV - "
72 << highEnergyLimit / MeV << " MeV"
73 << G4endl;
74 }
75}
76
77//////////////////////////////////////////////////////////////////////////////////
78
80
81//////////////////////////////////////////////////////////////////////////////////
82
84 const G4DataVector& cuts)
85{
86 if (verboseLevel > 3)
87 {
88 G4cout << "Calling G4XrayRayleighModel::Initialise()" << G4endl;
89 }
90
91 InitialiseElementSelectors(particle,cuts);
92
93
94 if(isInitialised) return;
95 fParticleChange = GetParticleChangeForGamma();
96 isInitialised = true;
97
98}
99
100//////////////////////////////////////////////////////////////////////////////////
101
104 G4double gammaEnergy,
107{
108 if (verboseLevel > 3)
109 {
110 G4cout << "Calling CrossSectionPerAtom() of G4XrayRayleighModel" << G4endl;
111 }
112 if (gammaEnergy < lowEnergyLimit || gammaEnergy > highEnergyLimit)
113 {
114 return 0.0;
115 }
116 G4double k = gammaEnergy/hbarc;
117 k *= Bohr_radius;
118 G4double p0 = 0.680654;
119 G4double p1 = -0.0224188;
120 G4double lnZ = std::log(Z);
121
122 G4double lna = p0 + p1*lnZ;
123
124 G4double alpha = std::exp(lna);
125
126 G4double fo = std::pow(k, alpha);
127
128 p0 = 3.68455;
129 p1 = -0.464806;
130 lna = p0 + p1*lnZ;
131
132 fo *= 0.01*std::exp(lna);
133
134 fFormFactor = fo;
135
136 G4double b = 1. + 2.*fo;
137 G4double b2 = b*b;
138 G4double b3 = b*b2;
139
140 G4double xsc = fCofR*Z*Z/b3;
141 xsc *= fo*fo + (1. + fo)*(1. + fo);
142
143
144 return xsc;
145
146}
147
148//////////////////////////////////////////////////////////////////////////////////
149
150void G4XrayRayleighModel::SampleSecondaries(std::vector<G4DynamicParticle*>* /*fvect*/,
151 const G4MaterialCutsCouple* couple,
152 const G4DynamicParticle* aDPGamma,
153 G4double,
154 G4double)
155{
156 if ( verboseLevel > 3)
157 {
158 G4cout << "Calling SampleSecondaries() of G4XrayRayleighModel" << G4endl;
159 }
160 G4double photonEnergy0 = aDPGamma->GetKineticEnergy();
161
162 G4ParticleMomentum photonDirection0 = aDPGamma->GetMomentumDirection();
163
164
165 // Sample the angle of the scattered photon
166 // according to 1 + cosTheta*cosTheta distribution
167
168 G4double cosDipole, cosTheta, sinTheta;
169 G4double c, delta, cofA, signc = 1., a, power = 1./3.;
170
171 c = 4. - 8.*G4UniformRand();
172 a = c;
173
174 if( c < 0. )
175 {
176 signc = -1.;
177 a = -c;
178 }
179 delta = std::sqrt(a*a+4.);
180 delta += a;
181 delta *= 0.5;
182 cofA = -signc*std::pow(delta, power);
183 cosDipole = cofA - 1./cofA;
184
185 // select atom
186 const G4Element* elm = SelectTargetAtom(couple, aDPGamma->GetParticleDefinition(),
187 photonEnergy0,aDPGamma->GetLogKineticEnergy());
188 G4double Z = elm->GetZ();
189
190 G4double k = photonEnergy0/hbarc;
191 k *= Bohr_radius;
192 G4double p0 = 0.680654;
193 G4double p1 = -0.0224188;
194 G4double lnZ = std::log(Z);
195
196 G4double lna = p0 + p1*lnZ;
197
198 G4double alpha = std::exp(lna);
199
200 G4double fo = std::pow(k, alpha);
201
202 p0 = 3.68455;
203 p1 = -0.464806;
204 lna = p0 + p1*lnZ;
205
206 fo *= 0.01*pi*std::exp(lna);
207
208
209 G4double beta = fo/(1 + fo);
210
211 cosTheta = (cosDipole + beta)/(1. + cosDipole*beta);
212
213
214 if( cosTheta > 1.) cosTheta = 1.;
215 if( cosTheta < -1.) cosTheta = -1.;
216
217 sinTheta = std::sqrt( (1. - cosTheta)*(1. + cosTheta) );
218
219 // Scattered photon angles. ( Z - axis along the parent photon)
220
221 G4double phi = twopi * G4UniformRand() ;
222 G4double dirX = sinTheta*std::cos(phi);
223 G4double dirY = sinTheta*std::sin(phi);
224 G4double dirZ = cosTheta;
225
226 // Update G4VParticleChange for the scattered photon
227
228 G4ThreeVector photonDirection1(dirX, dirY, dirZ);
229 photonDirection1.rotateUz(photonDirection0);
230 fParticleChange->ProposeMomentumDirection(photonDirection1);
231
232 fParticleChange->SetProposedKineticEnergy(photonEnergy0);
233}
234
235
double G4double
Definition: G4Types.hh:83
const G4int Z[17]
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:52
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:33
const G4ThreeVector & GetMomentumDirection() const
G4double GetLogKineticEnergy() const
const G4ParticleDefinition * GetParticleDefinition() const
G4double GetKineticEnergy() const
G4double GetZ() const
Definition: G4Element.hh:131
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeMomentumDirection(const G4ThreeVector &Pfinal)
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:746
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:124
const G4Element * SelectTargetAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double logKineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:577
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:139
~G4XrayRayleighModel() override
G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0, G4double cut=0, G4double emax=DBL_MAX) override
G4XrayRayleighModel(const G4ParticleDefinition *p=nullptr, const G4String &nam="XrayRayleigh")
void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
const G4double pi