Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4WentzelVIRelXSection Class Reference

#include <G4WentzelVIRelXSection.hh>

+ Inheritance diagram for G4WentzelVIRelXSection:

Public Member Functions

 G4WentzelVIRelXSection ()
 
 ~G4WentzelVIRelXSection () override
 
G4double SetupKinematic (G4double kinEnergy, const G4Material *mat) override
 
G4WentzelVIRelXSectionoperator= (const G4WentzelVIRelXSection &right)=delete
 
 G4WentzelVIRelXSection (const G4WentzelVIRelXSection &)=delete
 
- Public Member Functions inherited from G4WentzelOKandVIxSection
 G4WentzelOKandVIxSection (G4bool comb=true)
 
virtual ~G4WentzelOKandVIxSection ()
 
void Initialise (const G4ParticleDefinition *, G4double CosThetaLim)
 
void SetupParticle (const G4ParticleDefinition *)
 
virtual G4double SetupKinematic (G4double kinEnergy, const G4Material *mat)
 
G4double SetupTarget (G4int Z, G4double cut)
 
G4double ComputeTransportCrossSectionPerAtom (G4double CosThetaMax)
 
G4ThreeVectorSampleSingleScattering (G4double CosThetaMin, G4double CosThetaMax, G4double elecRatio)
 
G4double ComputeSecondTransportMoment (G4double CosThetaMax)
 
G4double ComputeNuclearCrossSection (G4double CosThetaMin, G4double CosThetaMax)
 
G4double ComputeElectronCrossSection (G4double CosThetaMin, G4double CosThetaMax)
 
void SetTargetMass (G4double value)
 
G4double GetMomentumSquare () const
 
G4double GetCosThetaNuc () const
 
G4double GetCosThetaElec () const
 
G4WentzelOKandVIxSectionoperator= (const G4WentzelOKandVIxSection &right)=delete
 
 G4WentzelOKandVIxSection (const G4WentzelOKandVIxSection &)=delete
 

Additional Inherited Members

- Protected Member Functions inherited from G4WentzelOKandVIxSection
void ComputeMaxElectronScattering (G4double cut)
 
void InitialiseA ()
 
G4double FlatFormfactor (G4double x)
 
- Protected Attributes inherited from G4WentzelOKandVIxSection
const G4ParticleDefinitiontheProton
 
const G4ParticleDefinitiontheElectron
 
const G4ParticleDefinitionthePositron
 
const G4ParticleDefinitionparticle = nullptr
 
const G4MaterialcurrentMaterial = nullptr
 
G4NistManagerfNistManager
 
G4PowfG4pow
 
G4ScreeningMottCrossSectionfMottXSection = nullptr
 
G4ThreeVector temp
 
G4double coeff
 
G4double cosTetMaxElec = 1.0
 
G4double cosTetMaxNuc = 1.0
 
G4double cosThetaMax = -1.0
 
G4double chargeSquare = 0.0
 
G4double charge3 = 0.0
 
G4double spin = 0.0
 
G4double mass = 0.0
 
G4double tkin = 0.0
 
G4double mom2 = 0.0
 
G4double momCM2 = 0.0
 
G4double invbeta2 = 1.0
 
G4double kinFactor = 1.0
 
G4double etag = DBL_MAX
 
G4double ecut = DBL_MAX
 
G4double targetMass
 
G4double screenZ = 0.0
 
G4double formfactA = 0.0
 
G4double factorA2 = 0.0
 
G4double factB = 0.0
 
G4double factD = 0.0
 
G4double fMottFactor = 1.0
 
G4double gam0pcmp = 1.0
 
G4double pcmp2 = 1.0
 
G4int targetZ = 0
 
G4int nwarnings = 0
 
G4NuclearFormfactorType fNucFormfactor = fExponentialNF
 
G4bool isCombined
 
- Static Protected Attributes inherited from G4WentzelOKandVIxSection
static G4double ScreenRSquareElec [100] = {0.0}
 
static G4double ScreenRSquare [100] = {0.0}
 
static G4double FormFactor [100] = {0.0}
 

Detailed Description

Definition at line 63 of file G4WentzelVIRelXSection.hh.

Constructor & Destructor Documentation

◆ G4WentzelVIRelXSection() [1/2]

G4WentzelVIRelXSection::G4WentzelVIRelXSection ( )
explicitdefault

◆ ~G4WentzelVIRelXSection()

G4WentzelVIRelXSection::~G4WentzelVIRelXSection ( )
overridedefault

◆ G4WentzelVIRelXSection() [2/2]

G4WentzelVIRelXSection::G4WentzelVIRelXSection ( const G4WentzelVIRelXSection )
delete

Member Function Documentation

◆ operator=()

G4WentzelVIRelXSection & G4WentzelVIRelXSection::operator= ( const G4WentzelVIRelXSection right)
delete

◆ SetupKinematic()

G4double G4WentzelVIRelXSection::SetupKinematic ( G4double  kinEnergy,
const G4Material mat 
)
overridevirtual

Reimplemented from G4WentzelOKandVIxSection.

Definition at line 60 of file G4WentzelVIRelXSection.cc.

62{
63 if(kinEnergy != tkin || mat != currentMaterial) {
64
65 currentMaterial = mat;
66 tkin = kinEnergy;
67 G4double momLab2 = tkin*(tkin + 2.0*mass);
68
69 G4double etot = tkin + mass;
70 G4double ptot = std::sqrt(momLab2);
71 G4double m12 = mass*mass;
72
73 // relativistic reduced mass from publucation
74 // A.P. Martynenko, R.N. Faustov, Teoret. mat. Fiz. 64 (1985) 179
75
76 //incident particle & target nucleus
77 G4double Ecm = std::sqrt(m12 + targetMass*targetMass + 2.0*etot*targetMass);
78 G4double mu_rel = mass*targetMass/Ecm;
79 G4double momCM = ptot*targetMass/Ecm;
80 // relative system
81 mom2 = momCM*momCM;
82 invbeta2 = 1.0 + mu_rel*mu_rel/mom2;
83
85 factD = std::sqrt(mom2)/targetMass;
87 std::max(cosThetaMax, 1.-factorA2*mat->GetIonisation()->GetInvA23()/mom2)
89 }
90 return cosTetMaxNuc;
91}
double G4double
Definition: G4Types.hh:83
G4double GetInvA23() const
G4IonisParamMat * GetIonisation() const
Definition: G4Material.hh:221

Referenced by G4hCoulombScatteringModel::ComputeCrossSectionPerAtom(), and G4hCoulombScatteringModel::SampleSecondaries().


The documentation for this class was generated from the following files: