source: trunk/source/processes/hadronic/cross_sections/src/G4IonsKoxCrossSection.cc @ 969

Last change on this file since 969 was 819, checked in by garnier, 16 years ago

import all except CVS

File size: 5.4 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// 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
31#include "G4IonsKoxCrossSection.hh"
32#include "G4ParticleTable.hh"
33#include "G4IonTable.hh"
34
35G4double G4IonsKoxCrossSection::
36GetIsoZACrossSection(const G4DynamicParticle* aParticle, G4double ZZ, 
37                     G4double AA, G4double /*temperature*/)
38{
39   G4double xsection = 0.0;
40
41   G4int Ap = aParticle->GetDefinition()->GetBaryonNumber();
42   G4int Zp = int ( aParticle->GetDefinition()->GetPDGCharge() / eplus + 0.5); 
43   G4double ke_per_N = aParticle->GetKineticEnergy() / Ap;
44
45// Apply energy check, if less than lower limit then 0 value is returned
46   if (  ke_per_N < lowerLimit )
47      return xsection;
48
49   G4int At = int (AA + 0.5);
50   G4int Zt = int (ZZ + 0.5 ); 
51
52   G4double one_third = 1.0 / 3.0;
53
54   G4double cubicrAt = std::pow ( G4double(At) , G4double(one_third) ); 
55   G4double cubicrAp = std::pow ( G4double(Ap) , G4double(one_third) ); 
56
57
58   G4double Bc = Zt * Zp / ( ( rc / fermi ) * (  cubicrAp + cubicrAt ) );   // rc divide fermi
59
60   G4double Rvol = r0 * (  cubicrAp + cubicrAt );
61
62//   G4double ke_per_N = aParticle->GetKineticEnergy() / Ap;
63   G4double c = calCeValue ( ke_per_N / MeV  ); 
64
65   G4double a = 1.85;
66   G4double Rsurf = r0 * ( a * cubicrAp * cubicrAt / ( cubicrAp + cubicrAt ) - c); 
67   G4double D = 5.0 * ( At - 2 * Zt ) * Zp / ( Ap * At );
68   Rsurf = Rsurf + D * fermi;  // multiply D by fermi
69
70   G4double Rint = Rvol + Rsurf;
71
72   G4double targ_mass = G4ParticleTable::GetParticleTable()->GetIonTable()->GetIonMass( Zt , At ); 
73   G4double proj_mass = aParticle->GetMass(); 
74   G4double proj_momentum = aParticle->GetMomentum().mag(); 
75
76   G4double Ecm = calEcm ( proj_mass , targ_mass , proj_momentum ); 
77
78   xsection = pi * Rint * Rint * ( 1 - Bc / ( Ecm / MeV ) );
79 
80   return xsection; 
81}
82
83G4double G4IonsKoxCrossSection::
84GetCrossSection(const G4DynamicParticle* aParticle, 
85                const G4Element* anElement, G4double temperature)
86{
87  G4int nIso = anElement->GetNumberOfIsotopes();
88  G4double xsection = 0;
89
90  if (nIso) {
91    G4double sig;
92    G4IsotopeVector* isoVector = anElement->GetIsotopeVector();
93    G4double* abundVector = anElement->GetRelativeAbundanceVector();
94    G4double ZZ;
95    G4double AA;
96     
97    for (G4int i = 0; i < nIso; i++) {
98      ZZ = G4double( (*isoVector)[i]->GetZ() );
99      AA = G4double( (*isoVector)[i]->GetN() );
100      sig = GetIsoZACrossSection(aParticle, ZZ, AA, temperature);
101      xsection += sig*abundVector[i];
102    }
103   
104  } else {
105    xsection =
106      GetIsoZACrossSection(aParticle, anElement->GetZ(), anElement->GetN(),
107                           temperature);
108  }
109   
110  return xsection;
111}
112
113
114G4double G4IonsKoxCrossSection::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 ) / std::pow ( G4double(1.5) , G4double(3) ) * std::pow ( G4double(log10_ke) , G4double(3) );
141
142   }
143   return Ce;
144}
Note: See TracBrowser for help on using the repository browser.