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