Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4LivermoreGammaConversion5DModel.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// Author: Zhuxin Li@CENBG
27// 11 March 2020
28// on the base of G4LivermoreGammaConversionModel
29// derives from G4BetheHeitler5DModel
30// -------------------------------------------------------------------
31
33#include "G4Electron.hh"
34#include "G4Positron.hh"
35#include "G4EmParameters.hh"
38#include "G4PhysicsLogVector.hh"
41#include "G4SystemOfUnits.hh"
42#include "G4Exp.hh"
43
44//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
45//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
46
47G4double G4LivermoreGammaConversion5DModel::lowEnergyLimit = 2.*CLHEP::electron_mass_c2;
48G4int G4LivermoreGammaConversion5DModel::verboseLevel = 0;
49constexpr G4int G4LivermoreGammaConversion5DModel::maxZ;
50G4LPhysicsFreeVector* G4LivermoreGammaConversion5DModel::data[] = {nullptr};
51
53const G4String& nam)
54: G4BetheHeitler5DModel(p, nam), fParticleChange(nullptr)
55{
56 // Verbosity scale for debugging purposes:
57 // 0 = nothing
58 // 1 = calculation of cross sections, file openings...
59 // 2 = entering in methods
60 if(verboseLevel > 0)
61 {
62 G4cout << "G4LivermoreGammaConversion5DModel is constructed " << G4endl;
63 }
64}
65
66//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
67
69{
70 if(IsMaster()) {
71 for(G4int i=0; i<maxZ; ++i) {
72 if(data[i]) {
73 delete data[i];
74 data[i] = nullptr;
75 }
76 }
77 }
78}
79
80//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
81
82void
84 const G4DataVector& cuts)
85{
87
88 if (verboseLevel > 1)
89 {
90 G4cout << "Calling Initialise() of G4LivermoreGammaConversion5DModel."
91 << G4endl
92 << "Energy range: "
93 << LowEnergyLimit() / MeV << " MeV - "
94 << HighEnergyLimit() / GeV << " GeV isMater: " << IsMaster()
95 << G4endl;
96 }
97
98 if(!fParticleChange) {
99 fParticleChange = GetParticleChangeForGamma();
100 }
101
102 if(IsMaster())
103 {
104 // Initialise element selector
105 InitialiseElementSelectors(particle, cuts);
106 // Access to elements
107 char* path = std::getenv("G4LEDATA");
108 G4ProductionCutsTable* theCoupleTable =
110 G4int numOfCouples = theCoupleTable->GetTableSize();
111 for(G4int i=0; i<numOfCouples; ++i)
112 {
113 const G4MaterialCutsCouple* couple = theCoupleTable->GetMaterialCutsCouple(i);
114 SetCurrentCouple(couple);
115 const G4Material* mat = couple->GetMaterial();
116 const G4ElementVector* theElementVector = mat->GetElementVector();
117 G4int nelm = mat->GetNumberOfElements();
118 for (G4int j=0; j<nelm; ++j)
119 {
120 G4int Z = std::max(1, std::min((*theElementVector)[j]->GetZasInt(), maxZ));
121 if(!data[Z]) { ReadData(Z, path); }
122 }
123 }
124 }
125}
126
127//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
128
129void G4LivermoreGammaConversion5DModel::ReadData(size_t Z, const char* path)
130{
131 if (verboseLevel > 1)
132 {
133 G4cout << "Calling ReadData() of G4LivermoreGammaConversion5DModel"
134 << G4endl;
135 }
136
137 if(data[Z]) { return; }
138 const char* datadir = path;
139 if(!datadir)
140 {
141 datadir = std::getenv("G4LEDATA");
142 if(!datadir)
143 {
144 G4Exception("G4LivermoreGammaConversion5DModel::ReadData()",
145 "em0006",FatalException,
146 "Environment variable G4LEDATA not defined");
147 return;
148 }
149 }
150 data[Z] = new G4LPhysicsFreeVector();
151 std::ostringstream ost;
152 ost << datadir << "/epics2017/pair/pp-cs-" << Z <<".dat";
153 std::ifstream fin(ost.str().c_str());
154
155 if( !fin.is_open())
156 {
158 ed << "G4LivermoreGammaConversion5DModel data file <" << ost.str().c_str()
159 << "> is not opened!" << G4endl;
160 G4Exception("G4LivermoreGammaConversion5DModel::ReadData()",
161 "em0003",FatalException,
162 ed,"G4LEDATA version should be G4EMLOW6.27 or later.");
163 return;
164 }
165 else
166 {
167 if(verboseLevel > 1) { G4cout << "File " << ost.str()
168 << " is opened by G4LivermoreGammaConversion5DModel" << G4endl;}
169 data[Z]->Retrieve(fin, true);
170 }
171 // Activation of linear interpolation
172 data[Z]->SetSpline(false); // EPICS2017 has more points -> linear is fine
173}
174
175//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
176
179 const G4ParticleDefinition* particle, G4double GammaEnergy, G4double Z,
181{
182 if (verboseLevel > 1)
183 {
184 G4cout << "G4LivermoreGammaConversion5DModel::ComputeCrossSectionPerAtom() Z= "
185 << Z << G4endl;
186 }
187 G4double xs = 0.0;
188 if (GammaEnergy < lowEnergyLimit) { return xs; }
189
190 G4int intZ = std::max(1, std::min(G4lrint(Z), maxZ));
191 G4LPhysicsFreeVector* pv = data[intZ];
192 // if element was not initialised
193 // do initialisation safely for MT mode
194 if(!pv)
195 {
196 InitialiseForElement(particle, intZ);
197 pv = data[intZ];
198 if(!pv) { return xs; }
199 }
200 // x-section is taken from the table
201 xs = pv->Value(GammaEnergy);
202 if(verboseLevel > 0)
203 {
204 G4cout << "*** Gamma conversion xs for Z=" << Z << " at energy E(MeV)="
205 << GammaEnergy/MeV << " cs=" << xs/millibarn << " mb" << G4endl;
206 }
207 return xs;
208}
209
210//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
211
212#include "G4AutoLock.hh"
213namespace { G4Mutex LivermoreGammaConversion5DModelMutex = G4MUTEX_INITIALIZER; }
214
216 const G4ParticleDefinition*,
217 G4int Z)
218{
219 G4AutoLock l(&LivermoreGammaConversion5DModelMutex);
220 // G4cout << "G4LivermoreGammaConversion5DModel::InitialiseForElement Z= "
221 // << Z << G4endl;
222 if(!data[Z]) { ReadData(Z); }
223 l.unlock();
224}
225
226//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
std::vector< G4Element * > G4ElementVector
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
#define G4MUTEX_INITIALIZER
Definition: G4Threading.hh:85
std::mutex G4Mutex
Definition: G4Threading.hh:81
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
G4LivermoreGammaConversion5DModel(const G4ParticleDefinition *p=nullptr, const G4String &nam="Livermore5DConversion")
void InitialiseForElement(const G4ParticleDefinition *, G4int Z) override
G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0.0, G4double cut=0.0, G4double emax=DBL_MAX) override
const G4Material * GetMaterial() const
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:188
size_t GetNumberOfElements() const
Definition: G4Material.hh:184
G4double Value(G4double theEnergy, std::size_t &lastidx) const
virtual G4bool Retrieve(std::ifstream &fIn, G4bool ascii=false)
void SetSpline(G4bool)
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
std::size_t GetTableSize() const
static G4ProductionCutsTable * GetProductionCutsTable()
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:133
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:652
G4bool IsMaster() const
Definition: G4VEmModel.hh:736
void SetCurrentCouple(const G4MaterialCutsCouple *)
Definition: G4VEmModel.hh:465
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:645
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:148
int G4lrint(double ad)
Definition: templates.hh:134