Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4eBremParametrizedModel.cc File Reference
#include "G4eBremParametrizedModel.hh"
#include "G4PhysicalConstants.hh"
#include "G4SystemOfUnits.hh"
#include "G4Electron.hh"
#include "G4Positron.hh"
#include "G4Gamma.hh"
#include "Randomize.hh"
#include "G4Material.hh"
#include "G4Element.hh"
#include "G4ElementVector.hh"
#include "G4ProductionCutsTable.hh"
#include "G4ParticleChangeForLoss.hh"
#include "G4LossTableManager.hh"
#include "G4ModifiedTsai.hh"

Go to the source code of this file.

Functions

G4double ScreenFunction1 (G4double ScreenVariable)
 
G4double ScreenFunction2 (G4double ScreenVariable)
 
G4double ComputeParametrizedDXSectionPerAtom (G4double kineticEnergy, G4double gammaEnergy, G4double Z)
 

Function Documentation

◆ ComputeParametrizedDXSectionPerAtom()

G4double ComputeParametrizedDXSectionPerAtom ( G4double  kineticEnergy,
G4double  gammaEnergy,
G4double  Z 
)

Definition at line 357 of file G4eBremParametrizedModel.cc.

357 {
358 G4double lnZ = std::log(Z); // 3.*(anElement->GetIonisation()->GetlogZ3());
359 G4double FZ = lnZ* (4.- 0.55*lnZ);
360 G4double ZZ = std::pow (Z*(Z+1.),1./3.); // anElement->GetIonisation()->GetZZ3();
361 G4double Z3 = std::pow (Z,1./3.); // (anElement->GetIonisation()->GetZ3())
362
363 G4double totalEnergy = kineticEnergy + electron_mass_c2;
364
365 // G4double x, epsil, greject, migdal, grejmax, q;
366 G4double epsil, greject;
367 G4double U = log(kineticEnergy/electron_mass_c2);
368 G4double U2 = U*U;
369
370 // precalculated parameters
371 G4double ah, bh;
372
373if (kineticEnergy > tlow) {
374
375 G4double ah1 = ah10 + ZZ* (ah11 + ZZ* ah12);
376 G4double ah2 = ah20 + ZZ* (ah21 + ZZ* ah22);
377 G4double ah3 = ah30 + ZZ* (ah31 + ZZ* ah32);
378
379 G4double bh1 = bh10 + ZZ* (bh11 + ZZ* bh12);
380 G4double bh2 = bh20 + ZZ* (bh21 + ZZ* bh22);
381 G4double bh3 = bh30 + ZZ* (bh31 + ZZ* bh32);
382
383 ah = 1. + (ah1*U2 + ah2*U + ah3) / (U2*U);
384 bh = 0.75 + (bh1*U2 + bh2*U + bh3) / (U2*U);
385
386 // limit of the screening variable
387 G4double screenfac =
388 136.*electron_mass_c2/(Z3*totalEnergy);
389
390 epsil = gammaEnergy/totalEnergy; // epsil = x*kineticEnergy/totalEnergy;
391 G4double screenvar = screenfac*epsil/(1.0-epsil);
392 G4double F1 = max(ScreenFunction1(screenvar) - FZ ,0.);
393 G4double F2 = max(ScreenFunction2(screenvar) - FZ ,0.);
394
395
396 greject = (F1 - epsil* (ah*F1 - bh*epsil*F2))/8.; // 1./(42.392 - FZ);
397
398 std::cout << " yy = "<<epsil<<std::endl;
399 std::cout << " F1/(...) "<<F1/(42.392 - FZ)<<std::endl;
400 std::cout << " F2/(...) "<<F2/(42.392 - FZ)<<std::endl;
401 std::cout << " (42.392 - FZ) " << (42.392 - FZ) <<std::endl;
402
403 } else {
404
405 G4double al0 = al00 + ZZ* (al01 + ZZ* al02);
406 G4double al1 = al10 + ZZ* (al11 + ZZ* al12);
407 G4double al2 = al20 + ZZ* (al21 + ZZ* al22);
408
409 G4double bl0 = bl00 + ZZ* (bl01 + ZZ* bl02);
410 G4double bl1 = bl10 + ZZ* (bl11 + ZZ* bl12);
411 G4double bl2 = bl20 + ZZ* (bl21 + ZZ* bl22);
412
413 ah = al0 + al1*U + al2*U2;
414 bh = bl0 + bl1*U + bl2*U2;
415
416 G4double x=gammaEnergy/kineticEnergy;
417 greject=(1. + x* (ah + bh*x));
418
419 /*
420 // Compute the maximum of the rejection function
421 grejmax = max(1. + xmin* (ah + bh*xmin), 1.+ah+bh);
422 G4double xm = -ah/(2.*bh);
423 if ( xmin < xm && xm < xmax) grejmax = max(grejmax, 1.+ xm* (ah + bh*xm));
424 */
425 }
426
427 return greject;
428}
double G4double
Definition: G4Types.hh:64
G4double ScreenFunction2(G4double ScreenVariable)
G4double ScreenFunction1(G4double ScreenVariable)

◆ ScreenFunction1()

G4double ScreenFunction1 ( G4double  ScreenVariable)

Definition at line 323 of file G4eBremParametrizedModel.cc.

327{
328 G4double screenVal;
329
330 if (ScreenVariable > 1.)
331 screenVal = 42.24 - 8.368*std::log(ScreenVariable+0.952);
332 else
333 screenVal = 42.392 - ScreenVariable* (7.796 - 1.961*ScreenVariable);
334
335 return screenVal;
336}

Referenced by ComputeParametrizedDXSectionPerAtom().

◆ ScreenFunction2()

G4double ScreenFunction2 ( G4double  ScreenVariable)

Definition at line 340 of file G4eBremParametrizedModel.cc.

344{
345 G4double screenVal;
346
347 if (ScreenVariable > 1.)
348 screenVal = 42.24 - 8.368*std::log(ScreenVariable+0.952);
349 else
350 screenVal = 41.734 - ScreenVariable* (6.484 - 1.250*ScreenVariable);
351
352 return screenVal;
353}

Referenced by ComputeParametrizedDXSectionPerAtom().