Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions
G4IonsKoxCrossSection Class Reference

#include <G4IonsKoxCrossSection.hh>

Inheritance diagram for G4IonsKoxCrossSection:
G4VCrossSectionDataSet

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 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 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 G4int GetVerboseLevel () const
 
virtual 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 ( )

Definition at line 41 of file G4IonsKoxCrossSection.cc.

42  : G4VCrossSectionDataSet("IonsKox"),
43 // lowerLimit ( 10*MeV ),
44  r0 ( 1.1*fermi ), rc ( 1.3*fermi )
45 {}
G4VCrossSectionDataSet(const G4String &nam="")
G4IonsKoxCrossSection::~G4IonsKoxCrossSection ( )

Definition at line 47 of file G4IonsKoxCrossSection.cc.

48 {}

Member Function Documentation

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 51 of file G4IonsKoxCrossSection.cc.

52 {
53  outFile << "G4IonsKoxCrossSection calculates the total reaction cross\n"
54  << "section for nucleus-nucleus scattering using the Kox\n"
55  << "parameterization. It is valid for projectiles and targets\n"
56  << "of all Z, at projectile energies up to 10 GeV/n. If the\n"
57  << "projectile energy is less than 10 MeV/n, a zero cross section\n"
58  << "is returned.\n";
59 }
std::ofstream outFile
Definition: GammaRayTel.cc:68
G4double G4IonsKoxCrossSection::GetElementCrossSection ( const G4DynamicParticle aParticle,
G4int  Z,
const G4Material  
)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 68 of file G4IonsKoxCrossSection.cc.

References test::a, test::c, python.hepunit::eplus, python.hepunit::fermi, G4lrint(), G4ParticleDefinition::GetBaryonNumber(), G4DynamicParticle::GetDefinition(), G4DynamicParticle::GetKineticEnergy(), G4DynamicParticle::GetMass(), G4DynamicParticle::GetMomentum(), G4NucleiProperties::GetNuclearMass(), G4ParticleDefinition::GetPDGCharge(), G4NistManager::Instance(), CLHEP::Hep3Vector::mag(), python.hepunit::MeV, and python.hepunit::pi.

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 61 of file G4IonsKoxCrossSection.cc.

References G4ParticleDefinition::GetBaryonNumber(), and G4DynamicParticle::GetDefinition().

63 {
64  return (1 <= aDP->GetDefinition()->GetBaryonNumber());
65 }
G4ParticleDefinition * GetDefinition() const

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