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

Last change on this file since 1315 was 1228, checked in by garnier, 16 years ago

update geant4.9.3 tag

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