#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 |
Definition at line 51 of file G4IonsShenCrossSection.hh.
G4IonsShenCrossSection::G4IonsShenCrossSection | ( | ) |
Definition at line 43 of file G4IonsShenCrossSection.cc.
00044 : G4VCrossSectionDataSet("IonsShen"), 00045 upperLimit( 10*GeV ), lowerLimit( 10*MeV ), r0 ( 1.1 ) 00046 {}
G4IonsShenCrossSection::~G4IonsShenCrossSection | ( | ) | [virtual] |
void G4IonsShenCrossSection::CrossSectionDescription | ( | std::ostream & | ) | const [virtual] |
Reimplemented from G4VCrossSectionDataSet.
Definition at line 52 of file G4IonsShenCrossSection.cc.
00053 { 00054 outFile << "G4IonsShenCrossSection calculates the total reaction cross\n" 00055 << "section for nucleus-nucleus scattering using the Shen\n" 00056 << "parameterization. It is valid for projectiles and targets of\n" 00057 << "all Z, and projectile energies up to 1 TeV/n. Above 10 GeV/n" 00058 << "the cross section is constant. Below 10 MeV/n zero cross\n" 00059 << "is returned.\n"; 00060 }
G4double G4IonsShenCrossSection::GetElementCrossSection | ( | const G4DynamicParticle * | , | |
G4int | Z, | |||
const G4Material * | ||||
) | [virtual] |
Reimplemented from G4VCrossSectionDataSet.
Definition at line 69 of file G4IonsShenCrossSection.cc.
References G4lrint(), GetIsoCrossSection(), and G4NistManager::Instance().
00072 { 00073 G4int A = G4lrint(G4NistManager::Instance()->GetAtomicMassAmu(Z)); 00074 return GetIsoCrossSection(aParticle, Z, A); 00075 }
G4double G4IonsShenCrossSection::GetIsoCrossSection | ( | const G4DynamicParticle * | , | |
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.
References G4Pow::A13(), G4lrint(), G4ParticleDefinition::GetBaryonNumber(), G4DynamicParticle::GetDefinition(), G4Pow::GetInstance(), G4DynamicParticle::GetKineticEnergy(), G4DynamicParticle::GetMass(), G4DynamicParticle::GetMomentum(), G4NucleiProperties::GetNuclearMass(), G4ParticleDefinition::GetPDGCharge(), G4INCL::Math::pi, and G4Pow::Z13().
Referenced by GetElementCrossSection().
00083 { 00084 G4double xsection = 0.0; 00085 00086 G4int Ap = aParticle->GetDefinition()->GetBaryonNumber(); 00087 G4int Zp = G4lrint(aParticle->GetDefinition()->GetPDGCharge()/eplus); 00088 G4double ke_per_N = aParticle->GetKineticEnergy() / Ap; 00089 if ( ke_per_N > upperLimit ) { ke_per_N = upperLimit; } 00090 00091 // Apply energy check, if less than lower limit then 0 value is returned 00092 //if ( ke_per_N < lowerLimit ) { return xsection; } 00093 00094 G4Pow* g4pow = G4Pow::GetInstance(); 00095 00096 G4double cubicrAt = g4pow->Z13(At); 00097 G4double cubicrAp = g4pow->Z13(Ap); 00098 00099 G4double Rt = 1.12 * cubicrAt - 0.94 * ( 1.0 / cubicrAt ); 00100 G4double Rp = 1.12 * cubicrAp - 0.94 * ( 1.0 / cubicrAp ); 00101 00102 G4double r = Rt + Rp + 3.2; // in fm 00103 G4double b = 1.0; // in MeV/fm 00104 G4double targ_mass = G4NucleiProperties::GetNuclearMass(At, Zt); 00105 00106 G4double proj_mass = aParticle->GetMass(); 00107 G4double proj_momentum = aParticle->GetMomentum().mag(); 00108 00109 G4double Ecm = calEcmValue (proj_mass, targ_mass, proj_momentum); 00110 00111 G4double B = 1.44 * Zt * Zp / r - b * Rt * Rp / ( Rt + Rp ); 00112 if(Ecm <= B) { return xsection; } 00113 00114 G4double c = calCeValue ( ke_per_N / MeV ); 00115 00116 G4double R1 = r0 * (cubicrAt + cubicrAp + 1.85*cubicrAt*cubicrAp/(cubicrAt + cubicrAp) - c); 00117 00118 G4double R2 = 1.0 * ( At - 2 * Zt ) * Zp / ( Ap * At ); 00119 00120 00121 G4double R3 = (0.176 / g4pow->A13(Ecm)) * cubicrAt * cubicrAp /(cubicrAt + cubicrAp); 00122 00123 G4double R = R1 + R2 + R3; 00124 00125 xsection = 10 * pi * R * R * ( 1 - B / Ecm ); 00126 xsection = xsection * millibarn; // mulitply xsection by millibarn 00127 00128 return xsection; 00129 }
G4bool G4IonsShenCrossSection::IsElementApplicable | ( | const G4DynamicParticle * | aDP, | |
G4int | Z, | |||
const G4Material * | ||||
) | [virtual] |
Reimplemented from G4VCrossSectionDataSet.
Definition at line 62 of file G4IonsShenCrossSection.cc.
References G4ParticleDefinition::GetBaryonNumber(), and G4DynamicParticle::GetDefinition().
00064 { 00065 return (1 <= aDP->GetDefinition()->GetBaryonNumber()); 00066 }