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 | // $Id: G4KleinNishinaCompton.cc,v 1.10 2009/05/15 17:12:33 vnivanch Exp $ |
---|
27 | // GEANT4 tag $Name: geant4-09-03-beta-cand-01 $ |
---|
28 | // |
---|
29 | // ------------------------------------------------------------------- |
---|
30 | // |
---|
31 | // GEANT4 Class file |
---|
32 | // |
---|
33 | // |
---|
34 | // File name: G4KleinNishinaCompton |
---|
35 | // |
---|
36 | // Author: Vladimir Ivanchenko on base of Michel Maire code |
---|
37 | // |
---|
38 | // Creation date: 15.03.2005 |
---|
39 | // |
---|
40 | // Modifications: |
---|
41 | // 18-04-05 Use G4ParticleChangeForGamma (V.Ivantchenko) |
---|
42 | // 27-03-06 Remove upper limit of cross section (V.Ivantchenko) |
---|
43 | // |
---|
44 | // Class Description: |
---|
45 | // |
---|
46 | // ------------------------------------------------------------------- |
---|
47 | // |
---|
48 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
49 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
50 | |
---|
51 | #include "G4KleinNishinaCompton.hh" |
---|
52 | #include "G4Electron.hh" |
---|
53 | #include "G4Gamma.hh" |
---|
54 | #include "Randomize.hh" |
---|
55 | #include "G4DataVector.hh" |
---|
56 | #include "G4ParticleChangeForGamma.hh" |
---|
57 | |
---|
58 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
59 | |
---|
60 | using namespace std; |
---|
61 | |
---|
62 | G4KleinNishinaCompton::G4KleinNishinaCompton(const G4ParticleDefinition*, |
---|
63 | const G4String& nam) |
---|
64 | : G4VEmModel(nam) |
---|
65 | { |
---|
66 | theGamma = G4Gamma::Gamma(); |
---|
67 | theElectron = G4Electron::Electron(); |
---|
68 | lowestGammaEnergy = 1.0*eV; |
---|
69 | } |
---|
70 | |
---|
71 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
72 | |
---|
73 | G4KleinNishinaCompton::~G4KleinNishinaCompton() |
---|
74 | {} |
---|
75 | |
---|
76 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
77 | |
---|
78 | void G4KleinNishinaCompton::Initialise(const G4ParticleDefinition*, |
---|
79 | const G4DataVector&) |
---|
80 | { |
---|
81 | fParticleChange = GetParticleChangeForGamma(); |
---|
82 | } |
---|
83 | |
---|
84 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
85 | |
---|
86 | G4double G4KleinNishinaCompton::ComputeCrossSectionPerAtom( |
---|
87 | const G4ParticleDefinition*, |
---|
88 | G4double GammaEnergy, |
---|
89 | G4double Z, G4double, |
---|
90 | G4double, G4double) |
---|
91 | { |
---|
92 | G4double CrossSection = 0.0 ; |
---|
93 | if ( Z < 0.9999 ) return CrossSection; |
---|
94 | if ( GammaEnergy < 0.1*keV ) return CrossSection; |
---|
95 | // if ( GammaEnergy > (100.*GeV/Z) ) return CrossSection; |
---|
96 | |
---|
97 | static const G4double a = 20.0 , b = 230.0 , c = 440.0; |
---|
98 | |
---|
99 | static const G4double |
---|
100 | d1= 2.7965e-1*barn, d2=-1.8300e-1*barn, d3= 6.7527 *barn, d4=-1.9798e+1*barn, |
---|
101 | e1= 1.9756e-5*barn, e2=-1.0205e-2*barn, e3=-7.3913e-2*barn, e4= 2.7079e-2*barn, |
---|
102 | f1=-3.9178e-7*barn, f2= 6.8241e-5*barn, f3= 6.0480e-5*barn, f4= 3.0274e-4*barn; |
---|
103 | |
---|
104 | G4double p1Z = Z*(d1 + e1*Z + f1*Z*Z), p2Z = Z*(d2 + e2*Z + f2*Z*Z), |
---|
105 | p3Z = Z*(d3 + e3*Z + f3*Z*Z), p4Z = Z*(d4 + e4*Z + f4*Z*Z); |
---|
106 | |
---|
107 | G4double T0 = 15.0*keV; |
---|
108 | if (Z < 1.5) T0 = 40.0*keV; |
---|
109 | |
---|
110 | G4double X = max(GammaEnergy, T0) / electron_mass_c2; |
---|
111 | CrossSection = p1Z*std::log(1.+2.*X)/X |
---|
112 | + (p2Z + p3Z*X + p4Z*X*X)/(1. + a*X + b*X*X + c*X*X*X); |
---|
113 | |
---|
114 | // modification for low energy. (special case for Hydrogen) |
---|
115 | if (GammaEnergy < T0) { |
---|
116 | G4double dT0 = 1.*keV; |
---|
117 | X = (T0+dT0) / electron_mass_c2 ; |
---|
118 | G4double sigma = p1Z*log(1.+2*X)/X |
---|
119 | + (p2Z + p3Z*X + p4Z*X*X)/(1. + a*X + b*X*X + c*X*X*X); |
---|
120 | G4double c1 = -T0*(sigma-CrossSection)/(CrossSection*dT0); |
---|
121 | G4double c2 = 0.150; |
---|
122 | if (Z > 1.5) c2 = 0.375-0.0556*log(Z); |
---|
123 | G4double y = log(GammaEnergy/T0); |
---|
124 | CrossSection *= exp(-y*(c1+c2*y)); |
---|
125 | } |
---|
126 | // G4cout << "e= " << GammaEnergy << " Z= " << Z << " cross= " << CrossSection << G4endl; |
---|
127 | return CrossSection; |
---|
128 | } |
---|
129 | |
---|
130 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
131 | |
---|
132 | void G4KleinNishinaCompton::SampleSecondaries(std::vector<G4DynamicParticle*>* fvect, |
---|
133 | const G4MaterialCutsCouple*, |
---|
134 | const G4DynamicParticle* aDynamicGamma, |
---|
135 | G4double, |
---|
136 | G4double) |
---|
137 | { |
---|
138 | // The scattered gamma energy is sampled according to Klein - Nishina formula. |
---|
139 | // The random number techniques of Butcher & Messel are used |
---|
140 | // (Nuc Phys 20(1960),15). |
---|
141 | // Note : Effects due to binding of atomic electrons are negliged. |
---|
142 | |
---|
143 | G4double gamEnergy0 = aDynamicGamma->GetKineticEnergy(); |
---|
144 | G4double E0_m = gamEnergy0 / electron_mass_c2 ; |
---|
145 | |
---|
146 | G4ThreeVector gamDirection0 = aDynamicGamma->GetMomentumDirection(); |
---|
147 | |
---|
148 | // |
---|
149 | // sample the energy rate of the scattered gamma |
---|
150 | // |
---|
151 | |
---|
152 | G4double epsilon, epsilonsq, onecost, sint2, greject ; |
---|
153 | |
---|
154 | G4double epsilon0 = 1./(1. + 2.*E0_m); |
---|
155 | G4double epsilon0sq = epsilon0*epsilon0; |
---|
156 | G4double alpha1 = - log(epsilon0); |
---|
157 | G4double alpha2 = 0.5*(1.- epsilon0sq); |
---|
158 | |
---|
159 | do { |
---|
160 | if ( alpha1/(alpha1+alpha2) > G4UniformRand() ) { |
---|
161 | epsilon = exp(-alpha1*G4UniformRand()); // epsilon0**r |
---|
162 | epsilonsq = epsilon*epsilon; |
---|
163 | |
---|
164 | } else { |
---|
165 | epsilonsq = epsilon0sq + (1.- epsilon0sq)*G4UniformRand(); |
---|
166 | epsilon = sqrt(epsilonsq); |
---|
167 | }; |
---|
168 | |
---|
169 | onecost = (1.- epsilon)/(epsilon*E0_m); |
---|
170 | sint2 = onecost*(2.-onecost); |
---|
171 | greject = 1. - epsilon*sint2/(1.+ epsilonsq); |
---|
172 | |
---|
173 | } while (greject < G4UniformRand()); |
---|
174 | |
---|
175 | // |
---|
176 | // scattered gamma angles. ( Z - axis along the parent gamma) |
---|
177 | // |
---|
178 | |
---|
179 | G4double cosTeta = 1. - onecost; |
---|
180 | G4double sinTeta = sqrt (sint2); |
---|
181 | G4double Phi = twopi * G4UniformRand(); |
---|
182 | G4double dirx = sinTeta*cos(Phi), diry = sinTeta*sin(Phi), dirz = cosTeta; |
---|
183 | |
---|
184 | // |
---|
185 | // update G4VParticleChange for the scattered gamma |
---|
186 | // |
---|
187 | |
---|
188 | G4ThreeVector gamDirection1 ( dirx,diry,dirz ); |
---|
189 | gamDirection1.rotateUz(gamDirection0); |
---|
190 | G4double gamEnergy1 = epsilon*gamEnergy0; |
---|
191 | fParticleChange->SetProposedKineticEnergy(gamEnergy1); |
---|
192 | if(gamEnergy1 > lowestGammaEnergy) { |
---|
193 | fParticleChange->ProposeMomentumDirection(gamDirection1); |
---|
194 | } else { |
---|
195 | fParticleChange->ProposeTrackStatus(fStopAndKill); |
---|
196 | gamEnergy1 += fParticleChange->GetLocalEnergyDeposit(); |
---|
197 | fParticleChange->ProposeLocalEnergyDeposit(gamEnergy1); |
---|
198 | } |
---|
199 | |
---|
200 | // |
---|
201 | // kinematic of the scattered electron |
---|
202 | // |
---|
203 | |
---|
204 | G4double eKinEnergy = gamEnergy0 - gamEnergy1; |
---|
205 | |
---|
206 | if(eKinEnergy > DBL_MIN) { |
---|
207 | G4ThreeVector eDirection = gamEnergy0*gamDirection0 - gamEnergy1*gamDirection1; |
---|
208 | eDirection = eDirection.unit(); |
---|
209 | |
---|
210 | // create G4DynamicParticle object for the electron. |
---|
211 | G4DynamicParticle* dp = new G4DynamicParticle(theElectron,eDirection,eKinEnergy); |
---|
212 | fvect->push_back(dp); |
---|
213 | } |
---|
214 | } |
---|
215 | |
---|
216 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
217 | |
---|
218 | |
---|