Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4CrossSectionHP.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// V. Ivanchenko September 2023
27//
28// G4CrossSectionHP is a generic class implementing
29// cross sections for neutrons, protons and light ions
30// It is an alternative to code developed by J.P. Wellisch & T.Koi
31//
32
33#ifndef G4CrossSectionHP_h
34#define G4CrossSectionHP_h 1
35
37#include "globals.hh"
38#include "G4ElementData.hh"
39#include "G4PhysicsVector.hh"
40#include "G4LorentzVector.hh"
41#include "G4ThreeVector.hh"
42#include <vector>
43
47class G4Element;
48class G4Material;
49
51{
52 public:
54 const G4String& nameData,
55 const G4String& nameDir, G4double emaxHP,
56 G4int zmin, G4int zmax);
57
58 ~G4CrossSectionHP() override = default;
59
61 const G4Element*, const G4Material*) override;
62
64 const G4Isotope*, const G4Element*,
65 const G4Material*) override;
66
69 G4int Z, G4int A,
70 const G4Isotope* iso,
71 const G4Element* elm,
72 const G4Material* mat) override;
73
74 const G4Isotope* SelectIsotope(const G4Element*, G4double kinEnergy,
75 G4double logE) override;
76
77 void BuildPhysicsTable(const G4ParticleDefinition&) override;
78
79 void DumpPhysicsTable(const G4ParticleDefinition&) override;
80
83
84 protected:
85
86 inline G4double GetMaxHPEnergy() const;
87
88 inline void SetBinSearch(G4int n);
89
90 private:
91
92 void Initialise(const G4int Z);
93
94 void PrepareCache(const G4Material*);
95
96 G4double IsoCrossSection(const G4double kinE, const G4double loge,
97 const G4int Z, const G4int A,
98 const G4double temperature);
99
100 inline G4double GetCrossSection(const G4int Z, const G4int A);
101
102 inline G4bool CheckCache(const G4int Z);
103
104 const G4ParticleDefinition* fParticle;
105 const G4ParticleDefinition* fNeutron;
106 G4ParticleHPManager* fManagerHP;
107
108 const G4double emax;
109 const G4double emaxT;
110 const G4double elimit;
111 const G4double logElimit;
112 G4LorentzVector fLV;
113 G4ThreeVector fBoost;
114
115 const G4int minZ;
116 const G4int maxZ;
117
118 G4int binSearch{2};
119 std::size_t index{0};
120 G4bool fPrinted{false};
121
122 // cache
123 const G4Material* fCurrentMat{nullptr};
124 std::vector<std::pair<G4int, G4int> > fZA;
125 std::vector<G4double> fIsoXS;
126 std::vector<G4double> fTemp;
127
128 const G4String fDataName;
129 const G4String fDataDirectory;
130 G4ElementData* fData{nullptr};
131};
132
133inline G4double
134G4CrossSectionHP::GetCrossSection(const G4int Z, const G4int A)
135{
136 G4double res = 0.0;
137 for (std::size_t i=0; i<fZA.size(); ++i) {
138 if (Z == fZA[i].first && A == fZA[i].second) {
139 res = fIsoXS[i];
140 break;
141 }
142 }
143 return res;
144}
145
146inline G4bool G4CrossSectionHP::CheckCache(const G4int Z)
147{
148 G4bool res = false;
149 for (auto const & p : fZA) {
150 if (Z == p.first) {
151 res = true;
152 break;
153 }
154 }
155 return res;
156}
157
159{
160 return emax;
161}
162
164{
165 if (n > 0) { binSearch = n; }
166}
167
168#endif
CLHEP::HepLorentzVector G4LorentzVector
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
const G4double A[17]
G4CrossSectionHP & operator=(const G4CrossSectionHP &right)=delete
G4double ComputeIsoCrossSection(G4double kinEnergy, G4double loge, const G4ParticleDefinition *, G4int Z, G4int A, const G4Isotope *iso, const G4Element *elm, const G4Material *mat) override
G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *, const G4Element *, const G4Material *) override
~G4CrossSectionHP() override=default
G4CrossSectionHP(const G4CrossSectionHP &)=delete
void BuildPhysicsTable(const G4ParticleDefinition &) override
G4CrossSectionHP(const G4ParticleDefinition *, const G4String &nameData, const G4String &nameDir, G4double emaxHP, G4int zmin, G4int zmax)
void SetBinSearch(G4int n)
void DumpPhysicsTable(const G4ParticleDefinition &) override
G4bool IsIsoApplicable(const G4DynamicParticle *, G4int Z, G4int A, const G4Element *, const G4Material *) override
const G4Isotope * SelectIsotope(const G4Element *, G4double kinEnergy, G4double logE) override
G4double GetMaxHPEnergy() const
G4VCrossSectionDataSet(const G4String &nam="")