Geant4-11
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#include "G4AutoLock.hh"
44
46
47//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
48//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
49
53
55 const G4String& nam)
56 : G4BetheHeitler5DModel(p, nam), fParticleChange(nullptr)
57{
58 verboseLevel = 0;
59 // Verbosity scale for debugging purposes:
60 // 0 = nothing
61 // 1 = calculation of cross sections, file openings...
62 // 2 = entering in methods
63 if(verboseLevel > 0)
64 {
65 G4cout << "G4LivermoreGammaConversion5DModel is constructed " << G4endl;
66 }
67}
68
69//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
70
72{
73 if(IsMaster()) {
74 for(G4int i=0; i<maxZ; ++i) {
75 if(data[i]) {
76 delete data[i];
77 data[i] = nullptr;
78 }
79 }
80 }
81}
82
83//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
84
85void
87 const G4DataVector& cuts)
88{
90
91 if (verboseLevel > 1)
92 {
93 G4cout << "Calling Initialise() of G4LivermoreGammaConversion5DModel."
94 << G4endl
95 << "Energy range: "
96 << LowEnergyLimit() / MeV << " MeV - "
97 << HighEnergyLimit() / GeV << " GeV isMater: " << IsMaster()
98 << G4endl;
99 }
100
101 if(!fParticleChange) {
103 }
104
105 if(IsMaster())
106 {
107 // Initialise element selector
108 InitialiseElementSelectors(particle, cuts);
109 // Access to elements
110 char* path = std::getenv("G4LEDATA");
111 G4ProductionCutsTable* theCoupleTable =
113 G4int numOfCouples = theCoupleTable->GetTableSize();
114 for(G4int i=0; i<numOfCouples; ++i)
115 {
116 const G4MaterialCutsCouple* couple = theCoupleTable->GetMaterialCutsCouple(i);
117 SetCurrentCouple(couple);
118 const G4Material* mat = couple->GetMaterial();
119 const G4ElementVector* theElementVector = mat->GetElementVector();
120 G4int nelm = mat->GetNumberOfElements();
121 for (G4int j=0; j<nelm; ++j)
122 {
123 G4int Z = std::max(1, std::min((*theElementVector)[j]->GetZasInt(), maxZ));
124 if(!data[Z]) { ReadData(Z, path); }
125 }
126 }
127 }
128}
129
130//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
131
132void G4LivermoreGammaConversion5DModel::ReadData(size_t Z, const char* path)
133{
134 if (verboseLevel > 1)
135 {
136 G4cout << "Calling ReadData() of G4LivermoreGammaConversion5DModel"
137 << G4endl;
138 }
139
140 if(data[Z]) { return; }
141 const char* datadir = path;
142 if(!datadir)
143 {
144 datadir = std::getenv("G4LEDATA");
145 if(!datadir)
146 {
147 G4Exception("G4LivermoreGammaConversion5DModel::ReadData()",
148 "em0006",FatalException,
149 "Environment variable G4LEDATA not defined");
150 return;
151 }
152 }
153 data[Z] = new G4PhysicsFreeVector();
154 std::ostringstream ost;
155 ost << datadir << "/epics2017/pair/pp-cs-" << Z <<".dat";
156 std::ifstream fin(ost.str().c_str());
157
158 if( !fin.is_open())
159 {
161 ed << "G4LivermoreGammaConversion5DModel data file <" << ost.str().c_str()
162 << "> is not opened!" << G4endl;
163 G4Exception("G4LivermoreGammaConversion5DModel::ReadData()",
164 "em0003",FatalException,
165 ed,"G4LEDATA version should be G4EMLOW6.27 or later.");
166 return;
167 }
168 else
169 {
170 if(verboseLevel > 1) { G4cout << "File " << ost.str()
171 << " is opened by G4LivermoreGammaConversion5DModel" << G4endl;}
172 data[Z]->Retrieve(fin, true);
173 }
174}
175
176//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
177
180 const G4ParticleDefinition* particle, G4double GammaEnergy, G4double Z,
182{
183 if (verboseLevel > 1)
184 {
185 G4cout << "G4LivermoreGammaConversion5DModel::ComputeCrossSectionPerAtom() Z= "
186 << Z << G4endl;
187 }
188 G4double xs = 0.0;
189 if (GammaEnergy < lowEnergyLimit) { return xs; }
190
191 G4int intZ = std::max(1, std::min(G4lrint(Z), maxZ));
192 G4PhysicsFreeVector* pv = data[intZ];
193 // if element was not initialised
194 // do initialisation safely for MT mode
195 if(!pv)
196 {
197 InitialiseForElement(particle, intZ);
198 pv = data[intZ];
199 if(!pv) { return xs; }
200 }
201 // x-section is taken from the table
202 xs = pv->Value(GammaEnergy);
203 if(verboseLevel > 0)
204 {
205 G4cout << "*** Gamma conversion xs for Z=" << Z << " at energy E(MeV)="
206 << GammaEnergy/MeV << " cs=" << xs/millibarn << " mb" << G4endl;
207 }
208 return xs;
209}
210
211//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
212
214 const G4ParticleDefinition*,
215 G4int Z)
216{
218 if(!data[Z]) { ReadData(Z); }
219 l.unlock();
220}
221
222//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
std::vector< const 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
static constexpr double millibarn
Definition: G4SIunits.hh:86
static constexpr double GeV
Definition: G4SIunits.hh:203
static constexpr double MeV
Definition: G4SIunits.hh:200
#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
const G4int Z[17]
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
void ReadData(size_t Z, const char *path=nullptr)
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:186
size_t GetNumberOfElements() const
Definition: G4Material.hh:182
G4bool Retrieve(std::ifstream &fIn, G4bool ascii=false)
G4double Value(const G4double energy, std::size_t &lastidx) const
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
std::size_t GetTableSize() const
static G4ProductionCutsTable * GetProductionCutsTable()
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:123
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:662
G4bool IsMaster() const
Definition: G4VEmModel.hh:746
void SetCurrentCouple(const G4MaterialCutsCouple *)
Definition: G4VEmModel.hh:472
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:655
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:138
static constexpr double electron_mass_c2
T max(const T t1, const T t2)
brief Return the largest of the two arguments
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
int G4lrint(double ad)
Definition: templates.hh:134