Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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// Author V.Ivantchenko, Geant4, 22 October 2020
35//
36
37// Class Description:
38// This is a base class for gamma-nuclear cross section based on
39// data files from G4PARTICLEXSDATA data set
40// Class Description - End
41
42#ifndef G4GammaNuclearXS_h
43#define G4GammaNuclearXS_h 1
44
46#include "globals.hh"
47#include "G4PhysicsVector.hh"
48#include "G4Threading.hh"
49#include <vector>
50
51const G4int MAXZEL = 93;
52
55class G4Element;
56
58{
59public:
60
61 explicit G4GammaNuclearXS();
62
63 ~G4GammaNuclearXS() final;
64
65 static const char* Default_Name() {return "G4GammaNuclearXS";}
66
68 G4int Z, const G4Material*) final;
69
71 const G4Element*, const G4Material* mat) final;
72
74 G4int Z, const G4Material* mat) final;
75
77 const G4Isotope* iso,
78 const G4Element* elm,
79 const G4Material* mat) final;
80
81 const G4Isotope* SelectIsotope(const G4Element*,
82 G4double kinEnergy, G4double logE) final;
83
84 void BuildPhysicsTable(const G4ParticleDefinition&) final;
85
86 void CrossSectionDescription(std::ostream&) const final;
87
90
91private:
92
93 void Initialise(G4int Z);
94
95 void InitialiseOnFly(G4int Z);
96
97 const G4String& FindDirectoryPath();
98
99 inline G4PhysicsVector* GetPhysicsVector(G4int Z);
100
101 G4VCrossSectionDataSet* ggXsection;
102 const G4ParticleDefinition* gamma;
103
104 static G4PhysicsVector* data[MAXZEL];
105 static G4double coeff[MAXZEL];
106 static G4String gDataDirectory;
107
108 G4bool isMaster;
109
110#ifdef G4MULTITHREADED
111 static G4Mutex gNuclearXSMutex;
112#endif
113};
114
115inline
116G4PhysicsVector* G4GammaNuclearXS::GetPhysicsVector(G4int Z)
117{
118 if(nullptr == data[Z]) { InitialiseOnFly(Z); }
119 return data[Z];
120}
121
122#endif
const G4int MAXZEL
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
static const char * Default_Name()
G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso, const G4Element *elm, const G4Material *mat) final
G4GammaNuclearXS(const G4GammaNuclearXS &)=delete
void BuildPhysicsTable(const G4ParticleDefinition &) final
G4GammaNuclearXS & operator=(const G4GammaNuclearXS &right)=delete
G4bool IsIsoApplicable(const G4DynamicParticle *, G4int Z, G4int A, const G4Element *, const G4Material *mat) final
void CrossSectionDescription(std::ostream &) const final
const G4Isotope * SelectIsotope(const G4Element *, G4double kinEnergy, G4double logE) final
G4double GetElementCrossSection(const G4DynamicParticle *, G4int Z, const G4Material *mat) final
G4bool IsElementApplicable(const G4DynamicParticle *, G4int Z, const G4Material *) final