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

#include <G4IonsSihverCrossSection.hh>

+ Inheritance diagram for G4IonsSihverCrossSection:

Public Member Functions

 G4IonsSihverCrossSection ()
 
virtual ~G4IonsSihverCrossSection ()
 
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 G4IonsSihverCrossSection.hh.

Constructor & Destructor Documentation

◆ G4IonsSihverCrossSection()

G4IonsSihverCrossSection::G4IonsSihverCrossSection ( )

Definition at line 40 of file G4IonsSihverCrossSection.cc.

41 : G4VCrossSectionDataSet("IonsSihver"), square_r0 ( (1.36*fermi) * (1.36*fermi) )
42{}

◆ ~G4IonsSihverCrossSection()

G4IonsSihverCrossSection::~G4IonsSihverCrossSection ( )
virtual

Definition at line 44 of file G4IonsSihverCrossSection.cc.

45{}

Member Function Documentation

◆ CrossSectionDescription()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 48 of file G4IonsSihverCrossSection.cc.

49{
50 outFile << "G4IonsSihverCrossSection calculates the total reaction cross\n"
51 << "section for nucleus-nucleus scattering using the Sihver\n"
52 << "parameterization. It is valid for projectiles and targets of\n"
53 << "all Z, and for all projectile energies above 100 MeV/n.\n";
54}

◆ GetElementCrossSection()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 67 of file G4IonsSihverCrossSection.cc.

69{
70 G4double xsection = 0.0;
71 G4int At = G4lrint(G4NistManager::Instance()->GetAtomicMassAmu(Z));
72
73 G4int Ap = aParticle->GetDefinition()->GetBaryonNumber();
74
75 G4Pow* g4pow = G4Pow::GetInstance();
76
77 G4double cubicrAt = g4pow->Z13(At);
78 G4double cubicrAp = g4pow->Z13(Ap);
79
80 G4double b0 = 1.581 - 0.876 * (1.0/cubicrAp + 1.0/cubicrAt);
81
82 G4double xr = cubicrAp + cubicrAt - b0 * (1.0/cubicrAp + 1.0/cubicrAt);
83 xsection = pi * square_r0 * xr * xr;
84
85 return xsection;
86}
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
G4ParticleDefinition * GetDefinition() const
static G4NistManager * Instance()
Definition: G4Pow.hh:54
static G4Pow * GetInstance()
Definition: G4Pow.cc:50
G4double Z13(G4int Z)
Definition: G4Pow.hh:110
const G4double pi
int G4lrint(double ad)
Definition: templates.hh:163

◆ IsElementApplicable()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 57 of file G4IonsSihverCrossSection.cc.

59{
60 G4int BaryonNumber = aDP->GetDefinition()->GetBaryonNumber();
61 G4double KineticEnergy = aDP->GetKineticEnergy();
62 if ( KineticEnergy / BaryonNumber >= 100*MeV && BaryonNumber > 1 ) { return true; }
63 return false;
64}
G4double GetKineticEnergy() const

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