source: trunk/source/processes/hadronic/cross_sections/src/G4IonsShenCrossSection.cc@ 1007

Last change on this file since 1007 was 819, checked in by garnier, 17 years ago

import all except CVS

File size: 5.7 KB
Line 
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//
31
32#include "G4IonsShenCrossSection.hh"
33#include "G4ParticleTable.hh"
34#include "G4IonTable.hh"
35
36
37G4double G4IonsShenCrossSection::
38GetIsoZACrossSection(const G4DynamicParticle* aParticle, G4double ZZ,
39 G4double AA, G4double /*temperature*/)
40{
41 G4double xsection = 0.0;
42
43 G4int Ap = aParticle->GetDefinition()->GetBaryonNumber();
44 G4int Zp = int ( aParticle->GetDefinition()->GetPDGCharge() / eplus + 0.5 );
45 G4double ke_per_N = aParticle->GetKineticEnergy() / Ap;
46 if ( ke_per_N > 10*GeV ) ke_per_N = 10*GeV;
47
48// Apply energy check, if less than lower limit then 0 value is returned
49 if ( ke_per_N < lowerLimit )
50 return xsection;
51
52 G4int At = G4int(AA);
53 G4int Zt = G4int(ZZ);
54
55 G4double one_third = 1.0 / 3.0;
56
57 G4double cubicrAt = std::pow ( G4double(At) , G4double(one_third) );
58 G4double cubicrAp = std::pow ( G4double(Ap) , G4double(one_third) );
59
60 G4double Rt = 1.12 * cubicrAt - 0.94 * ( 1.0 / cubicrAt );
61 G4double Rp = 1.12 * cubicrAp - 0.94 * ( 1.0 / cubicrAp );
62
63 G4double r = Rt + Rp + 3.2; // in fm
64 G4double b = 1.0; // in MeV/fm
65
66 G4double B = 1.44 * Zt * Zp / r - b * Rt * Rp / ( Rt + Rp );
67
68 //G4double ke_per_N = aParticle->GetKineticEnergy() / Ap;
69
70 G4double c = calCeValue ( ke_per_N / MeV );
71
72 G4double R1 = r0 * ( cubicrAt + cubicrAp + 1.85 * cubicrAt * cubicrAp / ( cubicrAt + cubicrAp ) - c);
73
74 G4double R2 = 1.0 * ( At - 2 * Zt ) * Zp / ( Ap * At );
75
76 G4double targ_mass = G4ParticleTable::GetParticleTable()->GetIonTable()->GetIonMass( Zt , At );
77 G4double proj_mass = aParticle->GetMass();
78 G4double proj_momentum = aParticle->GetMomentum().mag();
79
80 G4double Ecm = calEcmValue ( proj_mass , targ_mass , proj_momentum );
81
82 G4double R3 = 0.176 / std::pow ( G4double(Ecm) , G4double(one_third) ) * cubicrAt * cubicrAp / ( cubicrAt + cubicrAp );
83
84 G4double R = R1 + R2 + R3;
85
86 xsection = 10 * pi * R * R * ( 1 - B / Ecm );
87 xsection = xsection * millibarn; // mulitply xsection by millibarn
88
89 return xsection;
90}
91
92
93G4double G4IonsShenCrossSection::
94GetCrossSection(const G4DynamicParticle* aParticle, const G4Element* anElement,
95 G4double temperature)
96{
97 G4int nIso = anElement->GetNumberOfIsotopes();
98 G4double xsection = 0;
99
100 if (nIso) {
101 G4double sig;
102 G4IsotopeVector* isoVector = anElement->GetIsotopeVector();
103 G4double* abundVector = anElement->GetRelativeAbundanceVector();
104 G4double ZZ;
105 G4double AA;
106
107 for (G4int i = 0; i < nIso; i++) {
108 ZZ = G4double( (*isoVector)[i]->GetZ() );
109 AA = G4double( (*isoVector)[i]->GetN() );
110 sig = GetIsoZACrossSection(aParticle, ZZ, AA, temperature);
111 xsection += sig*abundVector[i];
112 }
113
114 } else {
115 xsection =
116 GetIsoZACrossSection(aParticle, anElement->GetZ(), anElement->GetN(),
117 temperature);
118 }
119
120 return xsection;
121}
122
123
124G4double G4IonsShenCrossSection::calEcmValue( const G4double mp , const G4double mt , const G4double Plab )
125{
126 G4double Elab = std::sqrt ( mp * mp + Plab * Plab );
127 G4double Ecm = std::sqrt ( mp * mp + mt * mt + 2 * Elab * mt );
128 G4double Pcm = Plab * mt / Ecm;
129 G4double KEcm = std::sqrt ( Pcm * Pcm + mp * mp ) - mp;
130 return KEcm;
131}
132
133
134G4double G4IonsShenCrossSection::calCeValue( const G4double ke )
135{
136 // Calculate c value
137 // This value is indepenent from projectile and target particle
138 // ke is projectile kinetic energy per nucleon in the Lab system with MeV unit
139 // fitting function is made by T. Koi
140 // There are no data below 30 MeV/n in Kox et al.,
141
142 G4double Ce;
143 G4double log10_ke = std::log10 ( ke );
144 if ( log10_ke > 1.5 )
145 {
146 Ce = - 10.0 / std::pow ( G4double(log10_ke) , G4double(5) ) + 2.0;
147 }
148 else
149 {
150 Ce = ( - 10.0 / std::pow ( G4double(1.5) , G4double(5) ) + 2.0 ) / std::pow ( G4double(1.5) , G4double(3) ) * std::pow ( G4double(log10_ke) , G4double(3) );
151 }
152 return Ce;
153}
Note: See TracBrowser for help on using the repository browser.