#include <G4IonsKoxCrossSection.hh>
|
| 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 |
|
| 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 G4Isotope * | SelectIsotope (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 G4String & | GetName () const |
|
Definition at line 46 of file G4IonsKoxCrossSection.hh.
◆ G4IonsKoxCrossSection()
G4IonsKoxCrossSection::G4IonsKoxCrossSection |
( |
| ) |
|
◆ ~G4IonsKoxCrossSection()
G4IonsKoxCrossSection::~G4IonsKoxCrossSection |
( |
| ) |
|
◆ 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()
Reimplemented from G4VCrossSectionDataSet.
Definition at line 67 of file G4IonsKoxCrossSection.cc.
69{
71
75
76
77
78
81
83
86
87
88 G4double Bc = Zt * Zp / ( (rc/fermi) * (cubicrAp+cubicrAt) );
89
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
100 G4double c = calCeValue ( ke_per_N / MeV );
101
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;
106
108 xsection =
pi * Rint * Rint * ( 1 - Bc / ( Ecm / MeV ) );
109
110 return xsection;
111}
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4ThreeVector GetMomentum() const
static G4NistManager * Instance()
static G4double GetNuclearMass(const G4double A, const G4double Z)
G4double GetPDGCharge() const
G4int GetBaryonNumber() const
◆ IsElementApplicable()
The documentation for this class was generated from the following files: