Geant4  10.02.p01
G4IonsShenCrossSection.cc
Go to the documentation of this file.
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the result of the scientific and *
19 // * technical work of the GEANT4 collaboration. *
20 // * By using, copying, modifying or distributing the software (or *
21 // * any work based on the software) you agree to acknowledge its *
22 // * use in resulting scientific publications, and indicate your *
23 // * acceptance of all terms of the Geant4 Software license. *
24 // ********************************************************************
25 //
26 // 18-Sep-2003 First version is written by T. Koi
27 // 12-Nov-2003 Add energy check at lower side T. Koi
28 // 15-Nov-2006 Above 10GeV/n Cross Section become constant T. Koi (SLAC/SCCS)
29 // 23-Dec-2006 Isotope dependence adde by D. Wright
30 // 14-Mar-2011 Moved constructor, destructor and virtual methods to source by V.Ivanchenko
31 // 19-Aug-2011 V.Ivanchenko move to new design and make x-section per element
32 //
33 
35 #include "G4PhysicalConstants.hh"
36 #include "G4SystemOfUnits.hh"
37 #include "G4DynamicParticle.hh"
38 #include "G4NucleiProperties.hh"
39 #include "G4HadTmpUtil.hh"
40 #include "G4NistManager.hh"
41 
43  : G4VCrossSectionDataSet("IonsShen"),
44  upperLimit( 10*GeV ),
45 // lowerLimit( 10*MeV ),
46  r0 ( 1.1 )
47 {}
48 
50 {}
51 
52 void
54 {
55  outFile << "G4IonsShenCrossSection calculates the total reaction cross\n"
56  << "section for nucleus-nucleus scattering using the Shen\n"
57  << "parameterization. It is valid for projectiles and targets of\n"
58  << "all Z, and projectile energies up to 1 TeV/n. Above 10 GeV/n"
59  << "the cross section is constant. Below 10 MeV/n zero cross\n"
60  << "is returned.\n";
61 }
62 
64  G4int, const G4Material*)
65 {
66  return (1 <= aDP->GetDefinition()->GetBaryonNumber());
67 }
68 
69 G4double
71  G4int Z,
72  const G4Material*)
73 {
74  G4int A = G4lrint(G4NistManager::Instance()->GetAtomicMassAmu(Z));
75  return GetIsoCrossSection(aParticle, Z, A);
76 }
77 
79  G4int Zt, G4int At,
80  const G4Isotope*,
81  const G4Element*,
82  const G4Material*)
83 
84 {
85  G4double xsection = 0.0;
86 
87  G4int Ap = aParticle->GetDefinition()->GetBaryonNumber();
88  G4int Zp = G4lrint(aParticle->GetDefinition()->GetPDGCharge()/eplus);
89  G4double ke_per_N = aParticle->GetKineticEnergy() / Ap;
90  if ( ke_per_N > upperLimit ) { ke_per_N = upperLimit; }
91 
92  // Apply energy check, if less than lower limit then 0 value is returned
93  //if ( ke_per_N < lowerLimit ) { return xsection; }
94 
95  G4Pow* g4pow = G4Pow::GetInstance();
96 
97  G4double cubicrAt = g4pow->Z13(At);
98  G4double cubicrAp = g4pow->Z13(Ap);
99 
100  G4double Rt = 1.12 * cubicrAt - 0.94 * ( 1.0 / cubicrAt );
101  G4double Rp = 1.12 * cubicrAp - 0.94 * ( 1.0 / cubicrAp );
102 
103  G4double r = Rt + Rp + 3.2; // in fm
104  G4double b = 1.0; // in MeV/fm
105  G4double targ_mass = G4NucleiProperties::GetNuclearMass(At, Zt);
106 
107  G4double proj_mass = aParticle->GetMass();
108  G4double proj_momentum = aParticle->GetMomentum().mag();
109 
110  G4double Ecm = calEcmValue (proj_mass, targ_mass, proj_momentum);
111 
112  G4double B = 1.44 * Zt * Zp / r - b * Rt * Rp / ( Rt + Rp );
113  if(Ecm <= B) { return xsection; }
114 
115  G4double c = calCeValue ( ke_per_N / MeV );
116 
117  G4double R1 = r0 * (cubicrAt + cubicrAp + 1.85*cubicrAt*cubicrAp/(cubicrAt + cubicrAp) - c);
118 
119  G4double R2 = 1.0 * ( At - 2 * Zt ) * Zp / ( Ap * At );
120 
121 
122  G4double R3 = (0.176 / g4pow->A13(Ecm)) * cubicrAt * cubicrAp /(cubicrAt + cubicrAp);
123 
124  G4double R = R1 + R2 + R3;
125 
126  xsection = 10 * pi * R * R * ( 1 - B / Ecm );
127  xsection = xsection * millibarn; // mulitply xsection by millibarn
128 
129  return xsection;
130 }
131 
132 G4double
134  const G4double Plab)
135 {
136  G4double Elab = std::sqrt ( mp * mp + Plab * Plab );
137  G4double Ecm = std::sqrt ( mp * mp + mt * mt + 2 * Elab * mt );
138  G4double Pcm = Plab * mt / Ecm;
139  G4double KEcm = std::sqrt ( Pcm * Pcm + mp * mp ) - mp;
140  return KEcm;
141 }
142 
143 
145 {
146  // Calculate c value
147  // This value is indepenent from projectile and target particle
148  // ke is projectile kinetic energy per nucleon in the Lab system
149  // with MeV unit
150  // fitting function is made by T. Koi
151  // There are no data below 30 MeV/n in Kox et al.,
152 
153  G4double Ce;
154  G4double log10_ke = std::log10 ( ke );
155  if (log10_ke > 1.5)
156  {
157  Ce = -10.0/std::pow(G4double(log10_ke), G4double(5)) + 2.0;
158  }
159  else
160  {
161  Ce = (-10.0/std::pow(G4double(1.5), G4double(5) ) + 2.0) /
162  std::pow(G4double(1.5) , G4double(3)) * std::pow(G4double(log10_ke), G4double(3));
163  }
164  return Ce;
165 }
166 
static G4Pow * GetInstance()
Definition: G4Pow.cc:55
static const double MeV
Definition: G4SIunits.hh:211
static G4double GetNuclearMass(const G4double A, const G4double Z)
G4double GetKineticEnergy() const
Definition: G4Pow.hh:56
virtual void CrossSectionDescription(std::ostream &) const
G4ParticleDefinition * GetDefinition() const
double B(double temperature)
int G4int
Definition: G4Types.hh:78
static G4NistManager * Instance()
G4double calCeValue(const G4double)
G4double Z13(G4int Z) const
Definition: G4Pow.hh:127
double A(double temperature)
G4double GetMass() const
bool G4bool
Definition: G4Types.hh:79
virtual G4double GetElementCrossSection(const G4DynamicParticle *, G4int Z, const G4Material *)
static const double GeV
Definition: G4SIunits.hh:214
static const double pi
Definition: G4SIunits.hh:74
G4double A13(G4double A) const
Definition: G4Pow.hh:132
int G4lrint(double ad)
Definition: templates.hh:163
G4double calEcmValue(const G4double, const G4double, const G4double)
virtual G4bool IsElementApplicable(const G4DynamicParticle *aDP, G4int Z, const G4Material *)
static const double millibarn
Definition: G4SIunits.hh:105
virtual G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
double G4double
Definition: G4Types.hh:76
static const double eplus
Definition: G4SIunits.hh:196
G4double GetPDGCharge() const
const G4double r0
G4ThreeVector GetMomentum() const