Geant4-11
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
40
42
44
46
48 const G4String& nam)
49 :G4VEmModel(nam),isInitialised(false)
50{
52 lowEnergyLimit = 250*eV;
53 highEnergyLimit = 10.*MeV;
54 fFormFactor = 0.0;
55
56 // SetLowEnergyLimit(lowEnergyLimit);
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
78
80{
81
82}
83
85
87 const G4DataVector& cuts)
88{
89 if (verboseLevel > 3)
90 {
91 G4cout << "Calling G4XrayRayleighModel::Initialise()" << G4endl;
92 }
93
94 InitialiseElementSelectors(particle,cuts);
95
96
97 if(isInitialised) return;
99 isInitialised = true;
100
101}
102
104
107 G4double gammaEnergy,
110{
111 if (verboseLevel > 3)
112 {
113 G4cout << "Calling CrossSectionPerAtom() of G4XrayRayleighModel" << G4endl;
114 }
115 if (gammaEnergy < lowEnergyLimit || gammaEnergy > highEnergyLimit)
116 {
117 return 0.0;
118 }
119 G4double k = gammaEnergy/hbarc;
120 k *= Bohr_radius;
121 G4double p0 = 0.680654;
122 G4double p1 = -0.0224188;
123 G4double lnZ = std::log(Z);
124
125 G4double lna = p0 + p1*lnZ;
126
127 G4double alpha = std::exp(lna);
128
129 G4double fo = std::pow(k, alpha);
130
131 p0 = 3.68455;
132 p1 = -0.464806;
133 lna = p0 + p1*lnZ;
134
135 fo *= 0.01*std::exp(lna);
136
137 fFormFactor = fo;
138
139 G4double b = 1. + 2.*fo;
140 G4double b2 = b*b;
141 G4double b3 = b*b2;
142
143 G4double xsc = fCofR*Z*Z/b3;
144 xsc *= fo*fo + (1. + fo)*(1. + fo);
145
146
147 return xsc;
148
149}
150
152
153void G4XrayRayleighModel::SampleSecondaries(std::vector<G4DynamicParticle*>* /*fvect*/,
154 const G4MaterialCutsCouple* couple,
155 const G4DynamicParticle* aDPGamma,
156 G4double,
157 G4double)
158{
159 if ( verboseLevel > 3)
160 {
161 G4cout << "Calling SampleSecondaries() of G4XrayRayleighModel" << G4endl;
162 }
163 G4double photonEnergy0 = aDPGamma->GetKineticEnergy();
164
165 G4ParticleMomentum photonDirection0 = aDPGamma->GetMomentumDirection();
166
167
168 // Sample the angle of the scattered photon
169 // according to 1 + cosTheta*cosTheta distribution
170
171 G4double cosDipole, cosTheta, sinTheta;
172 G4double c, delta, cofA, signc = 1., a, power = 1./3.;
173
174 c = 4. - 8.*G4UniformRand();
175 a = c;
176
177 if( c < 0. )
178 {
179 signc = -1.;
180 a = -c;
181 }
182 delta = std::sqrt(a*a+4.);
183 delta += a;
184 delta *= 0.5;
185 cofA = -signc*std::pow(delta, power);
186 cosDipole = cofA - 1./cofA;
187
188 // select atom
189 const G4Element* elm = SelectTargetAtom(couple, aDPGamma->GetParticleDefinition(),
190 photonEnergy0,aDPGamma->GetLogKineticEnergy());
191 G4double Z = elm->GetZ();
192
193 G4double k = photonEnergy0/hbarc;
194 k *= Bohr_radius;
195 G4double p0 = 0.680654;
196 G4double p1 = -0.0224188;
197 G4double lnZ = std::log(Z);
198
199 G4double lna = p0 + p1*lnZ;
200
201 G4double alpha = std::exp(lna);
202
203 G4double fo = std::pow(k, alpha);
204
205 p0 = 3.68455;
206 p1 = -0.464806;
207 lna = p0 + p1*lnZ;
208
209 fo *= 0.01*pi*std::exp(lna);
210
211
212 G4double beta = fo/(1 + fo);
213
214 cosTheta = (cosDipole + beta)/(1. + cosDipole*beta);
215
216
217 if( cosTheta > 1.) cosTheta = 1.;
218 if( cosTheta < -1.) cosTheta = -1.;
219
220 sinTheta = std::sqrt( (1. - cosTheta)*(1. + cosTheta) );
221
222 // Scattered photon angles. ( Z - axis along the parent photon)
223
224 G4double phi = twopi * G4UniformRand() ;
225 G4double dirX = sinTheta*std::cos(phi);
226 G4double dirY = sinTheta*std::sin(phi);
227 G4double dirZ = cosTheta;
228
229 // Update G4VParticleChange for the scattered photon
230
231 G4ThreeVector photonDirection1(dirX, dirY, dirZ);
232 photonDirection1.rotateUz(photonDirection0);
233 fParticleChange->ProposeMomentumDirection(photonDirection1);
234
236}
237
238
static const G4double alpha
static constexpr double pi2
Definition: G4SIunits.hh:58
static constexpr double twopi
Definition: G4SIunits.hh:56
static constexpr double eV
Definition: G4SIunits.hh:201
static constexpr double MeV
Definition: G4SIunits.hh:200
static constexpr double pi
Definition: G4SIunits.hh:55
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(G4double Px, G4double Py, G4double Pz)
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:767
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:123
const G4Element * SelectTargetAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double logKineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:597
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:138
static const G4double fCofR
G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0, G4double cut=0, G4double emax=DBL_MAX) override
G4XrayRayleighModel(const G4ParticleDefinition *p=0, 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
static const G4double fCofA
G4ParticleChangeForGamma * fParticleChange
float Bohr_radius
Definition: hepunit.py:289
int classic_electr_radius
Definition: hepunit.py:287
float hbarc
Definition: hepunit.py:264