Geant4
9.6.p02
Main Page
Related Pages
Modules
Namespaces
Classes
Files
File List
File Members
All
Classes
Namespaces
Files
Functions
Variables
Typedefs
Enumerations
Enumerator
Friends
Macros
Groups
Pages
geant4_9_6_p02
source
processes
hadronic
cross_sections
include
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
// $Id: $
27
// GEANT4 tag $Name: not supported by cvs2svn $
28
//
29
// -------------------------------------------------------------------
30
// File name: G4CrossSectionDataStore
31
//
32
// Modifications:
33
// 23.01.2009 V.Ivanchenko move constructor and destructor to source,
34
// use STL vector instead of C-array
35
//
36
// August 2011 Re-designed
37
// by G. Folger, V. Ivantchenko, T. Koi and D.H. Wright
38
39
// Class Description
40
// This is the class to which cross section data sets may be registered.
41
// An instance of it is contained in each hadronic process, allowing
42
// the use of the AddDataSet() method to tailor the cross sections to
43
// your application.
44
// Class Description - End
45
46
#ifndef G4CrossSectionDataStore_h
47
#define G4CrossSectionDataStore_h 1
48
49
#include "
globals.hh
"
50
#include "
G4VCrossSectionDataSet.hh
"
51
#include <vector>
52
53
class
G4Nucleus
;
54
class
G4DynamicParticle
;
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
();
68
69
// Cross section per unit volume is computed (inverse mean free path)
70
G4double
GetCrossSection
(
const
G4DynamicParticle
*,
const
G4Material
*);
71
72
// Cross section per element is computed
73
G4double
GetCrossSection
(
const
G4DynamicParticle
*,
74
const
G4Element
*,
const
G4Material
*);
75
76
// Cross section per isotope is computed
77
G4double
GetCrossSection
(
const
G4DynamicParticle
*,
G4int
Z
,
G4int
A,
78
const
G4Isotope
*,
79
const
G4Element
*,
const
G4Material
*);
80
81
// Sample Z and A of a target nucleus and upload into G4Nucleus
82
G4Element
*
SampleZandA
(
const
G4DynamicParticle
*,
const
G4Material
*,
83
G4Nucleus
&
target
);
84
85
// Initialisation before run
86
void
BuildPhysicsTable
(
const
G4ParticleDefinition
&);
87
88
// Dump store to G4cout
89
void
DumpPhysicsTable
(
const
G4ParticleDefinition
&);
90
91
// Dump store as html
92
void
DumpHtml
(
const
G4ParticleDefinition
&, std::ofstream&);
93
94
inline
void
AddDataSet
(
G4VCrossSectionDataSet
*);
95
96
inline
void
SetVerboseLevel
(
G4int
value
);
97
98
private
:
99
100
G4double
GetIsoCrossSection(
const
G4DynamicParticle
*,
G4int
Z,
G4int
A,
101
const
G4Isotope
*,
102
const
G4Element
*,
const
G4Material
* aMaterial,
103
G4int
index
);
104
105
G4CrossSectionDataStore
& operator=(
const
G4CrossSectionDataStore
&
right
);
106
G4CrossSectionDataStore
(
const
G4CrossSectionDataStore
&);
107
108
G4NistManager
* nist;
109
110
std::vector<G4VCrossSectionDataSet*> dataSetList;
111
std::vector<G4double> xsecelm;
112
std::vector<G4double> xseciso;
113
114
const
G4Material
* currentMaterial;
115
const
G4ParticleDefinition
* matParticle;
116
G4double
matKinEnergy;
117
G4double
matCrossSection;
118
119
const
G4Material
* elmMaterial;
120
const
G4Element
* currentElement;
121
const
G4ParticleDefinition
* elmParticle;
122
G4double
elmKinEnergy;
123
G4double
elmCrossSection;
124
125
G4int
nDataSetList;
126
G4int
verboseLevel;
127
};
128
129
inline
void
G4CrossSectionDataStore::AddDataSet
(
G4VCrossSectionDataSet
*
p
)
130
{
131
dataSetList.push_back(p);
132
++nDataSetList;
133
}
134
135
inline
void
G4CrossSectionDataStore::SetVerboseLevel
(
G4int
value
)
136
{
137
verboseLevel =
value
;
138
}
139
140
#endif
Generated on Sat May 25 2013 14:33:40 for Geant4 by
1.8.4