Geant4
11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4UPiNuclearCrossSection.hh
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
// Calculation of the total, elastic and inelastic cross-sections
27
// based on Barashenkov parametrisations of pion data
28
//
29
// 16.08.06 V.Ivanchenko - first implementation
30
// 22.01.07 V.Ivanchenko - add cross section interfaces with Z and A
31
// 05.03.07 V.Ivanchenko - add IfZAApplicable
32
//
33
34
#ifndef G4UPiNuclearCrossSection_h
35
#define G4UPiNuclearCrossSection_h
36
37
#include "
G4VCrossSectionDataSet.hh
"
38
#include "
G4DynamicParticle.hh
"
39
#include "
G4ParticleDefinition.hh
"
40
#include "
globals.hh
"
41
#include "
G4Threading.hh
"
42
43
class
G4PhysicsTable
;
44
45
class
G4UPiNuclearCrossSection
:
public
G4VCrossSectionDataSet
46
{
47
public
:
48
49
explicit
G4UPiNuclearCrossSection
();
50
51
~G4UPiNuclearCrossSection
()
override
;
52
53
G4bool
IsElementApplicable
(
const
G4DynamicParticle
* aParticle,
54
G4int
Z
,
const
G4Material
*)
final
;
55
56
inline
57
G4double
GetElasticCrossSection
(
const
G4DynamicParticle
* aParticle,
58
G4int
Z
,
G4int
A
)
const
;
59
60
inline
61
G4double
GetInelasticCrossSection
(
const
G4DynamicParticle
* aParticle,
62
G4int
Z
,
G4int
A
)
const
;
63
64
void
BuildPhysicsTable
(
const
G4ParticleDefinition
&)
final
;
65
66
void
DumpPhysicsTable
(
const
G4ParticleDefinition
&)
final
;
67
68
void
CrossSectionDescription
(std::ostream&)
const
final
;
69
70
private
:
71
72
G4double
Interpolate(
G4int
Z
,
G4int
A
,
G4double
ekin,
73
const
G4PhysicsTable
*)
const
;
74
75
void
AddDataSet(
const
G4String
& p,
const
G4double
* tot,
76
const
G4double
* in,
const
G4double
* e,
G4int
n);
77
78
void
LoadData();
79
80
const
G4ParticleDefinition
* piPlus;
81
const
G4ParticleDefinition
* piMinus;
82
83
static
const
G4int
NZ = 16;
84
static
G4int
theZ[NZ];
85
static
G4int
idxZ[93];
86
87
static
G4double
theA[NZ];
88
static
G4double
APower[93];
89
90
static
G4PhysicsTable
* piPlusElastic;
91
static
G4PhysicsTable
* piPlusInelastic;
92
static
G4PhysicsTable
* piMinusElastic;
93
static
G4PhysicsTable
* piMinusInelastic;
94
95
G4double
aPower;
96
G4double
elow;
97
98
G4bool
isMaster;
99
G4bool
spline;
100
101
#ifdef G4MULTITHREADED
102
static
G4Mutex
pionUXSMutex;
103
#endif
104
};
105
106
inline
G4double
107
G4UPiNuclearCrossSection::GetElasticCrossSection
(
108
const
G4DynamicParticle
* dp,
G4int
Z
,
G4int
A
)
const
109
{
110
const
G4PhysicsTable
* table =
111
(dp->
GetDefinition
() == piPlus) ? piPlusElastic : piMinusElastic;
112
return
Interpolate(
Z
,
A
, dp->
GetKineticEnergy
(), table);
113
}
114
115
inline
G4double
116
G4UPiNuclearCrossSection::GetInelasticCrossSection
(
117
const
G4DynamicParticle
* dp,
G4int
Z
,
G4int
A
)
const
118
{
119
const
G4PhysicsTable
* table =
120
(dp->
GetDefinition
() == piPlus) ? piPlusInelastic : piMinusInelastic;
121
return
Interpolate(
Z
,
A
, dp->
GetKineticEnergy
(), table);
122
}
123
124
#endif
G4DynamicParticle.hh
G4ParticleDefinition.hh
G4Threading.hh
G4Mutex
std::mutex G4Mutex
Definition:
G4Threading.hh:81
G4double
double G4double
Definition:
G4Types.hh:83
G4bool
bool G4bool
Definition:
G4Types.hh:86
G4int
int G4int
Definition:
G4Types.hh:85
G4VCrossSectionDataSet.hh
Z
const G4int Z[17]
Definition:
G4WaterStopping.cc:51
A
const G4double A[17]
Definition:
G4WaterStopping.cc:53
G4DynamicParticle
Definition:
G4DynamicParticle.hh:65
G4DynamicParticle::GetDefinition
G4ParticleDefinition * GetDefinition() const
G4DynamicParticle::GetKineticEnergy
G4double GetKineticEnergy() const
G4Material
Definition:
G4Material.hh:117
G4ParticleDefinition
Definition:
G4ParticleDefinition.hh:61
G4PhysicsTable
Definition:
G4PhysicsTable.hh:56
G4String
Definition:
G4String.hh:62
G4UPiNuclearCrossSection
Definition:
G4UPiNuclearCrossSection.hh:46
G4UPiNuclearCrossSection::GetElasticCrossSection
G4double GetElasticCrossSection(const G4DynamicParticle *aParticle, G4int Z, G4int A) const
Definition:
G4UPiNuclearCrossSection.hh:107
G4UPiNuclearCrossSection::~G4UPiNuclearCrossSection
~G4UPiNuclearCrossSection() override
Definition:
G4UPiNuclearCrossSection.cc:74
G4UPiNuclearCrossSection::G4UPiNuclearCrossSection
G4UPiNuclearCrossSection()
Definition:
G4UPiNuclearCrossSection.cc:61
G4UPiNuclearCrossSection::IsElementApplicable
G4bool IsElementApplicable(const G4DynamicParticle *aParticle, G4int Z, const G4Material *) final
Definition:
G4UPiNuclearCrossSection.cc:101
G4UPiNuclearCrossSection::BuildPhysicsTable
void BuildPhysicsTable(const G4ParticleDefinition &) final
Definition:
G4UPiNuclearCrossSection.cc:173
G4UPiNuclearCrossSection::CrossSectionDescription
void CrossSectionDescription(std::ostream &) const final
Definition:
G4UPiNuclearCrossSection.cc:594
G4UPiNuclearCrossSection::DumpPhysicsTable
void DumpPhysicsTable(const G4ParticleDefinition &) final
Definition:
G4UPiNuclearCrossSection.cc:158
G4UPiNuclearCrossSection::GetInelasticCrossSection
G4double GetInelasticCrossSection(const G4DynamicParticle *aParticle, G4int Z, G4int A) const
Definition:
G4UPiNuclearCrossSection.hh:116
G4VCrossSectionDataSet
Definition:
G4VCrossSectionDataSet.hh:70
globals.hh
geant4-v11.1.1
source
processes
hadronic
cross_sections
include
G4UPiNuclearCrossSection.hh
Generated by
1.9.6