Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4WentzelOKandVIxSection.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// $Id$
27//
28// -------------------------------------------------------------------
29//
30// GEANT4 Class file
31//
32//
33// File name: G4WentzelOKandVIxSection
34//
35// Author: V.Ivanchenko
36//
37// Creation date: 09.04.2008 from G4MuMscModel
38//
39// Modifications:
40//
41//
42
43// -------------------------------------------------------------------
44//
45
46//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
47//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
48
51#include "G4SystemOfUnits.hh"
52#include "Randomize.hh"
53#include "G4Electron.hh"
54#include "G4Positron.hh"
55#include "G4Proton.hh"
56#include "G4LossTableManager.hh"
57
58//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
59
60G4double G4WentzelOKandVIxSection::ScreenRSquare[] = {0.0};
61G4double G4WentzelOKandVIxSection::FormFactor[] = {0.0};
62
63using namespace std;
64
66 numlimit(0.1),
67 nwarnings(0),
68 nwarnlimit(50),
69 alpha2(fine_structure_const*fine_structure_const)
70{
71 fNistManager = G4NistManager::Instance();
72 fG4pow = G4Pow::GetInstance();
73 theElectron = G4Electron::Electron();
74 thePositron = G4Positron::Positron();
75 theProton = G4Proton::Proton();
76 lowEnergyLimit = 1.0*eV;
77 G4double p0 = electron_mass_c2*classic_electr_radius;
78 coeff = twopi*p0*p0;
79 particle = 0;
80
81 // Thomas-Fermi screening radii
82 // Formfactors from A.V. Butkevich et al., NIM A 488 (2002) 282
83
84 if(0.0 == ScreenRSquare[0]) {
85 G4double a0 = electron_mass_c2/0.88534;
86 G4double constn = 6.937e-6/(MeV*MeV);
87
88 ScreenRSquare[0] = alpha2*a0*a0;
89 for(G4int j=1; j<100; ++j) {
90 G4double x = a0*fG4pow->Z13(j);
91 if(1 == j) { ScreenRSquare[j] = 0.5*alpha2*a0*a0; }
92 else {
93 ScreenRSquare[j] = 0.5*(1 + exp(-j*j*0.001))*alpha2*x*x;
94 //ScreenRSquare[j] = (0.5 + 0.5/G4double(j))*alpha2*x*x;
95 //ScreenRSquare[j] = 0.5*alpha2*x*x;
96 }
97 x = fNistManager->GetA27(j);
98 FormFactor[j] = constn*x*x;
99 }
100 }
101 currentMaterial = 0;
102 elecXSRatio = factB = factD = formfactA = screenZ = 0.0;
103 cosTetMaxElec = cosTetMaxNuc = invbeta2 = kinFactor = gam0pcmp = pcmp2 = 1.0;
104
105 factB1= 0.5*CLHEP::pi*fine_structure_const;
106
107 Initialise(theElectron, 1.0);
108 targetMass = proton_mass_c2;
109}
110
111//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
112
114{}
115
116//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
117
119 G4double CosThetaLim)
120{
121 SetupParticle(p);
122 tkin = mom2 = momCM2 = 0.0;
123 ecut = etag = DBL_MAX;
124 targetZ = 0;
125 cosThetaMax = CosThetaLim;
126 G4double a =
127 G4LossTableManager::Instance()->FactorForAngleLimit()*CLHEP::hbarc/CLHEP::fermi;
128 factorA2 = 0.5*a*a;
129 currentMaterial = 0;
130}
131
132//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
133
135{
136 particle = p;
137 mass = particle->GetPDGMass();
138 spin = particle->GetPDGSpin();
139 if(0.0 != spin) { spin = 0.5; }
140 G4double q = std::fabs(particle->GetPDGCharge()/eplus);
141 chargeSquare = q*q;
142 charge3 = chargeSquare*q;
143 tkin = 0.0;
144 currentMaterial = 0;
145 targetZ = 0;
146}
147
148//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
149
152{
153 G4double cosTetMaxNuc2 = cosTetMaxNuc;
154 if(Z != targetZ || tkin != etag) {
155 etag = tkin;
156 targetZ = Z;
157 if(targetZ > 99) { targetZ = 99; }
158 SetTargetMass(fNistManager->GetAtomicMassAmu(targetZ)*CLHEP::amu_c2);
159 //G4double tmass2 = targetMass*targetMass;
160 //G4double etot = tkin + mass;
161 //G4double invmass2 = mass*mass + tmass2 + 2*etot*targetMass;
162 //momCM2 = mom2*tmass2/invmass2;
163 //gam0pcmp = (etot + targetMass)*targetMass/invmass2;
164 //pcmp2 = tmass2/invmass2;
165
166 kinFactor = coeff*Z*chargeSquare*invbeta2/mom2;
167
168 screenZ = ScreenRSquare[targetZ]/mom2;
169 if(Z > 1) {
170 if(mass > MeV) {
171 screenZ *= std::min(Z*1.13,1.13 +3.76*Z*Z*invbeta2*alpha2*chargeSquare);
172 } else {
173 G4double tau = tkin/mass;
174 screenZ *= std::min(Z*1.13,(1.13 +3.76*Z*Z
175 *invbeta2*alpha2*std::sqrt(tau/(tau + fG4pow->Z23(targetZ)))));
176 }
177 }
178 if(targetZ == 1 && cosTetMaxNuc2 < 0.0 && particle == theProton) {
179 cosTetMaxNuc2 = 0.0;
180 }
181 formfactA = FormFactor[targetZ]*mom2;
182
183 cosTetMaxElec = 1.0;
184 ComputeMaxElectronScattering(cut);
185 }
186 return cosTetMaxNuc2;
187}
188
189//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
190
193{
194 G4double xsec = 0.0;
195 if(cosTMax >= 1.0) { return xsec; }
196
197 G4double xSection = 0.0;
198 G4double x = 0;
199 G4double y = 0;
200 G4double x1= 0;
201 G4double x2= 0;
202 G4double xlog = 0.0;
203
204 G4double costm = std::max(cosTMax,cosTetMaxElec);
205 G4double fb = screenZ*factB;
206
207 // scattering off electrons
208 if(costm < 1.0) {
209 x = (1.0 - costm)/screenZ;
210 if(x < numlimit) {
211 x2 = 0.5*x*x;
212 y = x2*(1.0 - 1.3333333*x + 3*x2);
213 if(0.0 < factB) { y -= fb*x2*x*(0.6666667 - x); }
214 } else {
215 x1= x/(1 + x);
216 xlog = log(1.0 + x);
217 y = xlog - x1;
218 if(0.0 < factB) { y -= fb*(x + x1 - 2*xlog); }
219 }
220
221 if(y < 0.0) {
222 ++nwarnings;
223 if(nwarnings < nwarnlimit) {
224 G4cout << "G4WentzelOKandVIxSection::ComputeTransportCrossSectionPerAtom scattering on e- <0"
225 << G4endl;
226 G4cout << "y= " << y
227 << " e(MeV)= " << tkin << " p(MeV/c)= " << sqrt(mom2)
228 << " Z= " << targetZ << " "
229 << particle->GetParticleName() << G4endl;
230 G4cout << " 1-costm= " << 1.0-costm << " screenZ= " << screenZ
231 << " x= " << x << G4endl;
232 }
233 y = 0.0;
234 }
235 xSection = y;
236 }
237 /*
238 G4cout << "G4WentzelVI:XS per A " << " Z= " << targetZ
239 << " e(MeV)= " << tkin/MeV << " XSel= " << xSection
240 << " cut(MeV)= " << ecut/MeV
241 << " zmaxE= " << (1.0 - cosTetMaxElec)/screenZ
242 << " zmaxN= " << (1.0 - cosThetaMax)/screenZ
243 << " 1-costm= " << 1.0 - cosThetaMax << G4endl;
244 */
245 // scattering off nucleus
246 if(cosTMax < 1.0) {
247 x = (1.0 - cosTMax)/screenZ;
248 if(x < numlimit) {
249 x2 = 0.5*x*x;
250 y = x2*(1.0 - 1.3333333*x + 3*x2);
251 if(0.0 < factB) { y -= fb*x2*x*(0.6666667 - x); }
252 } else {
253 x1= x/(1 + x);
254 xlog = log(1.0 + x);
255 y = xlog - x1;
256 if(0.0 < factB) { y -= fb*(x + x1 - 2*xlog); }
257 }
258
259 if(y < 0.0) {
260 ++nwarnings;
261 if(nwarnings < nwarnlimit) {
262 G4cout << "G4WentzelOKandVIxSection::ComputeTransportCrossSectionPerAtom scattering on e- <0"
263 << G4endl;
264 G4cout << "y= " << y
265 << " e(MeV)= " << tkin << " Z= " << targetZ << " "
266 << particle->GetParticleName() << G4endl;
267 G4cout << " formfactA= " << formfactA << " screenZ= " << screenZ
268 << " x= " << " x1= " << x1 <<G4endl;
269 }
270 y = 0.0;
271 }
272 xSection += y*targetZ;
273 }
274 xSection *= kinFactor;
275 /*
276 G4cout << "Z= " << targetZ << " XStot= " << xSection/barn
277 << " screenZ= " << screenZ << " formF= " << formfactA
278 << " for " << particle->GetParticleName()
279 << " m= " << mass << " 1/v= " << sqrt(invbeta2) << " p= " << sqrt(mom2)
280 << " x= " << x
281 << G4endl;
282 */
283 return xSection;
284}
285
286//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
287
290 G4double cosTMax,
291 G4double elecRatio)
292{
293 G4ThreeVector v(0.0,0.0,1.0);
294
295 G4double formf = formfactA;
296 G4double cost1 = cosTMin;
297 G4double cost2 = cosTMax;
298 if(elecRatio > 0.0) {
299 if(G4UniformRand() <= elecRatio) {
300 formf = 0.0;
301 cost1 = std::max(cost1,cosTetMaxElec);
302 cost2 = std::max(cost2,cosTetMaxElec);
303 }
304 }
305 if(cost1 < cost2) { return v; }
306
307 G4double w1 = 1. - cost1 + screenZ;
308 G4double w2 = 1. - cost2 + screenZ;
309 G4double z1 = w1*w2/(w1 + G4UniformRand()*(w2 - w1)) - screenZ;
310
311 //G4double fm = 1.0 + formf*z1/(1.0 + (mass + tkin)*z1/targetMass);
312 G4double fm = 1.0 + formf*z1;
313 //G4double grej = (1. - z1*factB)/( (1.0 + z1*factD)*fm*fm );
314 G4double grej = (1. - z1*factB + factB1*targetZ*sqrt(z1*factB)*(2 - z1))/( (1.0 + z1*factD)*fm*fm );
315 // "false" scattering
316 if( G4UniformRand() > grej ) { return v; }
317 // }
318 G4double cost = 1.0 - z1;
319
320 if(cost > 1.0) { cost = 1.0; }
321 else if(cost < -1.0) { cost =-1.0; }
322 G4double sint = sqrt((1.0 - cost)*(1.0 + cost));
323 //G4cout << "sint= " << sint << G4endl;
324 G4double phi = twopi*G4UniformRand();
325 G4double vx1 = sint*cos(phi);
326 G4double vy1 = sint*sin(phi);
327
328 // only direction is changed
329 v.set(vx1,vy1,cost);
330 return v;
331}
332
333//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
334
335void
336G4WentzelOKandVIxSection::ComputeMaxElectronScattering(G4double cutEnergy)
337{
338 if(mass > MeV) {
339 G4double ratio = electron_mass_c2/mass;
340 G4double tau = tkin/mass;
341 G4double tmax = 2.0*electron_mass_c2*tau*(tau + 2.)/
342 (1.0 + 2.0*ratio*(tau + 1.0) + ratio*ratio);
343 //tmax = std::min(tmax, targetZ*targetZ*10*eV);
344 cosTetMaxElec = 1.0 - std::min(cutEnergy, tmax)*electron_mass_c2/mom2;
345 } else {
346
347 G4double tmax = tkin;
348 if(particle == theElectron) { tmax *= 0.5; }
349 //tmax = std::min(tmax, targetZ*targetZ*10*eV);
350 G4double t = std::min(cutEnergy, tmax);
351 G4double mom21 = t*(t + 2.0*electron_mass_c2);
352 G4double t1 = tkin - t;
353 //G4cout <<"tkin=" <<tkin<<" tmax= "<<tmax<<" t= "
354 //<<t<< " t1= "<<t1<<" cut= "<<ecut<<G4endl;
355 if(t1 > 0.0) {
356 G4double mom22 = t1*(t1 + 2.0*mass);
357 G4double ctm = (mom2 + mom22 - mom21)*0.5/sqrt(mom2*mom22);
358 if(ctm < 1.0) { cosTetMaxElec = ctm; }
359 if(particle == theElectron && cosTetMaxElec < 0.0) { cosTetMaxElec = 0.0; }
360 }
361 }
362}
363
364//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:53
void set(double x, double y, double z)
static G4Electron * Electron()
Definition: G4Electron.cc:94
static G4LossTableManager * Instance()
G4double FactorForAngleLimit() const
G4double GetA27(G4int Z)
static G4NistManager * Instance()
G4double GetAtomicMassAmu(const G4String &symb) const
G4double GetPDGCharge() const
const G4String & GetParticleName() const
static G4Positron * Positron()
Definition: G4Positron.cc:94
static G4Pow * GetInstance()
Definition: G4Pow.cc:50
G4double Z23(G4int Z)
Definition: G4Pow.hh:134
G4double Z13(G4int Z)
Definition: G4Pow.hh:110
static G4Proton * Proton()
Definition: G4Proton.cc:93
G4ThreeVector SampleSingleScattering(G4double CosThetaMin, G4double CosThetaMax, G4double elecRatio=0.0)
void Initialise(const G4ParticleDefinition *, G4double CosThetaLim)
G4double ComputeTransportCrossSectionPerAtom(G4double CosThetaMax)
void SetupParticle(const G4ParticleDefinition *)
G4double SetupTarget(G4int Z, G4double cut=DBL_MAX)
#define DBL_MAX
Definition: templates.hh:83