Geant4  10.01.p01
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 85306 2014-10-27 14:17:47Z gcosmo $
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 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
45 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
46 
48 #include "G4PhysicalConstants.hh"
49 #include "G4SystemOfUnits.hh"
50 #include "Randomize.hh"
51 #include "G4Electron.hh"
52 #include "G4Positron.hh"
53 #include "G4Proton.hh"
54 #include "G4EmParameters.hh"
55 #include "G4Log.hh"
56 #include "G4Exp.hh"
57 
58 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
59 
63 
64 using namespace std;
65 
67  numlimit(0.1),
68  nwarnings(0),
69  nwarnlimit(50),
70  isCombined(combined),
71  cosThetaMax(-1.0),
72  alpha2(fine_structure_const*fine_structure_const)
73 {
79  lowEnergyLimit = 1.0*eV;
80  G4double p0 = electron_mass_c2*classic_electr_radius;
81  coeff = twopi*p0*p0;
82  particle = 0;
83 
84  // Thomas-Fermi screening radii
85  // Formfactors from A.V. Butkevich et al., NIM A 488 (2002) 282
86 
87  if(0.0 == ScreenRSquare[0]) {
88  G4double a0 = electron_mass_c2/0.88534;
89  G4double constn = 6.937e-6/(MeV*MeV);
90 
91  ScreenRSquare[0] = alpha2*a0*a0;
93  for(G4int j=1; j<100; ++j) {
94  G4double x = a0*fG4pow->Z13(j);
95  if(1 == j) { ScreenRSquare[j] = 0.5*alpha2*a0*a0; }
96  else {
97  ScreenRSquare[j] = 0.5*(1 + G4Exp(-j*j*0.001))*alpha2*x*x;
98  ScreenRSquareElec[j] = 0.5*alpha2*x*x;
99  }
100  x = fNistManager->GetA27(j);
101  FormFactor[j] = constn*x*x;
102  }
103  }
104  currentMaterial = 0;
105  factB = factD = formfactA = screenZ = 0.0;
107 
108  factB1= 0.5*CLHEP::pi*fine_structure_const;
109 
110  tkin = mom2 = momCM2 = factorA2 = mass = spin = chargeSquare = charge3 = 0.0;
111  ecut = etag = DBL_MAX;
112  targetZ = 0;
113  targetMass = proton_mass_c2;
114  particle = 0;
115 }
116 
117 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
118 
120 {}
121 
122 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
123 
125  G4double cosThetaLim)
126 {
127  SetupParticle(p);
128  tkin = mom2 = momCM2 = 0.0;
129  ecut = etag = DBL_MAX;
130  targetZ = 0;
131 
132  // cosThetaMax is below 1.0 only when MSC is combined with SS
133  if(isCombined) { cosThetaMax = cosThetaLim; }
134 
136  *CLHEP::hbarc/CLHEP::fermi;
137  factorA2 = 0.5*a*a;
138  currentMaterial = 0;
139 
140  //G4cout << "G4WentzelOKandVIxSection::Initialise mass= " << mass
141  // << " " << p->GetParticleName()
142  // << " cosThetaMax= " << cosThetaMax << G4endl;
143 
144 }
145 
146 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
147 
149 {
150  particle = p;
151  mass = particle->GetPDGMass();
152  spin = particle->GetPDGSpin();
153  if(0.0 != spin) { spin = 0.5; }
154  G4double q = std::fabs(particle->GetPDGCharge()/eplus);
155  chargeSquare = q*q;
156  charge3 = chargeSquare*q;
157  tkin = 0.0;
158  currentMaterial = 0;
159  targetZ = 0;
160 }
161 
162 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
163 
164 G4double
166 {
167  G4double cosTetMaxNuc2 = cosTetMaxNuc;
168  if(Z != targetZ || tkin != etag) {
169  etag = tkin;
170  targetZ = Z;
171  if(targetZ > 99) { targetZ = 99; }
172  G4double massT = proton_mass_c2;
173  if(targetZ > 1) {
174  massT = fNistManager->GetAtomicMassAmu(targetZ)*CLHEP::amu_c2;
175  }
176  SetTargetMass(massT);
177 
179 
180  if(1 == Z) {
182  } else if(mass > MeV) {
183  screenZ = std::min(Z*1.13,1.13 +3.76*Z*Z*invbeta2*alpha2*chargeSquare)*
185  } else {
186  G4double tau = tkin/mass;
187  screenZ = std::min(Z*1.13,(1.13 +3.76*Z*Z
188  *invbeta2*alpha2*std::sqrt(tau/(tau + fG4pow->Z23(targetZ)))))*
190  }
191  if(targetZ == 1 && cosTetMaxNuc2 < 0.0 && particle == theProton) {
192  cosTetMaxNuc2 = 0.0;
193  }
195 
196  cosTetMaxElec = 1.0;
198  }
199  return cosTetMaxNuc2;
200 }
201 
202 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
203 
204 G4double
206 {
207  G4double xSection = 0.0;
208  if(cosTMax >= 1.0) { return xSection; }
209 
210  G4double x = 0;
211  G4double y = 0;
212  G4double x1= 0;
213  G4double x2= 0;
214  G4double xlog = 0.0;
215 
216  G4double costm = std::max(cosTMax,cosTetMaxElec);
217  G4double fb = screenZ*factB;
218 
219  // scattering off electrons
220  if(costm < 1.0) {
221  x = (1.0 - costm)/screenZ;
222  if(x < numlimit) {
223  x2 = 0.5*x*x;
224  y = x2*(1.0 - 1.3333333*x + 3*x2);
225  if(0.0 < factB) { y -= fb*x2*x*(0.6666667 - x); }
226  } else {
227  x1= x/(1 + x);
228  xlog = G4Log(1.0 + x);
229  y = xlog - x1;
230  if(0.0 < factB) { y -= fb*(x + x1 - 2*xlog); }
231  }
232 
233  if(y < 0.0) {
234  ++nwarnings;
235  if(nwarnings < nwarnlimit) {
236  G4cout << "G4WentzelOKandVIxSection::ComputeTransportCrossSectionPerAtom"
237  << " scattering on e- <0"
238  << G4endl;
239  G4cout << "y= " << y
240  << " e(MeV)= " << tkin << " p(MeV/c)= " << sqrt(mom2)
241  << " Z= " << targetZ << " "
243  G4cout << " 1-costm= " << 1.0-costm << " screenZ= " << screenZ
244  << " x= " << x << G4endl;
245  }
246  y = 0.0;
247  }
248  xSection = y;
249  }
250  /*
251  G4cout << "G4WentzelVI:XS per A " << " Z= " << targetZ
252  << " e(MeV)= " << tkin/MeV << " XSel= " << xSection
253  << " cut(MeV)= " << ecut/MeV
254  << " zmaxE= " << (1.0 - cosTetMaxElec)/screenZ
255  << " zmaxN= " << (1.0 - cosThetaMax)/screenZ
256  << " 1-costm= " << 1.0 - cosThetaMax << G4endl;
257  */
258  // scattering off nucleus
259  if(cosTMax < 1.0) {
260  x = (1.0 - cosTMax)/screenZ;
261  if(x < numlimit) {
262  x2 = 0.5*x*x;
263  y = x2*(1.0 - 1.3333333*x + 3*x2);
264  if(0.0 < factB) { y -= fb*x2*x*(0.6666667 - x); }
265  } else {
266  x1= x/(1 + x);
267  xlog = G4Log(1.0 + x);
268  y = xlog - x1;
269  if(0.0 < factB) { y -= fb*(x + x1 - 2*xlog); }
270  }
271 
272  if(y < 0.0) {
273  ++nwarnings;
274  if(nwarnings < nwarnlimit) {
275  G4cout << "G4WentzelOKandVIxSection::ComputeTransportCrossSectionPerAtom"
276  << " scattering on nucleus <0"
277  << G4endl;
278  G4cout << "y= " << y
279  << " e(MeV)= " << tkin << " Z= " << targetZ << " "
281  G4cout << " formfactA= " << formfactA << " screenZ= " << screenZ
282  << " x= " << " x1= " << x1 <<G4endl;
283  }
284  y = 0.0;
285  }
286  xSection += y*targetZ;
287  }
288  xSection *= kinFactor;
289 
290  /*
291  G4cout << "Z= " << targetZ << " XStot= " << xSection/barn
292  << " screenZ= " << screenZ << " formF= " << formfactA
293  << " for " << particle->GetParticleName()
294  << " m= " << mass << " 1/v= " << sqrt(invbeta2) << " p= " << sqrt(mom2)
295  << " x= " << x
296  << G4endl;
297  */
298  return xSection;
299 }
300 
301 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
302 
305  G4double cosTMax,
306  G4double elecRatio)
307 {
308  temp.set(0.0,0.0,1.0);
309 
310  G4double formf = formfactA;
311  G4double cost1 = cosTMin;
312  G4double cost2 = cosTMax;
313  if(elecRatio > 0.0) {
314  if(G4UniformRand() <= elecRatio) {
315  formf = 0.0;
316  cost1 = std::max(cost1,cosTetMaxElec);
317  cost2 = std::max(cost2,cosTetMaxElec);
318  }
319  }
320  if(cost1 > cost2) {
321  G4double z1, z2;
322  G4double w1 = 1. - cost1 + screenZ;
323  G4double w2 = 1. - cost2 + screenZ;
324  do {
325  z1 = w1*w2/(w1 + G4UniformRand()*(w2 - w1)) - screenZ;
326  z2 = 1. - z1*factB + factB1*targetZ*sqrt(z1*factB)*(2 - z1);
327  } while(G4UniformRand() > z2);
328 
329  G4double fm = 1.0 + formf*z1;
330  G4double grej = (1.0 + z1*factD)*fm*fm;
331 
332  if(G4UniformRand()*grej < 1.0) {
333  // exclude "false" scattering due to formfactor
334  G4double cost = 1.0 - z1;
335  if(cost > 1.0) { cost = 1.0; }
336  else if(cost < -1.0) { cost =-1.0; }
337  G4double sint = sqrt((1.0 - cost)*(1.0 + cost));
338  //G4cout << "sint= " << sint << G4endl;
339  G4double phi = twopi*G4UniformRand();
340  temp.set(sint*cos(phi),sint*sin(phi),cost);
341  }
342  }
343  return temp;
344 }
345 
346 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
347 
348 void
350 {
351  if(mass > MeV) {
352  G4double ratio = electron_mass_c2/mass;
353  G4double tau = tkin/mass;
354  G4double tmax = 2.0*electron_mass_c2*tau*(tau + 2.)/
355  (1.0 + 2.0*ratio*(tau + 1.0) + ratio*ratio);
356  //tmax = std::min(tmax, targetZ*targetZ*10*eV);
357  cosTetMaxElec = 1.0 - std::min(cutEnergy, tmax)*electron_mass_c2/mom2;
358  } else {
359 
360  G4double tmax = tkin;
361  if(particle == theElectron) { tmax *= 0.5; }
362  //tmax = std::min(tmax, targetZ*targetZ*10*eV);
363  G4double t = std::min(cutEnergy, tmax);
364  G4double mom21 = t*(t + 2.0*electron_mass_c2);
365  G4double t1 = tkin - t;
366  //G4cout <<"tkin=" <<tkin<<" tmax= "<<tmax<<" t= "
367  //<<t<< " t1= "<<t1<<" cut= "<<ecut<<G4endl;
368  if(t1 > 0.0) {
369  G4double mom22 = t1*(t1 + 2.0*mass);
370  G4double ctm = (mom2 + mom22 - mom21)*0.5/sqrt(mom2*mom22);
371  if(ctm < 1.0) { cosTetMaxElec = ctm; }
372  if(particle == theElectron && cosTetMaxElec < 0.0) {
373  cosTetMaxElec = 0.0;
374  }
375  }
376  }
377 }
378 
379 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
380 
381 G4double
383 {
384  return 0.0;
385 }
386 
387 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
const G4ParticleDefinition * theElectron
const G4ParticleDefinition * theProton
const G4double a0
static G4Pow * GetInstance()
Definition: G4Pow.cc:55
void SetupParticle(const G4ParticleDefinition *)
static const double MeV
Definition: G4SIunits.hh:193
static G4double ScreenRSquareElec[100]
G4WentzelOKandVIxSection(G4bool combined=true)
CLHEP::Hep3Vector G4ThreeVector
G4double ComputeSecondTransportMoment(G4double CosThetaMax)
const G4double pi
const G4ParticleDefinition * particle
void ComputeMaxElectronScattering(G4double cut)
G4double GetA27(G4int Z)
G4double a
Definition: TRTMaterials.hh:39
int G4int
Definition: G4Types.hh:78
static G4NistManager * Instance()
G4double SetupTarget(G4int Z, G4double cut=DBL_MAX)
const G4String & GetParticleName() const
G4double ComputeTransportCrossSectionPerAtom(G4double CosThetaMax)
#define G4UniformRand()
Definition: Randomize.hh:95
G4GLOB_DLL std::ostream G4cout
G4double Z13(G4int Z) const
Definition: G4Pow.hh:127
bool G4bool
Definition: G4Types.hh:79
static G4Proton * Proton()
Definition: G4Proton.cc:93
G4double G4Log(G4double x)
Definition: G4Log.hh:230
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
static const double eV
Definition: G4SIunits.hh:194
static G4Positron * Positron()
Definition: G4Positron.cc:94
G4double GetPDGMass() const
const G4ParticleDefinition * thePositron
T max(const T t1, const T t2)
brief Return the largest of the two arguments
static G4EmParameters * Instance()
G4double GetAtomicMassAmu(const G4String &symb) const
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
G4double Z23(G4int Z) const
Definition: G4Pow.hh:151
G4double GetPDGSpin() const
static G4Electron * Electron()
Definition: G4Electron.cc:94
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
static const double eplus
Definition: G4SIunits.hh:178
G4double GetPDGCharge() const
#define DBL_MAX
Definition: templates.hh:83
static const double fermi
Definition: G4SIunits.hh:93
void Initialise(const G4ParticleDefinition *, G4double CosThetaLim)
G4double FactorForAngleLimit() const
G4ThreeVector & SampleSingleScattering(G4double CosThetaMin, G4double CosThetaMax, G4double elecRatio=0.0)