Geant4
11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4CrossSectionDataStore.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
//
27
// -------------------------------------------------------------------
28
// File name: G4CrossSectionDataStore
29
//
30
// Modifications:
31
// 23.01.2009 V.Ivanchenko move constructor and destructor to source,
32
// use STL vector instead of C-array
33
//
34
// August 2011 Re-designed
35
// by G. Folger, V. Ivantchenko, T. Koi and D.H. Wright
36
37
// Class Description
38
// This is the class to which cross section data sets may be registered.
39
// An instance of it is contained in each hadronic process, allowing
40
// the use of the AddDataSet() method to tailor the cross sections to
41
// your application.
42
// Class Description - End
43
44
#ifndef G4CrossSectionDataStore_h
45
#define G4CrossSectionDataStore_h 1
46
47
#include "
globals.hh
"
48
#include "
G4VCrossSectionDataSet.hh
"
49
#include "
G4DynamicParticle.hh
"
50
#include "
G4PhysicsVector.hh
"
51
#include <vector>
52
#include <iostream>
53
54
class
G4Nucleus
;
55
class
G4ParticleDefinition
;
56
class
G4Isotope
;
57
class
G4Element
;
58
class
G4Material
;
59
class
G4NistManager
;
60
61
class
G4CrossSectionDataStore
62
{
63
public
:
64
65
G4CrossSectionDataStore
();
66
67
~G4CrossSectionDataStore
() =
default
;
68
69
// Run time cross section per unit volume
70
inline
G4double
GetCrossSection
(
const
G4DynamicParticle
*,
const
G4Material
*);
71
G4double
ComputeCrossSection
(
const
G4DynamicParticle
*,
const
G4Material
*);
72
73
// Cross section per element
74
G4double
GetCrossSection
(
const
G4DynamicParticle
*,
75
const
G4Element
*,
const
G4Material
*);
76
77
// Cross section per isotope
78
G4double
GetCrossSection
(
const
G4DynamicParticle
*,
G4int
Z
,
G4int
A
,
79
const
G4Isotope
*,
80
const
G4Element
*,
const
G4Material
*);
81
82
// Sample Z and A of a target nucleus and upload into G4Nucleus
83
const
G4Element
*
SampleZandA
(
const
G4DynamicParticle
*,
const
G4Material
*,
84
G4Nucleus
& target);
85
86
// Initialisation before run
87
void
BuildPhysicsTable
(
const
G4ParticleDefinition
&);
88
89
// Dump store to G4cout
90
void
DumpPhysicsTable
(
const
G4ParticleDefinition
&);
91
92
// Dump store as html
93
void
DumpHtml
(
const
G4ParticleDefinition
&, std::ofstream&)
const
;
94
void
PrintCrossSectionHtml
(
const
G4VCrossSectionDataSet
*cs)
const
;
95
96
void
AddDataSet
(
G4VCrossSectionDataSet
*);
97
void
AddDataSet
(
G4VCrossSectionDataSet
*, std::size_t);
98
inline
const
std::vector<G4VCrossSectionDataSet*>&
GetDataSetList
()
const
;
99
100
inline
void
SetVerboseLevel
(
G4int
value);
101
102
// may be used by special processes
103
inline
void
SetForcedElement
(
const
G4Element
*);
104
105
G4CrossSectionDataStore
&
operator
=
106
(
const
G4CrossSectionDataStore
&right) =
delete
;
107
G4CrossSectionDataStore
(
const
G4CrossSectionDataStore
&) =
delete
;
108
109
private
:
110
111
G4double
GetIsoCrossSection(
const
G4DynamicParticle
*,
G4int
Z
,
G4int
A
,
112
const
G4Isotope
*,
const
G4Element
*,
113
const
G4Material
*,
const
G4int
index);
114
115
G4String
HtmlFileName(
const
G4String
& in)
const
;
116
117
G4NistManager
* nist;
118
const
G4Material
* currentMaterial =
nullptr
;
119
const
G4ParticleDefinition
* matParticle =
nullptr
;
120
const
G4Element
* forcedElement =
nullptr
;
121
G4double
matKinEnergy = 0.0;
122
G4double
matCrossSection = 0.0;
123
124
G4int
nDataSetList = 0;
125
G4int
verboseLevel = 1;
126
127
std::vector<G4VCrossSectionDataSet*> dataSetList;
128
std::vector<G4double> xsecelm;
129
std::vector<G4double> xseciso;
130
};
131
132
inline
void
G4CrossSectionDataStore::SetVerboseLevel
(
G4int
value)
133
{
134
verboseLevel = value;
135
}
136
137
inline
void
G4CrossSectionDataStore::SetForcedElement
(
const
G4Element
* ptr)
138
{
139
forcedElement = ptr;
140
}
141
142
inline
const
std::vector<G4VCrossSectionDataSet*>&
143
G4CrossSectionDataStore::GetDataSetList
()
const
144
{
145
return
dataSetList;
146
}
147
148
inline
G4double
149
G4CrossSectionDataStore::GetCrossSection
(
const
G4DynamicParticle
* dp,
150
const
G4Material
* mat)
151
{
152
if
(dp->
GetKineticEnergy
() != matKinEnergy || mat != currentMaterial ||
153
dp->
GetDefinition
() != matParticle) {
154
ComputeCrossSection
(dp, mat);
155
}
156
return
matCrossSection;
157
}
158
159
#endif
G4DynamicParticle.hh
G4PhysicsVector.hh
G4double
double G4double
Definition:
G4Types.hh:83
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
G4CrossSectionDataStore
Definition:
G4CrossSectionDataStore.hh:62
G4CrossSectionDataStore::DumpHtml
void DumpHtml(const G4ParticleDefinition &, std::ofstream &) const
Definition:
G4CrossSectionDataStore.cc:331
G4CrossSectionDataStore::~G4CrossSectionDataStore
~G4CrossSectionDataStore()=default
G4CrossSectionDataStore::GetDataSetList
const std::vector< G4VCrossSectionDataSet * > & GetDataSetList() const
Definition:
G4CrossSectionDataStore.hh:143
G4CrossSectionDataStore::BuildPhysicsTable
void BuildPhysicsTable(const G4ParticleDefinition &)
Definition:
G4CrossSectionDataStore.cc:272
G4CrossSectionDataStore::AddDataSet
void AddDataSet(G4VCrossSectionDataSet *)
Definition:
G4CrossSectionDataStore.cc:398
G4CrossSectionDataStore::G4CrossSectionDataStore
G4CrossSectionDataStore(const G4CrossSectionDataStore &)=delete
G4CrossSectionDataStore::PrintCrossSectionHtml
void PrintCrossSectionHtml(const G4VCrossSectionDataSet *cs) const
Definition:
G4CrossSectionDataStore.cc:362
G4CrossSectionDataStore::DumpPhysicsTable
void DumpPhysicsTable(const G4ParticleDefinition &)
Definition:
G4CrossSectionDataStore.cc:304
G4CrossSectionDataStore::ComputeCrossSection
G4double ComputeCrossSection(const G4DynamicParticle *, const G4Material *)
Definition:
G4CrossSectionDataStore.cc:68
G4CrossSectionDataStore::SetVerboseLevel
void SetVerboseLevel(G4int value)
Definition:
G4CrossSectionDataStore.hh:132
G4CrossSectionDataStore::GetCrossSection
G4double GetCrossSection(const G4DynamicParticle *, const G4Material *)
Definition:
G4CrossSectionDataStore.hh:149
G4CrossSectionDataStore::G4CrossSectionDataStore
G4CrossSectionDataStore()
Definition:
G4CrossSectionDataStore.cc:61
G4CrossSectionDataStore::SetForcedElement
void SetForcedElement(const G4Element *)
Definition:
G4CrossSectionDataStore.hh:137
G4CrossSectionDataStore::SampleZandA
const G4Element * SampleZandA(const G4DynamicParticle *, const G4Material *, G4Nucleus &target)
Definition:
G4CrossSectionDataStore.cc:192
G4DynamicParticle
Definition:
G4DynamicParticle.hh:65
G4DynamicParticle::GetDefinition
G4ParticleDefinition * GetDefinition() const
G4DynamicParticle::GetKineticEnergy
G4double GetKineticEnergy() const
G4Element
Definition:
G4Element.hh:98
G4Isotope
Definition:
G4Isotope.hh:72
G4Material
Definition:
G4Material.hh:117
G4NistManager
Definition:
G4NistManager.hh:83
G4Nucleus
Definition:
G4Nucleus.hh:52
G4ParticleDefinition
Definition:
G4ParticleDefinition.hh:61
G4String
Definition:
G4String.hh:62
G4VCrossSectionDataSet
Definition:
G4VCrossSectionDataSet.hh:70
globals.hh
geant4-v11.1.1
source
processes
hadronic
cross_sections
include
G4CrossSectionDataStore.hh
Generated by
1.9.6