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