Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
G4IonsKoxCrossSection.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// 10-Nov-2003 Bug fix at Cal. ke_per_n and D T. Koi
28// 12-Nov-2003 Add energy check at lower side T. Koi
29// 26-Dec-2006 Add isotope dependence 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
35#include "G4SystemOfUnits.hh"
36#include "G4DynamicParticle.hh"
37#include "G4NucleiProperties.hh"
38#include "G4HadTmpUtil.hh"
39#include "G4NistManager.hh"
40
42 : G4VCrossSectionDataSet("IonsKox"), lowerLimit ( 10*MeV ),
43 r0 ( 1.1*fermi ), rc ( 1.3*fermi )
44{}
45
47{}
48
49void
51{
52 outFile << "G4IonsKoxCrossSection calculates the total reaction cross\n"
53 << "section for nucleus-nucleus scattering using the Kox\n"
54 << "parameterization. It is valid for projectiles and targets\n"
55 << "of all Z, at projectile energies up to 10 GeV/n. If the\n"
56 << "projectile energy is less than 10 MeV/n, a zero cross section\n"
57 << "is returned.\n";
58}
59
61 G4int, const G4Material*)
62{
63 return (1 <= aDP->GetDefinition()->GetBaryonNumber());
64}
65
68 const G4DynamicParticle* aParticle, G4int ZZ, const G4Material*)
69{
70 G4double xsection = 0.0;
71
72 G4int Ap = aParticle->GetDefinition()->GetBaryonNumber();
73 G4int Zp = G4int(aParticle->GetDefinition()->GetPDGCharge() / eplus + 0.5);
74 G4double ke_per_N = aParticle->GetKineticEnergy() / Ap;
75
76 // Apply energy check, if less than lower limit then 0 value is returned
77 // if ( ke_per_N < lowerLimit ) return xsection;
78
79 G4int At = G4lrint(G4NistManager::Instance()->GetAtomicMassAmu(ZZ));
80 G4int Zt = ZZ;
81
82 G4double one_third = 1.0 / 3.0;
83
84 G4double cubicrAt = std::pow ( G4double(At) , G4double(one_third) );
85 G4double cubicrAp = std::pow ( G4double(Ap) , G4double(one_third) );
86
87 // rc divide fermi
88 G4double Bc = Zt * Zp / ( (rc/fermi) * (cubicrAp+cubicrAt) );
89
91 G4double proj_mass = aParticle->GetMass();
92 G4double proj_momentum = aParticle->GetMomentum().mag();
93
94 G4double Ecm = calEcm ( proj_mass , targ_mass , proj_momentum );
95 if( Ecm <= Bc) return xsection;
96
97 G4double Rvol = r0 * ( cubicrAp + cubicrAt );
98
99// G4double ke_per_N = aParticle->GetKineticEnergy() / Ap;
100 G4double c = calCeValue ( ke_per_N / MeV );
101
102 G4double a = 1.85;
103 G4double Rsurf = r0 * (a*cubicrAp * cubicrAt/(cubicrAp + cubicrAt) - c);
104 G4double D = 5.0 * ( At - 2 * Zt ) * Zp / ( Ap * At );
105 Rsurf = Rsurf + D * fermi; // multiply D by fermi
106
107 G4double Rint = Rvol + Rsurf;
108 xsection = pi * Rint * Rint * ( 1 - Bc / ( Ecm / MeV ) );
109
110 return xsection;
111}
112
114G4IonsKoxCrossSection::calEcm(G4double mp, G4double mt, G4double Plab)
115{
116 G4double Elab = std::sqrt ( mp * mp + Plab * Plab );
117 G4double Ecm = std::sqrt ( mp * mp + mt * mt + 2 * Elab * mt );
118 G4double Pcm = Plab * mt / Ecm;
119 G4double KEcm = std::sqrt ( Pcm * Pcm + mp * mp ) - mp;
120 return KEcm;
121}
122
123
124G4double G4IonsKoxCrossSection::calCeValue(const G4double ke)
125{
126 // Calculate c value
127 // This value is indepenent from projectile and target particle
128 // ke is projectile kinetic energy per nucleon in the Lab system with MeV unit
129 // fitting function is made by T. Koi
130 // There are no data below 30 MeV/n in Kox et al.,
131
132 G4double Ce;
133 G4double log10_ke = std::log10 ( ke );
134 if (log10_ke > 1.5)
135 {
136 Ce = - 10.0 / std::pow ( G4double(log10_ke) , G4double(5) ) + 2.0;
137 }
138 else
139 {
140 Ce = (-10.0/std::pow(G4double(1.5), G4double(5) ) + 2.0) /
141 std::pow(G4double(1.5), G4double(3)) * std::pow(G4double(log10_ke), G4double(3) );
142
143 }
144 return Ce;
145}
146
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
double mag() const
G4double GetMass() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4ThreeVector GetMomentum() const
virtual G4bool IsElementApplicable(const G4DynamicParticle *aDP, G4int Z, const G4Material *)
virtual G4double GetElementCrossSection(const G4DynamicParticle *, G4int Z, const G4Material *)
virtual void CrossSectionDescription(std::ostream &) const
static G4NistManager * Instance()
static G4double GetNuclearMass(const G4double A, const G4double Z)
G4double GetPDGCharge() const
int G4lrint(double ad)
Definition: templates.hh:163