Geant4-11
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//
27// -------------------------------------------------------------------
28//
29// GEANT4 Class file
30//
31//
32// File name: G4WentzelOKandVIxSection
33//
34// Author: V.Ivanchenko
35//
36// Creation date: 09.04.2008 from G4MuMscModel
37//
38// Modifications:
39//
40// -------------------------------------------------------------------
41//
42
43//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
44//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
45
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#ifdef G4MULTITHREADED
65G4Mutex G4WentzelOKandVIxSection::WentzelOKandVIxSectionMutex = G4MUTEX_INITIALIZER;
66#endif
67
70const G4double numlimit = 0.1;
71const G4int nwarnlimit = 50;
72
73using namespace std;
74
76 temp(0.,0.,0.),
77 isCombined(comb)
78{
81
85
87 coeff = CLHEP::twopi*p0*p0;
89}
90
91//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
92
94{
95 delete fMottXSection;
96}
97
98//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
99
101 G4double cosThetaLim)
102{
103 SetupParticle(p);
104 tkin = mom2 = momCM2 = 0.0;
105 ecut = etag = DBL_MAX;
106 targetZ = 0;
107
108 // cosThetaMax is below 1.0 only when MSC is combined with SS
109 if(isCombined) { cosThetaMax = cosThetaLim; }
112 factorA2 = 0.5*a*a;
113 currentMaterial = nullptr;
114
116 if(0.0 == ScreenRSquare[0]) { InitialiseA(); }
117
118 // Mott corrections always added
119 if((p == theElectron || p == thePositron) && !fMottXSection) {
121 fMottXSection->Initialise(p, 1.0);
122 }
123 /*
124 G4cout << "G4WentzelOKandVIxSection::Initialise for "
125 << p->GetParticleName() << " cosThetaMax= " << cosThetaMax
126 << " " << ScreenRSquare[0] << " coeff= " << coeff << G4endl;
127 */
128}
129
130//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
131
133{
134 // Thomas-Fermi screening radii
135 // Formfactors from A.V. Butkevich et al., NIM A 488 (2002) 282
136#ifdef G4MULTITHREADED
137 G4MUTEXLOCK(&G4WentzelOKandVIxSection::WentzelOKandVIxSectionMutex);
138 if(0.0 == ScreenRSquare[0]) {
139#endif
140 const G4double invmev2 = 1./(CLHEP::MeV*CLHEP::MeV);
142 G4double constn = 6.937e-6*invmev2;
144
145 G4double afact = 0.5*fct*alpha2*a0*a0;
146 ScreenRSquare[0] = afact;
147 ScreenRSquare[1] = afact;
148 ScreenRSquareElec[1] = afact;
149 FormFactor[1] = 3.097e-6*invmev2;
150
151 for(G4int j=2; j<100; ++j) {
152 G4double x = fG4pow->Z13(j);
153 ScreenRSquare[j] = afact*(1 + G4Exp(-j*j*0.001))*x*x;
154 ScreenRSquareElec[j] = afact*x*x;
155 x = fNistManager->GetA27(j);
156 FormFactor[j] = constn*x*x;
157 }
158#ifdef G4MULTITHREADED
159 }
160 G4MUTEXUNLOCK(&G4WentzelOKandVIxSection::WentzelOKandVIxSectionMutex);
161#endif
162
163 //G4cout << "G4WentzelOKandVIxSection::Initialise mass= " << mass
164 // << " " << p->GetParticleName()
165 // << " cosThetaMax= " << cosThetaMax << G4endl;
166
167}
168
169//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
170
172{
173 particle = p;
176 if(0.0 != spin) { spin = 0.5; }
177 G4double q = std::abs(particle->GetPDGCharge()/eplus);
178 chargeSquare = q*q;
180 tkin = 0.0;
181 currentMaterial = nullptr;
182 targetZ = 0;
183}
184
185//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
186
189{
190 if(ekin != tkin || mat != currentMaterial) {
191 currentMaterial = mat;
192 tkin = ekin;
193 mom2 = tkin*(tkin + 2.0*mass);
194 invbeta2 = 1.0 + mass*mass/mom2;
198 : cosThetaMax;
199 }
200 return cosTetMaxNuc;
201}
202
203//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
204
207{
208 G4double cosTetMaxNuc2 = cosTetMaxNuc;
209 if(Z != targetZ || tkin != etag) {
210 etag = tkin;
211 targetZ = std::min(Z, 99);
212 G4double massT = (1 == Z) ? CLHEP::proton_mass_c2 :
214 SetTargetMass(massT);
215
218 fMottFactor = (1.0 + 2.0e-4*Z*Z);
219 }
220
221 if(1 == Z) {
223 } else if(mass > MeV) {
224 screenZ = std::min(Z*1.13,1.13 +3.76*Z*Z*invbeta2*alpha2*chargeSquare)*
226 } else {
227 G4double tau = tkin/mass;
228 screenZ = std::min(Z*1.13,(1.13 +3.76*Z*Z
229 *invbeta2*alpha2*std::sqrt(tau/(tau + fG4pow->Z23(targetZ)))))*
231 }
232 if(targetZ == 1 && particle == theProton && cosTetMaxNuc2 < 0.0) {
233 cosTetMaxNuc2 = 0.0;
234 }
236
237 cosTetMaxElec = 1.0;
239 }
240 //G4cout << "SetupTarget: Z= " << targetZ << " kinFactor= " << kinFactor
241 // << " fMottFactor= " << fMottFactor << " screenZ= " << screenZ <<G4endl;
242 return cosTetMaxNuc2;
243}
244
245//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
246
249{
250 G4double xSection = 0.0;
251 if(cosTMax >= 1.0) { return xSection; }
252
253 G4double costm = std::max(cosTMax,cosTetMaxElec);
255
256 // scattering off electrons
257 if(costm < 1.0) {
258 G4double x = (1.0 - costm)/screenZ;
259 if(x < numlimit) {
260 G4double x2 = 0.5*x*x;
261 xSection = x2*((1.0 - 1.3333333*x + 3*x2) - fb*x*(0.6666667 - x));
262 } else {
263 G4double x1= x/(1 + x);
264 G4double xlog = G4Log(1.0 + x);
265 xSection = xlog - x1 - fb*(x + x1 - 2*xlog);
266 }
267
268 if(xSection < 0.0) {
269 ++nwarnings;
270 if(nwarnings < nwarnlimit) {
271 G4cout << "G4WentzelOKandVIxSection::ComputeTransportCrossSectionPerAtom"
272 << " scattering on e- <0"
273 << G4endl;
274 G4cout << "cross= " << xSection
275 << " e(MeV)= " << tkin << " p(MeV/c)= " << sqrt(mom2)
276 << " Z= " << targetZ << " "
278 G4cout << " 1-costm= " << 1.0-costm << " screenZ= " << screenZ
279 << " x= " << x << G4endl;
280 }
281 xSection = 0.0;
282 }
283 }
284 /*
285 G4cout << "G4WentzelOKandVIxSection::ComputeTransportCrossSectionPerAtom: \n"
286 << " Z= " << targetZ
287 << " e(MeV)= " << tkin/MeV << " XSel= " << xSection
288 << " zmaxE= " << (1.0 - cosTetMaxElec)/screenZ
289 << " zmaxN= " << (1.0 - cosThetaMax)/screenZ
290 << " 1-costm= " << 1.0 - cosThetaMax << G4endl;
291 */
292 // scattering off nucleus
293 if(cosTMax < 1.0) {
294 G4double x = (1.0 - cosTMax)/screenZ;
295 G4double y;
296 if(x < numlimit) {
297 G4double x2 = 0.5*x*x;
298 y = x2*((1.0 - 1.3333333*x + 3*x2) - fb*x*(0.6666667 - x));
299 } else {
300 G4double x1= x/(1 + x);
301 G4double xlog = G4Log(1.0 + x);
302 y = xlog - x1 - fb*(x + x1 - 2*xlog);
303 }
304
305 if(y < 0.0) {
306 ++nwarnings;
307 if(nwarnings < nwarnlimit) {
308 G4cout << "G4WentzelOKandVIxSection::ComputeTransportCrossSectionPerAtom"
309 << " scattering on nucleus <0"
310 << G4endl;
311 G4cout << "y= " << y
312 << " e(MeV)= " << tkin << " Z= " << targetZ << " "
314 G4cout << " formfactA= " << formfactA << " screenZ= " << screenZ
315 << " x= " << x <<G4endl;
316 }
317 y = 0.0;
318 }
319 xSection += y*targetZ;
320 }
321 xSection *= kinFactor;
322
323 /*
324 G4cout << "Z= " << targetZ << " XStot= " << xSection/barn
325 << " screenZ= " << screenZ << " formF= " << formfactA
326 << " for " << particle->GetParticleName()
327 << " m= " << mass << " 1/v= " << sqrt(invbeta2)
328 << " p= " << sqrt(mom2)
329 << " x= " << x << G4endl;
330 */
331 return xSection;
332}
333
334//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
335
338 G4double cosTMax,
339 G4double elecRatio)
340{
341 temp.set(0.0,0.0,1.0);
342 CLHEP::HepRandomEngine* rndmEngineMod = G4Random::getTheEngine();
343
344 G4double formf = formfactA;
345 G4double cost1 = cosTMin;
346 G4double cost2 = cosTMax;
347 if(elecRatio > 0.0) {
348 if(rndmEngineMod->flat() <= elecRatio) {
349 formf = 0.0;
350 cost1 = std::max(cost1,cosTetMaxElec);
351 cost2 = std::max(cost2,cosTetMaxElec);
352 }
353 }
354 if(cost1 > cost2) {
355
356 G4double w1 = 1. - cost1 + screenZ;
357 G4double w2 = 1. - cost2 + screenZ;
358 G4double z1 = w1*w2/(w1 + rndmEngineMod->flat()*(w2 - w1)) - screenZ;
359
360 G4double fm = 1.0;
362 fm += formf*z1;
363 fm = 1.0/(fm*fm);
364 } else if(fNucFormfactor == fGaussianNF) {
365 fm = G4Exp(-2*formf*z1);
366 } else if(fNucFormfactor == fFlatNF) {
367 static const G4double ccoef = 0.00508/MeV;
368 G4double x = std::sqrt(2.*mom2*z1)*ccoef*2.;
369 fm = FlatFormfactor(x);
370 fm *= FlatFormfactor(x*0.6
372 }
373 G4double grej;
374 if(fMottXSection) {
376 grej = fMottXSection->RatioMottRutherfordCosT(std::sqrt(z1))*fm*fm;
377 } else {
378 grej = (1. - z1*factB + factB1*targetZ*sqrt(z1*factB)*(2. - z1))
379 *fm*fm/(1.0 + z1*factD);
380 }
381 // G4cout << "SampleSingleScattering: E= " << tkin << " z1= "
382 // << z1 << " grej= "<< grej << " mottFact= "<< fMottFactor<< G4endl;
383 if(fMottFactor*rndmEngineMod->flat() <= grej ) {
384 // exclude "false" scattering due to formfactor and spin effect
385 G4double cost = 1.0 - z1;
386 if(cost > 1.0) { cost = 1.0; }
387 else if(cost < -1.0) { cost =-1.0; }
388 G4double sint = sqrt((1.0 - cost)*(1.0 + cost));
389 //G4cout << "sint= " << sint << G4endl;
390 G4double phi = twopi*rndmEngineMod->flat();
391 temp.set(sint*cos(phi),sint*sin(phi),cost);
392 }
393 }
394 return temp;
395}
396
397//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
398
399void
401{
402 if(mass > MeV) {
404 G4double tau = tkin/mass;
405 G4double tmax = 2.0*electron_mass_c2*tau*(tau + 2.)/
406 (1.0 + 2.0*ratio*(tau + 1.0) + ratio*ratio);
407 cosTetMaxElec = 1.0 - std::min(cutEnergy, tmax)*electron_mass_c2/mom2;
408 } else {
409
410 G4double tmax = (particle == theElectron) ? 0.5*tkin : tkin;
411 G4double t = std::min(cutEnergy, tmax);
412 G4double mom21 = t*(t + 2.0*electron_mass_c2);
413 G4double t1 = tkin - t;
414 //G4cout <<"tkin=" <<tkin<<" tmax= "<<tmax<<" t= "
415 //<<t<< " t1= "<<t1<<" cut= "<<ecut<<G4endl;
416 if(t1 > 0.0) {
417 G4double mom22 = t1*(t1 + 2.0*mass);
418 G4double ctm = (mom2 + mom22 - mom21)*0.5/sqrt(mom2*mom22);
419 if(ctm < 1.0) { cosTetMaxElec = ctm; }
420 if(particle == theElectron && cosTetMaxElec < 0.0) {
421 cosTetMaxElec = 0.0;
422 }
423 }
424 }
425}
426
427//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
428
431{
432 return 0.0;
433}
434
435//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:179
const G4double a0
G4double G4Log(G4double x)
Definition: G4Log.hh:226
static constexpr double twopi
Definition: G4SIunits.hh:56
static constexpr double eplus
Definition: G4SIunits.hh:184
static constexpr double MeV
Definition: G4SIunits.hh:200
#define G4MUTEX_INITIALIZER
Definition: G4Threading.hh:85
#define G4MUTEXLOCK(mutex)
Definition: G4Threading.hh:251
#define G4MUTEXUNLOCK(mutex)
Definition: G4Threading.hh:254
std::mutex G4Mutex
Definition: G4Threading.hh:81
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
const G4double alpha2
const G4int nwarnlimit
const G4double factB1
const G4double numlimit
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
void set(double x, double y, double z)
virtual double flat()=0
static G4Electron * Electron()
Definition: G4Electron.cc:93
static G4EmParameters * Instance()
G4double ScreeningFactor() const
G4NuclearFormfactorType NuclearFormfactorType() const
G4double FactorForAngleLimit() const
G4double GetInvA23() const
G4IonisParamMat * GetIonisation() const
Definition: G4Material.hh:222
G4double GetA27(G4int Z) const
static G4NistManager * Instance()
G4double GetAtomicMassAmu(const G4String &symb) const
G4double GetPDGCharge() const
const G4String & GetParticleName() const
static G4Positron * Positron()
Definition: G4Positron.cc:93
static G4Pow * GetInstance()
Definition: G4Pow.cc:41
G4double A13(G4double A) const
Definition: G4Pow.cc:120
G4double Z13(G4int Z) const
Definition: G4Pow.hh:123
G4double Z23(G4int Z) const
Definition: G4Pow.hh:125
static G4Proton * Proton()
Definition: G4Proton.cc:92
void SetupKinematic(G4double kinEnergy, G4int Z)
void Initialise(const G4ParticleDefinition *, G4double cosThetaLim)
G4double RatioMottRutherfordCosT(G4double sin2t2)
G4double SetupTarget(G4int Z, G4double cut)
static G4double ScreenRSquareElec[100]
G4double ComputeSecondTransportMoment(G4double CosThetaMax)
void ComputeMaxElectronScattering(G4double cut)
const G4ParticleDefinition * theProton
void Initialise(const G4ParticleDefinition *, G4double CosThetaLim)
const G4ParticleDefinition * thePositron
const G4ParticleDefinition * particle
G4WentzelOKandVIxSection(G4bool comb=true)
G4double ComputeTransportCrossSectionPerAtom(G4double CosThetaMax)
void SetupParticle(const G4ParticleDefinition *)
virtual G4double SetupKinematic(G4double kinEnergy, const G4Material *mat)
G4ScreeningMottCrossSection * fMottXSection
G4NuclearFormfactorType fNucFormfactor
G4ThreeVector & SampleSingleScattering(G4double CosThetaMin, G4double CosThetaMax, G4double elecRatio)
const G4ParticleDefinition * theElectron
static constexpr double electron_mass_c2
static constexpr double amu_c2
static constexpr double proton_mass_c2
static constexpr double fine_structure_const
static constexpr double twopi
Definition: SystemOfUnits.h:56
static constexpr double MeV
static constexpr double classic_electr_radius
static constexpr double hbarc
static constexpr double pi
Definition: SystemOfUnits.h:55
static constexpr double fermi
Definition: SystemOfUnits.h:84
T max(const T t1, const T t2)
brief Return the largest of the two arguments
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
float electron_mass_c2
Definition: hepunit.py:273
#define DBL_MAX
Definition: templates.hh:62