Geant4 11.1.1
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//
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#include "G4AutoLock.hh"
58
59//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
60
64
65namespace
66{
67 G4Mutex theWOKVIMutex = G4MUTEX_INITIALIZER;
68}
69
70const G4double alpha2 = CLHEP::fine_structure_const*CLHEP::fine_structure_const;
71const G4double factB1= 0.5*CLHEP::pi*CLHEP::fine_structure_const;
72const G4double numlimit = 0.1;
73const G4int nwarnlimit = 50;
74
75using namespace std;
76
78 temp(0.,0.,0.),
79 isCombined(comb)
80{
83
87
88 G4double p0 = CLHEP::electron_mass_c2*CLHEP::classic_electr_radius;
89 coeff = CLHEP::twopi*p0*p0;
90 targetMass = CLHEP::proton_mass_c2;
91}
92
93//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
94
96{
97 delete fMottXSection;
98}
99
100//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
101
103 G4double cosThetaLim)
104{
105 SetupParticle(p);
106 tkin = mom2 = momCM2 = 0.0;
107 ecut = etag = DBL_MAX;
108 targetZ = 0;
109
110 // cosThetaMax is below 1.0 only when MSC is combined with SS
111 if(isCombined) { cosThetaMax = cosThetaLim; }
113 G4double a = param->FactorForAngleLimit()*CLHEP::hbarc/CLHEP::fermi;
114 factorA2 = 0.5*a*a;
115 currentMaterial = nullptr;
116
118 if(0.0 == ScreenRSquare[0]) { InitialiseA(); }
119
120 // Mott corrections always added
121 if((p == theElectron || p == thePositron) && !fMottXSection) {
123 fMottXSection->Initialise(p, 1.0);
124 }
125 /*
126 G4cout << "G4WentzelOKandVIxSection::Initialise for "
127 << p->GetParticleName() << " cosThetaMax= " << cosThetaMax
128 << " " << ScreenRSquare[0] << " coeff= " << coeff << G4endl;
129 */
130}
131
132//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
133
135{
136 // Thomas-Fermi screening radii
137 // Formfactors from A.V. Butkevich et al., NIM A 488 (2002) 282
138 if(0.0 != ScreenRSquare[0]) { return; }
139 G4AutoLock l(&theWOKVIMutex);
140 if(0.0 == ScreenRSquare[0]) {
141 const G4double invmev2 = 1./(CLHEP::MeV*CLHEP::MeV);
142 G4double a0 = CLHEP::electron_mass_c2/0.88534;
143 G4double constn = 6.937e-6*invmev2;
145
146 G4double afact = 0.5*fct*alpha2*a0*a0;
147 ScreenRSquare[0] = afact;
148 ScreenRSquare[1] = afact;
149 ScreenRSquareElec[1] = afact;
150 FormFactor[1] = 3.097e-6*invmev2;
151
152 for(G4int j=2; j<100; ++j) {
153 G4double x = fG4pow->Z13(j);
154 ScreenRSquare[j] = afact*(1 + G4Exp(-j*j*0.001))*x*x;
155 ScreenRSquareElec[j] = afact*x*x;
156 x = fNistManager->GetA27(j);
157 FormFactor[j] = constn*x*x;
158 }
159 }
160 l.unlock();
161}
162
163//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
164
166{
167 particle = p;
170 if(0.0 != spin) { spin = 0.5; }
171 G4double q = std::abs(particle->GetPDGCharge()/eplus);
172 chargeSquare = q*q;
174 tkin = 0.0;
175 currentMaterial = nullptr;
176 targetZ = 0;
177}
178
179//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
180
183{
184 if(ekin != tkin || mat != currentMaterial) {
185 currentMaterial = mat;
186 tkin = ekin;
187 mom2 = tkin*(tkin + 2.0*mass);
188 invbeta2 = 1.0 + mass*mass/mom2;
191 std::max(cosThetaMax, 1.-factorA2*mat->GetIonisation()->GetInvA23()/mom2)
192 : cosThetaMax;
193 }
194 return cosTetMaxNuc;
195}
196
197//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
198
201{
202 G4double cosTetMaxNuc2 = cosTetMaxNuc;
203 if(Z != targetZ || tkin != etag) {
204 etag = tkin;
205 targetZ = std::min(Z, 99);
206 G4double massT = (1 == Z) ? CLHEP::proton_mass_c2 :
207 fNistManager->GetAtomicMassAmu(Z)*CLHEP::amu_c2;
208 SetTargetMass(massT);
209
212 fMottFactor = (1.0 + 2.0e-4*Z*Z);
213 }
214
215 if(1 == Z) {
217 } else if(mass > MeV) {
218 screenZ = std::min(Z*1.13,1.13 +3.76*Z*Z*invbeta2*alpha2*chargeSquare)*
220 } else {
221 G4double tau = tkin/mass;
222 screenZ = std::min(Z*1.13,(1.13 +3.76*Z*Z
223 *invbeta2*alpha2*std::sqrt(tau/(tau + fG4pow->Z23(targetZ)))))*
225 }
226 if(targetZ == 1 && particle == theProton && cosTetMaxNuc2 < 0.0) {
227 cosTetMaxNuc2 = 0.0;
228 }
230
231 cosTetMaxElec = 1.0;
233 }
234 //G4cout << "SetupTarget: Z= " << targetZ << " kinFactor= " << kinFactor
235 // << " fMottFactor= " << fMottFactor << " screenZ= " << screenZ <<G4endl;
236 return cosTetMaxNuc2;
237}
238
239//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
240
243{
244 G4double xSection = 0.0;
245 if(cosTMax >= 1.0) { return xSection; }
246
247 G4double costm = std::max(cosTMax,cosTetMaxElec);
249
250 // scattering off electrons
251 if(costm < 1.0) {
252 G4double x = (1.0 - costm)/screenZ;
253 if(x < numlimit) {
254 G4double x2 = 0.5*x*x;
255 xSection = x2*((1.0 - 1.3333333*x + 3*x2) - fb*x*(0.6666667 - x));
256 } else {
257 G4double x1= x/(1 + x);
258 G4double xlog = G4Log(1.0 + x);
259 xSection = xlog - x1 - fb*(x + x1 - 2*xlog);
260 }
261
262 if(xSection < 0.0) {
263 ++nwarnings;
264 if(nwarnings < nwarnlimit) {
265 G4cout << "G4WentzelOKandVIxSection::ComputeTransportCrossSectionPerAtom"
266 << " scattering on e- <0"
267 << G4endl;
268 G4cout << "cross= " << xSection
269 << " e(MeV)= " << tkin << " p(MeV/c)= " << sqrt(mom2)
270 << " Z= " << targetZ << " "
272 G4cout << " 1-costm= " << 1.0-costm << " screenZ= " << screenZ
273 << " x= " << x << G4endl;
274 }
275 xSection = 0.0;
276 }
277 }
278 /*
279 G4cout << "G4WentzelOKandVIxSection::ComputeTransportCrossSectionPerAtom: \n"
280 << " Z= " << targetZ
281 << " e(MeV)= " << tkin/MeV << " XSel= " << xSection
282 << " zmaxE= " << (1.0 - cosTetMaxElec)/screenZ
283 << " zmaxN= " << (1.0 - cosThetaMax)/screenZ
284 << " 1-costm= " << 1.0 - cosThetaMax << G4endl;
285 */
286 // scattering off nucleus
287 if(cosTMax < 1.0) {
288 G4double x = (1.0 - cosTMax)/screenZ;
289 G4double y;
290 if(x < numlimit) {
291 G4double x2 = 0.5*x*x;
292 y = x2*((1.0 - 1.3333333*x + 3*x2) - fb*x*(0.6666667 - x));
293 } else {
294 G4double x1= x/(1 + x);
295 G4double xlog = G4Log(1.0 + x);
296 y = xlog - x1 - fb*(x + x1 - 2*xlog);
297 }
298
299 if(y < 0.0) {
300 ++nwarnings;
301 if(nwarnings < nwarnlimit) {
302 G4cout << "G4WentzelOKandVIxSection::ComputeTransportCrossSectionPerAtom"
303 << " scattering on nucleus <0"
304 << G4endl;
305 G4cout << "y= " << y
306 << " e(MeV)= " << tkin << " Z= " << targetZ << " "
308 G4cout << " formfactA= " << formfactA << " screenZ= " << screenZ
309 << " x= " << x <<G4endl;
310 }
311 y = 0.0;
312 }
313 xSection += y*targetZ;
314 }
315 xSection *= kinFactor;
316
317 /*
318 G4cout << "Z= " << targetZ << " XStot= " << xSection/barn
319 << " screenZ= " << screenZ << " formF= " << formfactA
320 << " for " << particle->GetParticleName()
321 << " m= " << mass << " 1/v= " << sqrt(invbeta2)
322 << " p= " << sqrt(mom2)
323 << " x= " << x << G4endl;
324 */
325 return xSection;
326}
327
328//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
329
332 G4double cosTMax,
333 G4double elecRatio)
334{
335 temp.set(0.0,0.0,1.0);
336 CLHEP::HepRandomEngine* rndmEngineMod = G4Random::getTheEngine();
337
338 G4double formf = formfactA;
339 G4double cost1 = cosTMin;
340 G4double cost2 = cosTMax;
341 if(elecRatio > 0.0) {
342 if(rndmEngineMod->flat() <= elecRatio) {
343 formf = 0.0;
344 cost1 = std::max(cost1,cosTetMaxElec);
345 cost2 = std::max(cost2,cosTetMaxElec);
346 }
347 }
348 if(cost1 > cost2) {
349
350 G4double w1 = 1. - cost1 + screenZ;
351 G4double w2 = 1. - cost2 + screenZ;
352 G4double z1 = w1*w2/(w1 + rndmEngineMod->flat()*(w2 - w1)) - screenZ;
353
354 G4double fm = 1.0;
356 fm += formf*z1;
357 fm = 1.0/(fm*fm);
358 } else if(fNucFormfactor == fGaussianNF) {
359 fm = G4Exp(-2*formf*z1);
360 } else if(fNucFormfactor == fFlatNF) {
361 static const G4double ccoef = 0.00508/MeV;
362 G4double x = std::sqrt(2.*mom2*z1)*ccoef*2.;
363 fm = FlatFormfactor(x);
364 fm *= FlatFormfactor(x*0.6
366 }
367 G4double grej;
368 if(fMottXSection) {
370 grej = fMottXSection->RatioMottRutherfordCosT(std::sqrt(z1))*fm*fm;
371 } else {
372 grej = (1. - z1*factB + factB1*targetZ*sqrt(z1*factB)*(2. - z1))
373 *fm*fm/(1.0 + z1*factD);
374 }
375 // G4cout << "SampleSingleScattering: E= " << tkin << " z1= "
376 // << z1 << " grej= "<< grej << " mottFact= "<< fMottFactor<< G4endl;
377 if(fMottFactor*rndmEngineMod->flat() <= grej ) {
378 // exclude "false" scattering due to formfactor and spin effect
379 G4double cost = 1.0 - z1;
380 if(cost > 1.0) { cost = 1.0; }
381 else if(cost < -1.0) { cost =-1.0; }
382 G4double sint = sqrt((1.0 - cost)*(1.0 + cost));
383 //G4cout << "sint= " << sint << G4endl;
384 G4double phi = twopi*rndmEngineMod->flat();
385 temp.set(sint*cos(phi),sint*sin(phi),cost);
386 }
387 }
388 return temp;
389}
390
391//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
392
393void
395{
396 if(mass > MeV) {
397 G4double ratio = electron_mass_c2/mass;
398 G4double tau = tkin/mass;
399 G4double tmax = 2.0*electron_mass_c2*tau*(tau + 2.)/
400 (1.0 + 2.0*ratio*(tau + 1.0) + ratio*ratio);
401 cosTetMaxElec = 1.0 - std::min(cutEnergy, tmax)*electron_mass_c2/mom2;
402 } else {
403
404 G4double tmax = (particle == theElectron) ? 0.5*tkin : tkin;
405 G4double t = std::min(cutEnergy, tmax);
406 G4double mom21 = t*(t + 2.0*electron_mass_c2);
407 G4double t1 = tkin - t;
408 //G4cout <<"tkin=" <<tkin<<" tmax= "<<tmax<<" t= "
409 //<<t<< " t1= "<<t1<<" cut= "<<ecut<<G4endl;
410 if(t1 > 0.0) {
411 G4double mom22 = t1*(t1 + 2.0*mass);
412 G4double ctm = (mom2 + mom22 - mom21)*0.5/sqrt(mom2*mom22);
413 if(ctm < 1.0) { cosTetMaxElec = ctm; }
414 if(particle == theElectron && cosTetMaxElec < 0.0) {
415 cosTetMaxElec = 0.0;
416 }
417 }
418 }
419}
420
421//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
422
425{
426 return 0.0;
427}
428
429//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:180
const G4double a0
G4double G4Log(G4double x)
Definition: G4Log.hh:227
#define G4MUTEX_INITIALIZER
Definition: G4Threading.hh:85
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:221
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:116
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
#define DBL_MAX
Definition: templates.hh:62