#include <G4WentzelOKandVIxSection.hh>
Definition at line 72 of file G4WentzelOKandVIxSection.hh.
◆ G4WentzelOKandVIxSection() [1/2]
G4WentzelOKandVIxSection::G4WentzelOKandVIxSection |
( |
G4bool | combined = true | ) |
|
|
explicit |
Definition at line 77 of file G4WentzelOKandVIxSection.cc.
77 :
80{
83
87
88 G4double p0 = CLHEP::electron_mass_c2*CLHEP::classic_electr_radius;
89 coeff = CLHEP::twopi*p0*p0;
91}
static G4Electron * Electron()
static G4NistManager * Instance()
static G4Positron * Positron()
static G4Pow * GetInstance()
static G4Proton * Proton()
const G4ParticleDefinition * theProton
G4NistManager * fNistManager
const G4ParticleDefinition * thePositron
const G4ParticleDefinition * theElectron
Referenced by G4WentzelOKandVIxSection(), and operator=().
◆ ~G4WentzelOKandVIxSection()
G4WentzelOKandVIxSection::~G4WentzelOKandVIxSection |
( |
| ) |
|
|
virtual |
◆ G4WentzelOKandVIxSection() [2/2]
◆ ComputeElectronCrossSection()
◆ ComputeMaxElectronScattering()
void G4WentzelOKandVIxSection::ComputeMaxElectronScattering |
( |
G4double | cut | ) |
|
|
protected |
Definition at line 392 of file G4WentzelOKandVIxSection.cc.
393{
397 G4double tmax = 2.0*electron_mass_c2*tau*(tau + 2.)/
398 (1.0 + 2.0*ratio*(tau + 1.0) + ratio*ratio);
400 } else {
401
403 G4double t = std::min(cutEnergy, tmax);
404 G4double mom21 = t*(t + 2.0*electron_mass_c2);
406
407
408 if(t1 > 0.0) {
414 }
415 }
416 }
417}
const G4ParticleDefinition * particle
Referenced by SetupTarget().
◆ ComputeNuclearCrossSection()
◆ ComputeSecondTransportMoment()
G4double G4WentzelOKandVIxSection::ComputeSecondTransportMoment |
( |
G4double | CosThetaMax | ) |
|
◆ ComputeTransportCrossSectionPerAtom()
G4double G4WentzelOKandVIxSection::ComputeTransportCrossSectionPerAtom |
( |
G4double | CosThetaMax | ) |
|
Definition at line 242 of file G4WentzelOKandVIxSection.cc.
243{
245 if(cosTMax >= 1.0) { return xSection; }
246
249
250
251 if(costm < 1.0) {
255 xSection = x2*((1.0 - 1.3333333*x + 3*x2) - fb*x*(0.6666667 - x));
256 } else {
259 xSection = xlog - x1 - fb*(x + x1 - 2*xlog);
260 }
261
262 if(xSection < 0.0) {
265 G4cout <<
"G4WentzelOKandVIxSection::ComputeTransportCrossSectionPerAtom"
266 << " scattering on e- <0"
268 G4cout <<
"cross= " << xSection
269 <<
" e(MeV)= " <<
tkin <<
" p(MeV/c)= " << sqrt(
mom2)
272 G4cout <<
" 1-costm= " << 1.0-costm <<
" screenZ= " <<
screenZ
274 }
275 xSection = 0.0;
276 }
277 }
278
279
280
281
282
283
284
285
286
287 if(cosTMax < 1.0) {
292 y = x2*((1.0 - 1.3333333*x + 3*x2) - fb*x*(0.6666667 - x));
293 } else {
296 y = xlog - x1 - fb*(x + x1 - 2*xlog);
297 }
298
299 if(y < 0.0) {
302 G4cout <<
"G4WentzelOKandVIxSection::ComputeTransportCrossSectionPerAtom"
303 << " scattering on nucleus <0"
310 }
311 y = 0.0;
312 }
314 }
316
317
318
319
320
321
322
323
324
325 return xSection;
326}
G4double G4Log(G4double x)
G4GLOB_DLL std::ostream G4cout
◆ FlatFormfactor()
◆ GetCosThetaElec()
G4double G4WentzelOKandVIxSection::GetCosThetaElec |
( |
| ) |
const |
|
inline |
◆ GetCosThetaNuc()
G4double G4WentzelOKandVIxSection::GetCosThetaNuc |
( |
| ) |
const |
|
inline |
◆ GetMomentumSquare()
G4double G4WentzelOKandVIxSection::GetMomentumSquare |
( |
| ) |
const |
|
inline |
◆ Initialise()
Definition at line 102 of file G4WentzelOKandVIxSection.cc.
104{
109
110
116
119
120
124 }
125
126
127
128
129
130}
static G4EmParameters * Instance()
G4NuclearFormfactorType NuclearFormfactorType() const
G4double FactorForAngleLimit() const
void SetupParticle(const G4ParticleDefinition *)
static G4double ScreenRSquare[100]
G4NuclearFormfactorType fNucFormfactor
const G4Material * currentMaterial
◆ InitialiseA()
void G4WentzelOKandVIxSection::InitialiseA |
( |
| ) |
|
|
protected |
Definition at line 134 of file G4WentzelOKandVIxSection.cc.
135{
136
137
141 const G4double invmev2 = 1./(CLHEP::MeV*CLHEP::MeV);
142 G4double a0 = CLHEP::electron_mass_c2/0.88534;
145
151
152 for(
G4int j=2; j<100; ++j) {
158 }
159 }
160 l.unlock();
161}
G4TemplateAutoLock< G4Mutex > G4AutoLock
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
G4double ScreeningFactor() const
static G4double ScreenRSquareElec[100]
static G4double FormFactor[100]
Referenced by Initialise().
◆ operator=()
◆ SampleSingleScattering()
Definition at line 331 of file G4WentzelOKandVIxSection.cc.
334{
335 temp.set(0.0,0.0,1.0);
336 CLHEP::HepRandomEngine* rndmEngineMod = G4Random::getTheEngine();
337
341 if(elecRatio > 0.0) {
342 if(rndmEngineMod->
flat() <= elecRatio) {
343 formf = 0.0;
346 }
347 }
348 if(cost1 > cost2) {
354
356 fm += formf*z1;
357 fm = 1.0/(fm*fm);
359 fm =
G4Exp(-2*formf*z1);
361 static const G4double ccoef = 0.00508/CLHEP::MeV;
365 }
366
370 grej =
fMottXSection->RatioMottRutherfordCosT(std::sqrt(z1))*fm;
371 } else {
373 *fm/(1.0 + z1*
factD);
374 }
376
378 if(cost > 1.0) { cost = 1.0; }
379 else if(cost < -1.0) { cost =-1.0; }
380 G4double sint = sqrt((1.0 - cost)*(1.0 + cost));
381
383 temp.set(sint*cos(phi),sint*sin(phi),cost);
384 }
385 }
387}
G4double FlatFormfactor(G4double x)
◆ SetTargetMass()
void G4WentzelOKandVIxSection::SetTargetMass |
( |
G4double | value | ) |
|
|
inline |
◆ SetupKinematic()
◆ SetupParticle()
◆ SetupTarget()
Definition at line 200 of file G4WentzelOKandVIxSection.cc.
201{
206 G4double massT = (1 == Z) ? CLHEP::proton_mass_c2 :
209
213 }
214
215 if(1 == Z) {
217 }
else if(
mass > MeV) {
220 } else {
222 screenZ = std::min(Z*1.13,(1.13 +3.76*Z*Z
225 }
227 cosTetMaxNuc2 = 0.0;
228 }
230
233 }
234
235
236 return cosTetMaxNuc2;
237}
void SetTargetMass(G4double value)
void ComputeMaxElectronScattering(G4double cut)
◆ charge3
G4double G4WentzelOKandVIxSection::charge3 = 0.0 |
|
protected |
◆ chargeSquare
G4double G4WentzelOKandVIxSection::chargeSquare = 0.0 |
|
protected |
◆ coeff
G4double G4WentzelOKandVIxSection::coeff |
|
protected |
◆ cosTetMaxElec
G4double G4WentzelOKandVIxSection::cosTetMaxElec = 1.0 |
|
protected |
◆ cosTetMaxNuc
G4double G4WentzelOKandVIxSection::cosTetMaxNuc = 1.0 |
|
protected |
◆ cosThetaMax
G4double G4WentzelOKandVIxSection::cosThetaMax = 1.0 |
|
protected |
◆ currentMaterial
const G4Material* G4WentzelOKandVIxSection::currentMaterial = nullptr |
|
protected |
◆ ecut
◆ etag
◆ factB
G4double G4WentzelOKandVIxSection::factB = 0.0 |
|
protected |
◆ factD
G4double G4WentzelOKandVIxSection::factD = 0.0 |
|
protected |
◆ factorA2
G4double G4WentzelOKandVIxSection::factorA2 = 0.0 |
|
protected |
◆ fG4pow
G4Pow* G4WentzelOKandVIxSection::fG4pow |
|
protected |
◆ fMottFactor
G4double G4WentzelOKandVIxSection::fMottFactor = 1.0 |
|
protected |
◆ fMottXSection
◆ fNistManager
◆ fNucFormfactor
◆ formfactA
G4double G4WentzelOKandVIxSection::formfactA = 0.0 |
|
protected |
◆ FormFactor
G4double G4WentzelOKandVIxSection::FormFactor = {0.0} |
|
staticprotected |
◆ gam0pcmp
G4double G4WentzelOKandVIxSection::gam0pcmp = 1.0 |
|
protected |
◆ invbeta2
G4double G4WentzelOKandVIxSection::invbeta2 = 1.0 |
|
protected |
◆ isCombined
G4bool G4WentzelOKandVIxSection::isCombined |
|
protected |
◆ kinFactor
G4double G4WentzelOKandVIxSection::kinFactor = 1.0 |
|
protected |
◆ mass
G4double G4WentzelOKandVIxSection::mass = 0.0 |
|
protected |
◆ mom2
G4double G4WentzelOKandVIxSection::mom2 = 0.0 |
|
protected |
◆ momCM2
G4double G4WentzelOKandVIxSection::momCM2 = 0.0 |
|
protected |
◆ nwarnings
G4int G4WentzelOKandVIxSection::nwarnings = 0 |
|
protected |
◆ particle
◆ pcmp2
G4double G4WentzelOKandVIxSection::pcmp2 = 1.0 |
|
protected |
◆ ScreenRSquare
G4double G4WentzelOKandVIxSection::ScreenRSquare = {0.0} |
|
staticprotected |
◆ ScreenRSquareElec
G4double G4WentzelOKandVIxSection::ScreenRSquareElec = {0.0} |
|
staticprotected |
◆ screenZ
G4double G4WentzelOKandVIxSection::screenZ = 0.0 |
|
protected |
◆ spin
G4double G4WentzelOKandVIxSection::spin = 0.0 |
|
protected |
◆ targetMass
G4double G4WentzelOKandVIxSection::targetMass |
|
protected |
◆ targetZ
G4int G4WentzelOKandVIxSection::targetZ = 0 |
|
protected |
◆ temp
◆ theElectron
◆ thePositron
◆ theProton
◆ tkin
G4double G4WentzelOKandVIxSection::tkin = 0.0 |
|
protected |
The documentation for this class was generated from the following files: