Geant4
10.00.p03
Main Page
Related Pages
Modules
Namespaces
Classes
Files
File List
File Members
G4ProtonInelasticCrossSection.cc
Go to the documentation of this file.
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
// V.Ivanchenko, 18-Aug-2011: migration to new design and cleanup;
31
// make it applicable for Z>1
32
//
33
34
#include "
G4ProtonInelasticCrossSection.hh
"
35
#include "
G4SystemOfUnits.hh
"
36
#include "
G4DynamicParticle.hh
"
37
#include "
G4Proton.hh
"
38
#include "
G4HadTmpUtil.hh
"
39
#include "
G4NistManager.hh
"
40
41
G4ProtonInelasticCrossSection::G4ProtonInelasticCrossSection
()
42
:
G4VCrossSectionDataSet
(
"Axen-Wellisch"
),thEnergy(19.8*
CLHEP
::
GeV
)
43
{
44
nist
=
G4NistManager::Instance
();
45
theProton
=
G4Proton::Proton
();
46
}
47
48
G4ProtonInelasticCrossSection::~G4ProtonInelasticCrossSection
()
49
{}
50
51
G4bool
52
G4ProtonInelasticCrossSection::IsElementApplicable
(
53
const
G4DynamicParticle
* aPart,
54
G4int
Z,
const
G4Material
*)
55
{
56
return
((1 < Z) && (aPart->
GetDefinition
() ==
theProton
));
57
}
58
59
G4double
G4ProtonInelasticCrossSection::GetElementCrossSection
(
60
const
G4DynamicParticle
* aPart,
61
G4int
Z,
const
G4Material
*)
62
{
63
return
GetProtonCrossSection
(aPart->
GetKineticEnergy
(), Z);
64
}
65
66
G4double
G4ProtonInelasticCrossSection::GetProtonCrossSection
(
67
G4double
kineticEnergy,
G4int
Z)
68
{
69
if
(kineticEnergy <= 0.0) {
return
0.0; }
70
71
// constant cross section above ~20GeV
72
if
(kineticEnergy >
thEnergy
) { kineticEnergy =
thEnergy
; }
73
74
G4double
a
=
nist
->
GetAtomicMassAmu
(Z);
75
G4double
a13 = std::pow(a,-0.3333333333);
76
G4int
nOfNeutrons =
G4lrint
(a) - Z;
77
kineticEnergy /=
GeV
;
78
G4double
alog10E = std::log10(kineticEnergy);
79
80
static
const
G4double
nuleonRadius=1.36e-15;
81
static
const
G4double
fac
=
CLHEP::pi
*nuleonRadius*nuleonRadius;
82
83
G4double
b0
= 2.247-0.915*(1 - a13);
84
G4double
fac1 = b0*(1 - a13);
85
G4double
fac2 = 1.;
86
if
(nOfNeutrons > 1) { fac2=std::log((
G4double
(nOfNeutrons))); }
87
G4double
crossSection = 1.0E31*fac*fac2*(1. + 1./a13 - fac1);
88
89
// high energy correction
90
crossSection *= (1 - 0.15*std::exp(-kineticEnergy))/(1.0 - 0.0007*a);
91
92
// first try on low energies: rise
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
97
G4double
ff4= 1.0 - (1.0/(1+std::exp(-8*ff1*(alog10E + 1.37*ff2))));
98
99
crossSection *= (1 + ff3*ff4);
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
106
ff4=-8.*ff1*(alog10E + 2.0*ff2);
107
108
crossSection *=
millibarn
/(1. + std::exp(ff4));
109
return
crossSection;
110
}
fac
static const G4double fac
Definition:
G4ASTARStopping.cc:56
G4ProtonInelasticCrossSection::G4ProtonInelasticCrossSection
G4ProtonInelasticCrossSection()
Definition:
G4ProtonInelasticCrossSection.cc:41
G4DynamicParticle::GetKineticEnergy
G4double GetKineticEnergy() const
G4NistManager.hh
G4INCL::Math::pi
const G4double pi
Definition:
G4INCLGlobals.hh:66
G4Material
Definition:
G4Material.hh:118
G4DynamicParticle
Definition:
G4DynamicParticle.hh:73
G4DynamicParticle.hh
G4ProtonInelasticCrossSection.hh
G4ProtonInelasticCrossSection::theProton
G4Proton * theProton
Definition:
G4ProtonInelasticCrossSection.hh:73
G4ProtonInelasticCrossSection::GetElementCrossSection
virtual G4double GetElementCrossSection(const G4DynamicParticle *, G4int Z, const G4Material *)
Definition:
G4ProtonInelasticCrossSection.cc:59
G4ProtonInelasticCrossSection::~G4ProtonInelasticCrossSection
~G4ProtonInelasticCrossSection()
Definition:
G4ProtonInelasticCrossSection.cc:48
G4DynamicParticle::GetDefinition
G4ParticleDefinition * GetDefinition() const
a
G4double a
Definition:
TRTMaterials.hh:39
G4VCrossSectionDataSet
Definition:
G4VCrossSectionDataSet.hh:70
G4int
int G4int
Definition:
G4Types.hh:78
G4NistManager::Instance
static G4NistManager * Instance()
Definition:
G4NistManager.cc:68
G4bool
bool G4bool
Definition:
G4Types.hh:79
G4ProtonInelasticCrossSection::IsElementApplicable
virtual G4bool IsElementApplicable(const G4DynamicParticle *aPart, G4int Z, const G4Material *)
Definition:
G4ProtonInelasticCrossSection.cc:52
G4Proton.hh
G4Proton::Proton
static G4Proton * Proton()
Definition:
G4Proton.cc:93
GeV
static const double GeV
Definition:
G4SIunits.hh:196
G4lrint
int G4lrint(double ad)
Definition:
templates.hh:163
CLHEP
Definition:
AxisAngle.cc:16
G4ProtonInelasticCrossSection::thEnergy
G4double thEnergy
Definition:
G4ProtonInelasticCrossSection.hh:74
G4NistManager::GetAtomicMassAmu
G4double GetAtomicMassAmu(const G4String &symb) const
Definition:
G4NistManager.hh:344
millibarn
static const double millibarn
Definition:
G4SIunits.hh:96
G4ProtonInelasticCrossSection::nist
G4NistManager * nist
Definition:
G4ProtonInelasticCrossSection.hh:72
G4double
double G4double
Definition:
G4Types.hh:76
G4SystemOfUnits.hh
b0
static const G4double b0
Definition:
G4BetheHeitlerModel.cc:76
G4ProtonInelasticCrossSection::GetProtonCrossSection
G4double GetProtonCrossSection(G4double kineticEnergy, G4int Z)
Definition:
G4ProtonInelasticCrossSection.cc:66
G4HadTmpUtil.hh
source
source
processes
hadronic
cross_sections
src
G4ProtonInelasticCrossSection.cc
Generated on Sat Dec 13 2014 18:02:27 for Geant4 by
1.8.8