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

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

import all except CVS

File size: 6.1 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// Implementation of formulas in analogy to NASA technical paper 3621 by
27// Tripathi, et al.
28//
29// 26-Dec-2006 Isotope dependence added by D. Wright
30//
31
32#include "G4TripathiCrossSection.hh"
33#include "G4ParticleTable.hh"
34#include "G4IonTable.hh"
35#include "G4HadTmpUtil.hh"
36
37G4double G4TripathiCrossSection::
38GetIsoZACrossSection(const G4DynamicParticle* aPart, G4double ZZ, G4double AA,
39 G4double /*temperature*/)
40{
41 G4double result = 0;
42
43 const G4double targetAtomicNumber = AA;
44 const G4double nTargetProtons = ZZ;
45
46 const G4double kineticEnergy = aPart->GetKineticEnergy()/MeV;
47 const G4double nProjProtons = aPart->GetDefinition()->GetPDGCharge();
48 const G4double projectileAtomicNumber =
49 aPart->GetDefinition()->GetBaryonNumber();
50
51 const G4double nuleonRadius=1.1E-15;
52 const G4double myNuleonRadius=1.36E-15;
53
54 // needs target mass
55 G4double targetMass =
56 G4ParticleTable::GetParticleTable()->GetIonTable()
57 ->GetIonMass(G4lrint(nTargetProtons), G4lrint(targetAtomicNumber));
58 G4LorentzVector pTarget(0,0,0,targetMass);
59 G4LorentzVector pProjectile(aPart->Get4Momentum());
60 pTarget = pTarget+pProjectile;
61 G4double E_cm = (pTarget.mag()-targetMass-pProjectile.m())/MeV;
62
63 // done
64 G4double r_rms_p = 0.6 * myNuleonRadius *
65 std::pow(projectileAtomicNumber, 1./3.);
66 G4double r_rms_t = 0.6 * myNuleonRadius *
67 std::pow(targetAtomicNumber, 1./3.);
68
69 // done
70 G4double r_p = 1.29*r_rms_p/nuleonRadius ;
71 G4double r_t = 1.29*r_rms_t/nuleonRadius;
72
73 // done
74 G4double Radius = r_p + r_t +
75 1.2*(std::pow(targetAtomicNumber, 1./3.) +
76 std::pow(projectileAtomicNumber, 1./3.))/std::pow(E_cm, 1./3.);
77
78 //done
79 G4double B = 1.44*nProjProtons*nTargetProtons/Radius;
80
81 // done
82 G4double Energy = kineticEnergy/projectileAtomicNumber;
83
84 // done
85 //
86 // Note that this correction to G4TripathiCrossSection is just to accurately
87 // reflect Tripathi's algorithm. However, if you're using alpha
88 // particles/protons consider using the more accurate
89 // G4TripathiLightCrossSection, which Tripathi developed specifically for
90 // light systems.
91 //
92
93 G4double D;
94 if (nProjProtons==1 && projectileAtomicNumber==1)
95 {
96 D = 2.05;
97 }
98 else if (nProjProtons==2 && projectileAtomicNumber==4)
99 {
100 D = 2.77-(8.0E-3*targetAtomicNumber)+
101 (1.8E-5*targetAtomicNumber*targetAtomicNumber)
102 - 0.8/(1+std::exp((250.-Energy)/75.));
103 }
104 else
105 {
106 //
107 // This is the original value used in the G4TripathiCrossSection
108 // implementation, and was used for all projectile/target conditions.
109 // I'm not touching this, although judging from Tripathi's paper, this is
110 // valid for cases where the nucleon density changes little with A.
111 //
112 D = 1.75;
113 }
114 // done
115 G4double C_E = D * (1-std::exp(-Energy/40.)) -
116 0.292*std::exp(-Energy/792.)*std::cos(0.229*std::pow(Energy, 0.453));
117
118 // done
119 G4double S = std::pow(projectileAtomicNumber, 1./3.)*
120 std::pow(targetAtomicNumber, 1./3.)/
121 (std::pow(projectileAtomicNumber, 1./3.) +
122 std::pow(targetAtomicNumber, 1./3.));
123
124 // done
125 G4double deltaE = 1.85*S + 0.16*S/std::pow(E_cm,1./3.) - C_E +
126 0.91*(targetAtomicNumber-2.*nTargetProtons)*nProjProtons/
127 (targetAtomicNumber*projectileAtomicNumber);
128
129 // done
130 result = pi * nuleonRadius*nuleonRadius *
131 std::pow(( std::pow(targetAtomicNumber, 1./3.) +
132 std::pow(projectileAtomicNumber, 1./3.) + deltaE),2.) *
133 (1-B/E_cm);
134
135 if(result < 0) result = 0;
136 return result*m2;
137
138}
139
140
141G4double G4TripathiCrossSection::
142GetCrossSection(const G4DynamicParticle* aPart, const G4Element* anEle,
143 G4double temperature)
144{
145 G4int nIso = anEle->GetNumberOfIsotopes();
146 G4double xsection = 0;
147
148 if (nIso) {
149 G4double sig;
150 G4IsotopeVector* isoVector = anEle->GetIsotopeVector();
151 G4double* abundVector = anEle->GetRelativeAbundanceVector();
152 G4double ZZ;
153 G4double AA;
154
155 for (G4int i = 0; i < nIso; i++) {
156 ZZ = G4double( (*isoVector)[i]->GetZ() );
157 AA = G4double( (*isoVector)[i]->GetN() );
158 sig = GetIsoZACrossSection(aPart, ZZ, AA, temperature);
159 xsection += sig*abundVector[i];
160 }
161
162 } else {
163 xsection =
164 GetIsoZACrossSection(aPart, anEle->GetZ(), anEle->GetN(),
165 temperature);
166 }
167
168 return xsection;
169}
Note: See TracBrowser for help on using the repository browser.