Geant4
11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4NucleonNuclearCrossSection.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
// author:
[email protected]
27
//
28
// Implements data from: Barashenkov V.S., Nucleon-Nucleus Cross Section,
29
// Preprint JINR P2-89-770, p. 12, Dubna 1989 (scanned version from KEK)
30
// Based on G. Folger version of G4PiNuclearCrossSection class
31
//
32
// Modified: 16.08.2018 V.Ivanchenko - major update
33
//
34
35
#include "
G4NucleonNuclearCrossSection.hh
"
36
#include "
G4DynamicParticle.hh
"
37
#include "
G4ParticleDefinition.hh
"
38
#include "
G4Neutron.hh
"
39
#include "
G4Proton.hh
"
40
#include "
G4ComponentBarNucleonNucleusXsc.hh
"
41
42
///////////////////////////////////////////////////////////////////////////////
43
44
G4NucleonNuclearCrossSection::G4NucleonNuclearCrossSection
()
45
:
G4VCrossSectionDataSet
(Default_Name()),
46
fTotalXsc(0.0), fInelasticXsc(0.0), fElasticXsc(0.0)
47
{
48
theNeutron =
G4Neutron::Neutron
();
49
theProton =
G4Proton::Proton
();
50
fBarash =
new
G4ComponentBarNucleonNucleusXsc
();
51
}
52
53
///////////////////////////////////////////////////////////////////////////////
54
//
55
56
G4NucleonNuclearCrossSection::~G4NucleonNuclearCrossSection
()
57
{}
58
59
////////////////////////////////////////////////////////////////////////////
60
61
G4bool
G4NucleonNuclearCrossSection::IsElementApplicable
(
62
const
G4DynamicParticle
*,
G4int
Z
,
const
G4Material
*)
63
{
64
return
(
Z
> 1);
65
}
66
67
////////////////////////////////////////////////////////////////////////////
68
69
G4double
G4NucleonNuclearCrossSection::GetElementCrossSection
(
70
const
G4DynamicParticle
* dp,
G4int
Z
,
const
G4Material
*)
71
{
72
ComputeCrossSections(dp->
GetDefinition
(), dp->
GetKineticEnergy
(),
Z
);
73
return
fInelasticXsc;
74
}
75
76
////////////////////////////////////////////////////////////////////////////
77
78
void
G4NucleonNuclearCrossSection::ComputeCrossSections(
79
const
G4ParticleDefinition
* pd,
80
G4double
kinEnergy,
G4int
Z
)
81
{
82
fBarash->
ComputeCrossSections
(pd, kinEnergy,
Z
);
83
fTotalXsc = fBarash->
GetTotalXsc
();
84
fInelasticXsc = fBarash->
GetInelasticXsc
();
85
fElasticXsc = fBarash->
GetElasticXsc
();
86
}
87
88
////////////////////////////////////////////////////////////////////////////
89
90
void
G4NucleonNuclearCrossSection::BuildPhysicsTable
(
const
G4ParticleDefinition
& part)
91
{
92
fBarash->
BuildPhysicsTable
(part);
93
}
94
95
////////////////////////////////////////////////////////////////////////////
96
97
void
98
G4NucleonNuclearCrossSection::CrossSectionDescription
(std::ostream& outFile)
const
99
{
100
outFile <<
"G4NucleonNuclearCrossSection is a variant of the Barashenkov\n"
101
<<
"cross section parameterization to be used of protons and\n"
102
<<
"nucleons on targets heavier than hydrogen. It is intended for\n"
103
<<
"use as a cross section component and is currently used by\n"
104
<<
"G4BGGNucleonInelasticXS. It is valid for incident energies up\n"
105
<<
"to 1 TeV.\n"
;
106
}
107
108
////////////////////////////////////////////////////////////////////////////
109
G4ComponentBarNucleonNucleusXsc.hh
G4DynamicParticle.hh
G4Neutron.hh
G4NucleonNuclearCrossSection.hh
G4ParticleDefinition.hh
G4Proton.hh
G4double
double G4double
Definition:
G4Types.hh:83
G4bool
bool G4bool
Definition:
G4Types.hh:86
G4int
int G4int
Definition:
G4Types.hh:85
Z
const G4int Z[17]
Definition:
G4WaterStopping.cc:51
G4ComponentBarNucleonNucleusXsc
Definition:
G4ComponentBarNucleonNucleusXsc.hh:50
G4ComponentBarNucleonNucleusXsc::GetElasticXsc
G4double GetElasticXsc()
Definition:
G4ComponentBarNucleonNucleusXsc.hh:92
G4ComponentBarNucleonNucleusXsc::BuildPhysicsTable
void BuildPhysicsTable(const G4ParticleDefinition &) final
Definition:
G4ComponentBarNucleonNucleusXsc.cc:202
G4ComponentBarNucleonNucleusXsc::GetInelasticXsc
G4double GetInelasticXsc()
Definition:
G4ComponentBarNucleonNucleusXsc.hh:93
G4ComponentBarNucleonNucleusXsc::GetTotalXsc
G4double GetTotalXsc()
Definition:
G4ComponentBarNucleonNucleusXsc.hh:91
G4ComponentBarNucleonNucleusXsc::ComputeCrossSections
void ComputeCrossSections(const G4ParticleDefinition *aParticle, G4double kinEnergy, G4int Z)
Definition:
G4ComponentBarNucleonNucleusXsc.cc:142
G4DynamicParticle
Definition:
G4DynamicParticle.hh:65
G4DynamicParticle::GetDefinition
G4ParticleDefinition * GetDefinition() const
G4DynamicParticle::GetKineticEnergy
G4double GetKineticEnergy() const
G4Material
Definition:
G4Material.hh:117
G4Neutron::Neutron
static G4Neutron * Neutron()
Definition:
G4Neutron.cc:103
G4NucleonNuclearCrossSection::BuildPhysicsTable
void BuildPhysicsTable(const G4ParticleDefinition &) final
Definition:
G4NucleonNuclearCrossSection.cc:90
G4NucleonNuclearCrossSection::GetElementCrossSection
G4double GetElementCrossSection(const G4DynamicParticle *aParticle, G4int Z, const G4Material *mat=nullptr) final
Definition:
G4NucleonNuclearCrossSection.cc:69
G4NucleonNuclearCrossSection::G4NucleonNuclearCrossSection
G4NucleonNuclearCrossSection()
Definition:
G4NucleonNuclearCrossSection.cc:44
G4NucleonNuclearCrossSection::CrossSectionDescription
void CrossSectionDescription(std::ostream &) const final
Definition:
G4NucleonNuclearCrossSection.cc:98
G4NucleonNuclearCrossSection::~G4NucleonNuclearCrossSection
~G4NucleonNuclearCrossSection() override
Definition:
G4NucleonNuclearCrossSection.cc:56
G4NucleonNuclearCrossSection::IsElementApplicable
G4bool IsElementApplicable(const G4DynamicParticle *aParticle, G4int Z, const G4Material *mat) final
Definition:
G4NucleonNuclearCrossSection.cc:61
G4ParticleDefinition
Definition:
G4ParticleDefinition.hh:61
G4Proton::Proton
static G4Proton * Proton()
Definition:
G4Proton.cc:92
G4VCrossSectionDataSet
Definition:
G4VCrossSectionDataSet.hh:70
geant4-v11.1.1
source
processes
hadronic
cross_sections
src
G4NucleonNuclearCrossSection.cc
Generated by
1.9.6