source: trunk/source/processes/hadronic/cross_sections/src/G4ProtonInelasticCrossSection.cc @ 1355

Last change on this file since 1355 was 1340, checked in by garnier, 14 years ago

update ti head

File size: 4.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// By JPW, working, but to be cleaned up. @@@
27// G.Folger, 29-sept-2006: extend to 1TeV, using a constant above 20GeV
28// 22 Dec 2006 - DHW added isotope dependence
29// G.Folger, 25-Nov-2009: extend to 100TeV, using a constant above 20GeV
30//
31
32#include "G4ProtonInelasticCrossSection.hh"
33#include "G4HadTmpUtil.hh"
34#include "globals.hh"
35
36
37G4double G4ProtonInelasticCrossSection::
38GetCrossSection(const G4DynamicParticle* aPart, 
39                const G4Element* anEle, G4double /*aTemperature*/)
40{
41  G4int nIso = anEle->GetNumberOfIsotopes();
42  G4double KE = aPart->GetKineticEnergy(); 
43  G4double cross_section = 0;
44 
45  if (nIso) {
46    G4double psig;
47    G4IsotopeVector* isoVector = anEle->GetIsotopeVector();
48    G4double* abundVector = anEle->GetRelativeAbundanceVector();
49    G4int ZZ;
50    G4int AA;
51 
52    for (G4int i = 0; i < nIso; i++) {
53      ZZ = (*isoVector)[i]->GetZ();
54      AA = (*isoVector)[i]->GetN();
55      psig = GetCrossSection(KE, AA, ZZ);
56      cross_section += psig*abundVector[i];
57    }
58 
59  } else {
60    G4int ZZ = G4lrint(anEle->GetZ());
61    G4int AA = G4lrint(anEle->GetN());
62    cross_section = GetCrossSection(KE, AA, ZZ);
63  }
64
65  return cross_section;
66}
67
68
69G4double G4ProtonInelasticCrossSection::
70GetCrossSection(G4double kineticEnergy, G4int atomicNumber, G4int nOfProtons)
71{   
72  if (kineticEnergy > 19.9*GeV ) 
73  { // constant cross section above ~20GeV.
74    return  GetCrossSection(19.8*GeV,atomicNumber,nOfProtons);
75  } 
76  G4int nOfNeutrons = atomicNumber-nOfProtons;
77  kineticEnergy /=GeV;
78  G4double a = atomicNumber;
79  const G4double nuleonRadius=1.36E-15;
80  const G4double pi=3.14159265;
81  G4double fac=pi*nuleonRadius*nuleonRadius;
82  G4double b0=2.247-0.915*(1-std::pow(a,-0.3333));
83  G4double fac1=b0*(1-std::pow(a,-0.3333));
84  G4double fac2=1.;
85  if(nOfNeutrons > 1) fac2=std::log((G4double(nOfNeutrons)));
86  G4double crossSection = 1E31*fac*fac2*(1+std::pow(a,0.3333)-fac1);
87
88  // high energy correction
89
90  crossSection = (1-0.15*std::exp(-kineticEnergy))*crossSection/(1.00-0.0007*a);
91  // first try on low energies: rise
92
93  G4double ff1= 0.70-0.002*a;  // slope of the drop at medium energies.
94  G4double ff2= 1.00+1/a;  // start of the slope.
95  G4double ff3= 0.8+18/a-0.002*a; // stephight
96  fac = 1.0;
97  if (kineticEnergy > DBL_MIN) 
98     fac= 1.0 - (1.0/(1+std::exp(-8*ff1*(std::log10(kineticEnergy)+1.37*ff2))));
99  crossSection = crossSection*(1+ff3*fac);
100
101  // low energy return to zero
102
103  ff1=1.-1/a-0.001*a; // slope of the rise
104  ff2=1.17-2.7/a-0.0014*a; // start of the rise
105  fac = 0.0;
106  if (kineticEnergy > DBL_MIN) {
107    fac=-8.*ff1*(std::log10(kineticEnergy)+2.0*ff2);
108    fac=1/(1+std::exp(fac));
109  }
110  crossSection = crossSection*fac;
111  return crossSection*millibarn;
112}
Note: See TracBrowser for help on using the repository browser.