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

Last change on this file since 1197 was 1196, checked in by garnier, 15 years ago

update CVS release candidate geant4.9.3.01

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