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

#include <G4IonsShenCrossSection.hh>

+ Inheritance diagram for G4IonsShenCrossSection:

Public Member Functions

 G4IonsShenCrossSection ()
 
virtual ~G4IonsShenCrossSection ()
 
virtual G4bool IsElementApplicable (const G4DynamicParticle *aDP, G4int Z, const G4Material *)
 
virtual G4double GetElementCrossSection (const G4DynamicParticle *, G4int Z, const G4Material *)
 
virtual G4double GetIsoCrossSection (const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
 
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 51 of file G4IonsShenCrossSection.hh.

Constructor & Destructor Documentation

◆ G4IonsShenCrossSection()

G4IonsShenCrossSection::G4IonsShenCrossSection ( )

Definition at line 43 of file G4IonsShenCrossSection.cc.

44 : G4VCrossSectionDataSet("IonsShen"),
45 upperLimit( 10*GeV ), lowerLimit( 10*MeV ), r0 ( 1.1 )
46{}

◆ ~G4IonsShenCrossSection()

G4IonsShenCrossSection::~G4IonsShenCrossSection ( )
virtual

Definition at line 48 of file G4IonsShenCrossSection.cc.

49{}

Member Function Documentation

◆ CrossSectionDescription()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 52 of file G4IonsShenCrossSection.cc.

53{
54 outFile << "G4IonsShenCrossSection calculates the total reaction cross\n"
55 << "section for nucleus-nucleus scattering using the Shen\n"
56 << "parameterization. It is valid for projectiles and targets of\n"
57 << "all Z, and projectile energies up to 1 TeV/n. Above 10 GeV/n"
58 << "the cross section is constant. Below 10 MeV/n zero cross\n"
59 << "is returned.\n";
60}

◆ GetElementCrossSection()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 69 of file G4IonsShenCrossSection.cc.

72{
73 G4int A = G4lrint(G4NistManager::Instance()->GetAtomicMassAmu(Z));
74 return GetIsoCrossSection(aParticle, Z, A);
75}
int G4int
Definition: G4Types.hh:66
virtual G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
static G4NistManager * Instance()
int G4lrint(double ad)
Definition: templates.hh:163

◆ GetIsoCrossSection()

G4double G4IonsShenCrossSection::GetIsoCrossSection ( const G4DynamicParticle aParticle,
G4int  Z,
G4int  A,
const G4Isotope iso = 0,
const G4Element elm = 0,
const G4Material mat = 0 
)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 77 of file G4IonsShenCrossSection.cc.

83{
84 G4double xsection = 0.0;
85
86 G4int Ap = aParticle->GetDefinition()->GetBaryonNumber();
87 G4int Zp = G4lrint(aParticle->GetDefinition()->GetPDGCharge()/eplus);
88 G4double ke_per_N = aParticle->GetKineticEnergy() / Ap;
89 if ( ke_per_N > upperLimit ) { ke_per_N = upperLimit; }
90
91 // Apply energy check, if less than lower limit then 0 value is returned
92 //if ( ke_per_N < lowerLimit ) { return xsection; }
93
94 G4Pow* g4pow = G4Pow::GetInstance();
95
96 G4double cubicrAt = g4pow->Z13(At);
97 G4double cubicrAp = g4pow->Z13(Ap);
98
99 G4double Rt = 1.12 * cubicrAt - 0.94 * ( 1.0 / cubicrAt );
100 G4double Rp = 1.12 * cubicrAp - 0.94 * ( 1.0 / cubicrAp );
101
102 G4double r = Rt + Rp + 3.2; // in fm
103 G4double b = 1.0; // in MeV/fm
105
106 G4double proj_mass = aParticle->GetMass();
107 G4double proj_momentum = aParticle->GetMomentum().mag();
108
109 G4double Ecm = calEcmValue (proj_mass, targ_mass, proj_momentum);
110
111 G4double B = 1.44 * Zt * Zp / r - b * Rt * Rp / ( Rt + Rp );
112 if(Ecm <= B) { return xsection; }
113
114 G4double c = calCeValue ( ke_per_N / MeV );
115
116 G4double R1 = r0 * (cubicrAt + cubicrAp + 1.85*cubicrAt*cubicrAp/(cubicrAt + cubicrAp) - c);
117
118 G4double R2 = 1.0 * ( At - 2 * Zt ) * Zp / ( Ap * At );
119
120
121 G4double R3 = (0.176 / g4pow->A13(Ecm)) * cubicrAt * cubicrAp /(cubicrAt + cubicrAp);
122
123 G4double R = R1 + R2 + R3;
124
125 xsection = 10 * pi * R * R * ( 1 - B / Ecm );
126 xsection = xsection * millibarn; // mulitply xsection by millibarn
127
128 return xsection;
129}
double G4double
Definition: G4Types.hh:64
double mag() const
G4double GetMass() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4ThreeVector GetMomentum() const
static G4double GetNuclearMass(const G4double A, const G4double Z)
G4double GetPDGCharge() const
Definition: G4Pow.hh:54
static G4Pow * GetInstance()
Definition: G4Pow.cc:50
G4double A13(G4double A)
Definition: G4Pow.hh:115
G4double Z13(G4int Z)
Definition: G4Pow.hh:110
const G4double pi

Referenced by GetElementCrossSection().

◆ IsElementApplicable()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 62 of file G4IonsShenCrossSection.cc.

64{
65 return (1 <= aDP->GetDefinition()->GetBaryonNumber());
66}

Referenced by G4GeneralSpaceNNCrossSection::GetElementCrossSection().


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