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

#include <G4IonsShenCrossSection.hh>

Inheritance diagram for G4IonsShenCrossSection:
G4VCrossSectionDataSet

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 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 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 51 of file G4IonsShenCrossSection.hh.

Constructor & Destructor Documentation

G4IonsShenCrossSection::G4IonsShenCrossSection ( )

Definition at line 43 of file G4IonsShenCrossSection.cc.

44  : G4VCrossSectionDataSet("IonsShen"),
45  upperLimit( 10*GeV ),
46 // lowerLimit( 10*MeV ),
47  r0 ( 1.1 )
48 {}
G4VCrossSectionDataSet(const G4String &nam="")
G4IonsShenCrossSection::~G4IonsShenCrossSection ( )
virtual

Definition at line 50 of file G4IonsShenCrossSection.cc.

51 {}

Member Function Documentation

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 54 of file G4IonsShenCrossSection.cc.

55 {
56  outFile << "G4IonsShenCrossSection calculates the total reaction cross\n"
57  << "section for nucleus-nucleus scattering using the Shen\n"
58  << "parameterization. It is valid for projectiles and targets of\n"
59  << "all Z, and projectile energies up to 1 TeV/n. Above 10 GeV/n"
60  << "the cross section is constant. Below 10 MeV/n zero cross\n"
61  << "is returned.\n";
62 }
std::ofstream outFile
Definition: GammaRayTel.cc:68
G4double G4IonsShenCrossSection::GetElementCrossSection ( const G4DynamicParticle aParticle,
G4int  Z,
const G4Material  
)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 71 of file G4IonsShenCrossSection.cc.

References G4lrint(), GetIsoCrossSection(), and G4NistManager::Instance().

74 {
75  G4int A = G4lrint(G4NistManager::Instance()->GetAtomicMassAmu(Z));
76  return GetIsoCrossSection(aParticle, Z, A);
77 }
int G4int
Definition: G4Types.hh:78
static G4NistManager * Instance()
int G4lrint(double ad)
Definition: templates.hh:163
virtual G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
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 79 of file G4IonsShenCrossSection.cc.

References G4Pow::A13(), test::b, test::c, python.hepunit::eplus, G4lrint(), G4ParticleDefinition::GetBaryonNumber(), G4DynamicParticle::GetDefinition(), G4Pow::GetInstance(), G4DynamicParticle::GetKineticEnergy(), G4DynamicParticle::GetMass(), G4DynamicParticle::GetMomentum(), G4NucleiProperties::GetNuclearMass(), G4ParticleDefinition::GetPDGCharge(), CLHEP::Hep3Vector::mag(), python.hepunit::MeV, python.hepunit::millibarn, python.hepunit::pi, and G4Pow::Z13().

Referenced by GetElementCrossSection().

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 64 of file G4IonsShenCrossSection.cc.

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

66 {
67  return (1 <= aDP->GetDefinition()->GetBaryonNumber());
68 }
G4ParticleDefinition * GetDefinition() const

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