Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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 #include "G4Pow.hh"
42 
44  : G4VCrossSectionDataSet("IonsShen"),
45  upperLimit( 10*GeV ),
46 // lowerLimit( 10*MeV ),
47  r0 ( 1.1 )
48 {}
49 
51 {}
52 
53 void
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 }
63 
65  G4int, const G4Material*)
66 {
67  return (1 <= aDP->GetDefinition()->GetBaryonNumber());
68 }
69 
70 G4double
72  G4int Z,
73  const G4Material*)
74 {
75  G4int A = G4lrint(G4NistManager::Instance()->GetAtomicMassAmu(Z));
76  return GetIsoCrossSection(aParticle, Z, A);
77 }
78 
80  G4int Zt, G4int At,
81  const G4Isotope*,
82  const G4Element*,
83  const G4Material*)
84 
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 }
132 
133 G4double
134 G4IonsShenCrossSection::calEcmValue(const G4double mp, const G4double mt,
135  const G4double Plab)
136 {
137  G4double Elab = std::sqrt ( mp * mp + Plab * Plab );
138  G4double Ecm = std::sqrt ( mp * mp + mt * mt + 2 * Elab * mt );
139  G4double Pcm = Plab * mt / Ecm;
140  G4double KEcm = std::sqrt ( Pcm * Pcm + mp * mp ) - mp;
141  return KEcm;
142 }
143 
144 
145 G4double G4IonsShenCrossSection::calCeValue(const G4double ke)
146 {
147  // Calculate c value
148  // This value is indepenent from projectile and target particle
149  // ke is projectile kinetic energy per nucleon in the Lab system
150  // with MeV unit
151  // fitting function is made by T. Koi
152  // There are no data below 30 MeV/n in Kox et al.,
153 
154  G4double Ce;
155  G4double log10_ke = std::log10 ( ke );
156  if (log10_ke > 1.5)
157  {
158  Ce = -10.0/std::pow(G4double(log10_ke), G4double(5)) + 2.0;
159  }
160  else
161  {
162  Ce = (-10.0/std::pow(G4double(1.5), G4double(5) ) + 2.0) /
163  std::pow(G4double(1.5) , G4double(3)) * std::pow(G4double(log10_ke), G4double(3));
164  }
165  return Ce;
166 }
167 
static G4Pow * GetInstance()
Definition: G4Pow.cc:53
static G4double GetNuclearMass(const G4double A, const G4double Z)
G4double GetKineticEnergy() const
Definition: G4Pow.hh:56
std::ofstream outFile
Definition: GammaRayTel.cc:68
virtual void CrossSectionDescription(std::ostream &) const
G4ParticleDefinition * GetDefinition() const
int G4int
Definition: G4Types.hh:78
int millibarn
Definition: hepunit.py:40
static G4NistManager * Instance()
G4double Z13(G4int Z) const
Definition: G4Pow.hh:129
G4double GetMass() const
bool G4bool
Definition: G4Types.hh:79
virtual G4double GetElementCrossSection(const G4DynamicParticle *, G4int Z, const G4Material *)
G4double A13(G4double A) const
Definition: G4Pow.hh:134
int G4lrint(double ad)
Definition: templates.hh:163
virtual G4bool IsElementApplicable(const G4DynamicParticle *aDP, 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)
double G4double
Definition: G4Types.hh:76
G4double GetPDGCharge() const
double mag() const
G4ThreeVector GetMomentum() const