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

#include <G4IonsKoxCrossSection.hh>

+ Inheritance diagram for G4IonsKoxCrossSection:

Public Member Functions

 G4IonsKoxCrossSection ()
 
 ~G4IonsKoxCrossSection ()
 
virtual G4bool IsElementApplicable (const G4DynamicParticle *aDP, G4int Z, const G4Material *)
 
virtual G4double GetElementCrossSection (const G4DynamicParticle *, G4int Z, const G4Material *)
 
virtual void CrossSectionDescription (std::ostream &) const
 
- Public Member Functions inherited from G4VCrossSectionDataSet
 G4VCrossSectionDataSet (const G4String &nam="")
 
virtual ~G4VCrossSectionDataSet ()
 
virtual G4bool IsElementApplicable (const G4DynamicParticle *, G4int Z, const G4Material *mat=0)
 
virtual G4bool IsIsoApplicable (const G4DynamicParticle *, G4int Z, G4int A, const G4Element *elm=0, const G4Material *mat=0)
 
G4double GetCrossSection (const G4DynamicParticle *, const G4Element *, const G4Material *mat=0)
 
G4double ComputeCrossSection (const G4DynamicParticle *, const G4Element *, const G4Material *mat=0)
 
virtual G4double GetElementCrossSection (const G4DynamicParticle *, G4int Z, const G4Material *mat=0)
 
virtual G4double GetIsoCrossSection (const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
 
virtual G4IsotopeSelectIsotope (const G4Element *, G4double kinEnergy)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void DumpPhysicsTable (const G4ParticleDefinition &)
 
virtual void CrossSectionDescription (std::ostream &) const
 
void SetVerboseLevel (G4int value)
 
G4double GetMinKinEnergy () const
 
void SetMinKinEnergy (G4double value)
 
G4double GetMaxKinEnergy () const
 
void SetMaxKinEnergy (G4double value)
 
const G4StringGetName () const
 

Additional Inherited Members

- Protected Member Functions inherited from G4VCrossSectionDataSet
void SetName (const G4String &)
 
- Protected Attributes inherited from G4VCrossSectionDataSet
G4int verboseLevel
 

Detailed Description

Definition at line 46 of file G4IonsKoxCrossSection.hh.

Constructor & Destructor Documentation

◆ G4IonsKoxCrossSection()

G4IonsKoxCrossSection::G4IonsKoxCrossSection ( )

Definition at line 41 of file G4IonsKoxCrossSection.cc.

42 : G4VCrossSectionDataSet("IonsKox"), lowerLimit ( 10*MeV ),
43 r0 ( 1.1*fermi ), rc ( 1.3*fermi )
44{}

◆ ~G4IonsKoxCrossSection()

G4IonsKoxCrossSection::~G4IonsKoxCrossSection ( )

Definition at line 46 of file G4IonsKoxCrossSection.cc.

47{}

Member Function Documentation

◆ CrossSectionDescription()

void G4IonsKoxCrossSection::CrossSectionDescription ( std::ostream &  outFile) const
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 50 of file G4IonsKoxCrossSection.cc.

51{
52 outFile << "G4IonsKoxCrossSection calculates the total reaction cross\n"
53 << "section for nucleus-nucleus scattering using the Kox\n"
54 << "parameterization. It is valid for projectiles and targets\n"
55 << "of all Z, at projectile energies up to 10 GeV/n. If the\n"
56 << "projectile energy is less than 10 MeV/n, a zero cross section\n"
57 << "is returned.\n";
58}

◆ GetElementCrossSection()

G4double G4IonsKoxCrossSection::GetElementCrossSection ( const G4DynamicParticle aParticle,
G4int  Z,
const G4Material  
)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 67 of file G4IonsKoxCrossSection.cc.

69{
70 G4double xsection = 0.0;
71
72 G4int Ap = aParticle->GetDefinition()->GetBaryonNumber();
73 G4int Zp = G4int(aParticle->GetDefinition()->GetPDGCharge() / eplus + 0.5);
74 G4double ke_per_N = aParticle->GetKineticEnergy() / Ap;
75
76 // Apply energy check, if less than lower limit then 0 value is returned
77 // if ( ke_per_N < lowerLimit ) return xsection;
78
79 G4int At = G4lrint(G4NistManager::Instance()->GetAtomicMassAmu(ZZ));
80 G4int Zt = ZZ;
81
82 G4double one_third = 1.0 / 3.0;
83
84 G4double cubicrAt = std::pow ( G4double(At) , G4double(one_third) );
85 G4double cubicrAp = std::pow ( G4double(Ap) , G4double(one_third) );
86
87 // rc divide fermi
88 G4double Bc = Zt * Zp / ( (rc/fermi) * (cubicrAp+cubicrAt) );
89
91 G4double proj_mass = aParticle->GetMass();
92 G4double proj_momentum = aParticle->GetMomentum().mag();
93
94 G4double Ecm = calEcm ( proj_mass , targ_mass , proj_momentum );
95 if( Ecm <= Bc) return xsection;
96
97 G4double Rvol = r0 * ( cubicrAp + cubicrAt );
98
99// G4double ke_per_N = aParticle->GetKineticEnergy() / Ap;
100 G4double c = calCeValue ( ke_per_N / MeV );
101
102 G4double a = 1.85;
103 G4double Rsurf = r0 * (a*cubicrAp * cubicrAt/(cubicrAp + cubicrAt) - c);
104 G4double D = 5.0 * ( At - 2 * Zt ) * Zp / ( Ap * At );
105 Rsurf = Rsurf + D * fermi; // multiply D by fermi
106
107 G4double Rint = Rvol + Rsurf;
108 xsection = pi * Rint * Rint * ( 1 - Bc / ( Ecm / MeV ) );
109
110 return xsection;
111}
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
double mag() const
G4double GetMass() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4ThreeVector GetMomentum() const
static G4NistManager * Instance()
static G4double GetNuclearMass(const G4double A, const G4double Z)
G4double GetPDGCharge() const
const G4double pi
int G4lrint(double ad)
Definition: templates.hh:163

◆ IsElementApplicable()

G4bool G4IonsKoxCrossSection::IsElementApplicable ( const G4DynamicParticle aDP,
G4int  Z,
const G4Material  
)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 60 of file G4IonsKoxCrossSection.cc.

62{
63 return (1 <= aDP->GetDefinition()->GetBaryonNumber());
64}

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