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
G4GammaNuclearXS.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//
29// GEANT4 Class header file
30//
31//
32// File name: G4GammaNuclearXS
33//
34// Authors V.Ivantchenko, Geant4, 20 October 2020
35// B.Kutsenko, BINP/NSU, 10 August 2021
36//
37// Modification
38
39
40// Class Description:
41// This is a base class for gamma-nuclear cross section based on
42// data files from IAEA Evaluated Photonuclear Data Library (IAEA/PD-2019)
43// https://www-nds.iaea.org/photonuclear/
44// Database review - https://www.sciencedirect.com/science/article/pii/S0090375219300699?via%3Dihub
45// Class Description - End
46
47#ifndef G4GammaNuclearXS_h
48#define G4GammaNuclearXS_h 1
49
51#include "globals.hh"
52#include "G4ElementData.hh"
53#include "G4PhysicsVector.hh"
54#include "G4Threading.hh"
55#include "G4IsotopeList.hh"
56#include <vector>
57
60class G4Element;
62
64{
65public:
66
67 explicit G4GammaNuclearXS();
68
69 ~G4GammaNuclearXS() final;
70
71 static const char* Default_Name() {return "GammaNuclearXS";}
72
74 G4int Z, const G4Material*) final;
75
77 const G4Element*, const G4Material* mat) final;
78
80 G4int Z,
81 const G4Material* mat = nullptr) final;
82
84 const G4Isotope* iso = nullptr,
85 const G4Element* elm = nullptr,
86 const G4Material* mat = nullptr) final;
87
88 const G4Isotope* SelectIsotope(const G4Element*,
89 G4double kinEnergy, G4double logE) final;
90
91 void BuildPhysicsTable(const G4ParticleDefinition&) final;
92
94
96
97 void CrossSectionDescription(std::ostream&) const final;
98
99 G4GammaNuclearXS & operator=(const G4GammaNuclearXS &right) = delete;
101
102private:
103
104 void Initialise(G4int Z);
105
106 void InitialiseOnFly(G4int Z);
107
108 const G4String& FindDirectoryPath();
109
110 inline G4PhysicsVector* GetPhysicsVector(G4int Z);
111
112 G4PhysicsVector* RetrieveVector(std::ostringstream& in, G4bool isElement, G4int Z);
113
114 G4VCrossSectionDataSet* ggXsection = nullptr;
115
116 std::vector<G4double> temp;
117 const G4ParticleDefinition* gamma;
118
119 G4bool isMaster = false;
120
121 static const G4int MAXZGAMMAXS = 95;
122 static G4ElementData* data;
123 // Upper limit of the linear transition between IAEA database and CHIPS model
124 static const G4int rTransitionBound = 150.*CLHEP::MeV;
125 // The list of elements with non-linear parametrisation for better precision
126 const G4int freeVectorException[11] = {4, 6, 7, 8, 27, 39, 45, 65, 67, 69, 73};
127 // CHIPS photonuclear model had a problem with high energy parametrisation
128 // for isotopes of H and He, coefficient is needed to connect isotope cross
129 // section with element cross-section on high energies.
130 static G4double coeff[3][3];
131 static G4double xs150[MAXZGAMMAXS];
132 static G4String gDataDirectory;
133
134#ifdef G4MULTITHREADED
135 static G4Mutex gNuclearXSMutex;
136#endif
137};
138
139inline
140G4PhysicsVector* G4GammaNuclearXS::GetPhysicsVector(G4int Z)
141{
142 G4PhysicsVector* pv = data->GetElementData(Z);
143 if(pv == nullptr) {
144 InitialiseOnFly(Z);
145 pv = data->GetElementData(Z);
146 }
147 return pv;
148}
149
150#endif
std::mutex G4Mutex
Definition: G4Threading.hh:81
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]
G4PhysicsVector * GetElementData(G4int Z)
G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=nullptr, const G4Element *elm=nullptr, const G4Material *mat=nullptr) final
static const char * Default_Name()
void BuildPhysicsTable(const G4ParticleDefinition &) final
G4double ElementCrossSection(G4double ekin, G4int Z)
G4bool IsIsoApplicable(const G4DynamicParticle *, G4int Z, G4int A, const G4Element *, const G4Material *mat) final
G4double GetElementCrossSection(const G4DynamicParticle *, G4int Z, const G4Material *mat=nullptr) final
void CrossSectionDescription(std::ostream &) const final
G4double IsoCrossSection(G4double ekin, G4int Z, G4int A)
const G4Isotope * SelectIsotope(const G4Element *, G4double kinEnergy, G4double logE) final
G4bool IsElementApplicable(const G4DynamicParticle *, G4int Z, const G4Material *) final
#define inline
Definition: internal.h:104
Definition: DoubConv.h:17