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

Detailed Description

Definition at line 46 of file G4IonsKoxCrossSection.hh.


Constructor & Destructor Documentation

G4IonsKoxCrossSection::G4IonsKoxCrossSection (  ) 

Definition at line 41 of file G4IonsKoxCrossSection.cc.

00042   : G4VCrossSectionDataSet("IonsKox"), lowerLimit ( 10*MeV ), 
00043     r0 ( 1.1*fermi ), rc ( 1.3*fermi )
00044 {}

G4IonsKoxCrossSection::~G4IonsKoxCrossSection (  ) 

Definition at line 46 of file G4IonsKoxCrossSection.cc.

00047 {}


Member Function Documentation

void G4IonsKoxCrossSection::CrossSectionDescription ( std::ostream &   )  const [virtual]

Reimplemented from G4VCrossSectionDataSet.

Definition at line 50 of file G4IonsKoxCrossSection.cc.

00051 {
00052   outFile << "G4IonsKoxCrossSection calculates the total reaction cross\n"
00053           << "section for nucleus-nucleus scattering using the Kox\n"
00054           << "parameterization.  It is valid for projectiles and targets\n"
00055           << "of all Z, at projectile energies up to 10 GeV/n.  If the\n"
00056           << "projectile energy is less than 10 MeV/n, a zero cross section\n"
00057           << "is returned.\n";
00058 }

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 67 of file G4IonsKoxCrossSection.cc.

References G4lrint(), G4ParticleDefinition::GetBaryonNumber(), G4DynamicParticle::GetDefinition(), G4DynamicParticle::GetKineticEnergy(), G4DynamicParticle::GetMass(), G4DynamicParticle::GetMomentum(), G4NucleiProperties::GetNuclearMass(), G4ParticleDefinition::GetPDGCharge(), G4NistManager::Instance(), and G4INCL::Math::pi.

00069 {
00070    G4double xsection = 0.0;
00071 
00072    G4int Ap = aParticle->GetDefinition()->GetBaryonNumber();
00073    G4int Zp = G4int(aParticle->GetDefinition()->GetPDGCharge() / eplus + 0.5); 
00074    G4double ke_per_N = aParticle->GetKineticEnergy() / Ap;
00075 
00076    // Apply energy check, if less than lower limit then 0 value is returned
00077    // if (  ke_per_N < lowerLimit ) return xsection;
00078 
00079    G4int At = G4lrint(G4NistManager::Instance()->GetAtomicMassAmu(ZZ));
00080    G4int Zt = ZZ;  
00081 
00082    G4double one_third = 1.0 / 3.0;
00083 
00084    G4double cubicrAt = std::pow ( G4double(At) , G4double(one_third) );  
00085    G4double cubicrAp = std::pow ( G4double(Ap) , G4double(one_third) );  
00086 
00087    // rc divide fermi
00088    G4double Bc = Zt * Zp / ( (rc/fermi) * (cubicrAp+cubicrAt) );
00089 
00090    G4double targ_mass = G4NucleiProperties::GetNuclearMass(At, Zt);
00091    G4double proj_mass = aParticle->GetMass();
00092    G4double proj_momentum = aParticle->GetMomentum().mag(); 
00093 
00094    G4double Ecm = calEcm ( proj_mass , targ_mass , proj_momentum ); 
00095    if( Ecm <= Bc) return xsection;
00096 
00097    G4double Rvol = r0 * (  cubicrAp + cubicrAt );
00098 
00099 //   G4double ke_per_N = aParticle->GetKineticEnergy() / Ap;
00100    G4double c = calCeValue ( ke_per_N / MeV  );  
00101 
00102    G4double a = 1.85;
00103    G4double Rsurf = r0 * (a*cubicrAp * cubicrAt/(cubicrAp + cubicrAt) - c); 
00104    G4double D = 5.0 * ( At - 2 * Zt ) * Zp / ( Ap * At );
00105    Rsurf = Rsurf + D * fermi;  // multiply D by fermi 
00106 
00107    G4double Rint = Rvol + Rsurf;
00108    xsection = pi * Rint * Rint * ( 1 - Bc / ( Ecm / MeV ) );
00109   
00110    return xsection; 
00111 }

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 60 of file G4IonsKoxCrossSection.cc.

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

00062 {
00063   return (1 <= aDP->GetDefinition()->GetBaryonNumber());
00064 }


The documentation for this class was generated from the following files:
Generated on Mon May 27 17:52:18 2013 for Geant4 by  doxygen 1.4.7