Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4VCrossSectionDataSet.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// -------------------------------------------------------------------
27//
28// GEANT4 Class file
29//
30// File name: G4VCrossSectionDataSet
31//
32// Author F.W. Jones, TRIUMF, 20-JAN-97
33//
34// Modifications:
35// 23.01.2009 V.Ivanchenko move constructor and destructor to source
36// 12.08.2011 G.Folger, V.Ivanchenko, T.Koi, D.Wright redesign the class
37//
38
40#include "G4SystemOfUnits.hh"
42#include "G4Material.hh"
43#include "G4Element.hh"
44#include "G4Isotope.hh"
45#include "G4NistManager.hh"
46#include "Randomize.hh"
48
50 verboseLevel(0),name(nam),minKinEnergy(0.0),
51 maxKinEnergy(G4HadronicParameters::Instance()->GetMaxEnergy()),
52 isForAllAtomsAndEnergies(false)
53{
55 registry->Register(this);
56}
57
59{
60 registry->DeRegister(this);
61}
62
63G4bool
65 G4int,
66 const G4Material*)
67{
68 return false;
69}
70
71G4bool
73 G4int, G4int,
74 const G4Element*,
75 const G4Material*)
76{
77 return false;
78}
79
82 const G4Element* elm,
83 const G4Material* mat)
84{
85 G4int Z = elm->GetZasInt();
86
87 if (IsElementApplicable(part, Z, mat)) {
88 return GetElementCrossSection(part, Z, mat);
89 }
90
91 // isotope-wise cross section making sum over available
92 // isotope cross sections, which may be incomplete, so
93 // the result is corrected
94 std::size_t nIso = elm->GetNumberOfIsotopes();
95 G4double fact = 0.0;
96 G4double xsec = 0.0;
97
98 // user-defined isotope abundances
99 const G4IsotopeVector* isoVector = elm->GetIsotopeVector();
100 const G4double* abundVector = elm->GetRelativeAbundanceVector();
101 for (std::size_t j=0; j<nIso; ++j) {
102 const G4Isotope* iso = (*isoVector)[j];
103 G4int A = iso->GetN();
104 if(abundVector[j] > 0.0 && IsIsoApplicable(part, Z, A, elm, mat)) {
105 fact += abundVector[j];
106 xsec += abundVector[j]*GetIsoCrossSection(part, Z, A, iso, elm, mat);
107 }
108 }
109 return (fact > 0.0) ? xsec/fact : 0.0;
110}
111
114 G4double kinEnergy, G4double loge,
115 const G4ParticleDefinition* pd,
116 const G4Element* elm, const G4Material* mat)
117{
118 G4int Z = elm->GetZasInt();
119 std::size_t nIso = elm->GetNumberOfIsotopes();
120 G4double xsec = 0.0;
121 const G4IsotopeVector* isoVector = elm->GetIsotopeVector();
122 const G4double* abundVector = elm->GetRelativeAbundanceVector();
123 for (std::size_t j=0; j<nIso; ++j) {
124 const G4Isotope* iso = (*isoVector)[j];
125 G4int A = iso->GetN();
126 xsec += abundVector[j]*
127 ComputeIsoCrossSection(kinEnergy, loge, pd, Z, A, iso, elm, mat);
128 }
129 return xsec;
130}
131
134 G4int Z,
135 const G4Material* mat)
136{
138 ed << "GetElementCrossSection is not implemented in <" << name << ">\n"
139 << "Particle: " << dynPart->GetDefinition()->GetParticleName()
140 << " Ekin(MeV)= " << dynPart->GetKineticEnergy()/MeV;
141 if(nullptr != mat) { ed << " material: " << mat->GetName(); }
142 ed << " target Z= " << Z << G4endl;
143 G4Exception("G4VCrossSectionDataSet::GetElementCrossSection", "had001",
144 FatalException, ed);
145 return 0.0;
146}
147
150 G4int Z, G4int A,
151 const G4Isotope*,
152 const G4Element* elm,
153 const G4Material* mat)
154{
156 ed << "GetIsoCrossSection is not implemented in <" << name << ">\n"
157 << "Particle: " << dynPart->GetDefinition()->GetParticleName()
158 << " Ekin(MeV)= " << dynPart->GetKineticEnergy()/MeV;
159 if(nullptr != mat) { ed << " material: " << mat->GetName(); }
160 if(nullptr != elm) { ed << " element: " << elm->GetName(); }
161 ed << " target Z= " << Z << " A= " << A << G4endl;
162 G4Exception("G4VCrossSectionDataSet::GetIsoCrossSection", "had001",
163 FatalException, ed);
164 return 0.0;
165}
166
169 const G4ParticleDefinition* pd,
170 G4int Z, G4int A,
171 const G4Isotope*,
172 const G4Element* elm,
173 const G4Material* mat)
174{
176 ed << "GetIsoCrossSection is not implemented in <" << name << ">\n"
177 << "Particle: " << pd->GetParticleName()
178 << " Ekin(MeV)= " << kinEnergy/CLHEP::MeV;
179 if(nullptr != mat) { ed << " material: " << mat->GetName(); }
180 if(nullptr != elm) { ed << " element: " << elm->GetName(); }
181 ed << " target Z= " << Z << " A= " << A << G4endl;
182 G4Exception("G4VCrossSectionDataSet::GetIsoCrossSection", "had001",
183 FatalException, ed);
184 return 0.0;
185}
186
187const G4Isotope*
190{
191 G4int nIso = (G4int)anElement->GetNumberOfIsotopes();
192 const G4Isotope* iso = anElement->GetIsotope(0);
193
194 // more than 1 isotope
195 if(1 < nIso) {
196 const G4double* abundVector = anElement->GetRelativeAbundanceVector();
197 G4double sum = 0.0;
199 for (G4int j=0; j<nIso; ++j) {
200 sum += abundVector[j];
201 if(q <= sum) {
202 iso = anElement->GetIsotope(j);
203 break;
204 }
205 }
206 }
207 return iso;
208}
209
211{}
212
214{}
215
216void G4VCrossSectionDataSet::CrossSectionDescription(std::ostream& outFile) const
217{
218 outFile << "The description for this cross section data set has not been written yet.\n";
219}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
std::vector< G4Isotope * > G4IsotopeVector
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
const G4double A[17]
#define G4endl
Definition: G4ios.hh:57
#define G4UniformRand()
Definition: Randomize.hh:52
void DeRegister(G4VCrossSectionDataSet *)
void Register(G4VCrossSectionDataSet *)
static G4CrossSectionDataSetRegistry * Instance()
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4double * GetRelativeAbundanceVector() const
Definition: G4Element.hh:167
const G4Isotope * GetIsotope(G4int iso) const
Definition: G4Element.hh:170
size_t GetNumberOfIsotopes() const
Definition: G4Element.hh:159
const G4String & GetName() const
Definition: G4Element.hh:127
G4int GetZasInt() const
Definition: G4Element.hh:132
G4IsotopeVector * GetIsotopeVector() const
Definition: G4Element.hh:163
G4int GetN() const
Definition: G4Isotope.hh:93
const G4String & GetName() const
Definition: G4Material.hh:172
const G4String & GetParticleName() const
G4VCrossSectionDataSet(const G4String &nam="")
virtual G4double GetElementCrossSection(const G4DynamicParticle *, G4int Z, const G4Material *mat=nullptr)
virtual G4double ComputeIsoCrossSection(G4double kinEnergy, G4double loge, const G4ParticleDefinition *, G4int Z, G4int A, const G4Isotope *iso=nullptr, const G4Element *elm=nullptr, const G4Material *mat=nullptr)
virtual G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=nullptr, const G4Element *elm=nullptr, const G4Material *mat=nullptr)
virtual void DumpPhysicsTable(const G4ParticleDefinition &)
virtual G4double ComputeCrossSectionPerElement(G4double kinEnergy, G4double loge, const G4ParticleDefinition *, const G4Element *, const G4Material *mat=nullptr)
G4double ComputeCrossSection(const G4DynamicParticle *, const G4Element *, const G4Material *mat=nullptr)
virtual G4bool IsIsoApplicable(const G4DynamicParticle *, G4int Z, G4int A, const G4Element *elm=nullptr, const G4Material *mat=nullptr)
virtual void CrossSectionDescription(std::ostream &) const
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
virtual G4bool IsElementApplicable(const G4DynamicParticle *, G4int Z, const G4Material *mat=nullptr)
virtual const G4Isotope * SelectIsotope(const G4Element *, G4double kinEnergy, G4double logE)