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

Last change on this file since 1347 was 1347, checked in by garnier, 13 years ago

geant4 tag 9.4

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