Geant4 11.1.1
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#include "G4AutoLock.hh"
44
45namespace { G4Mutex LivermoreGammaConversion5DModelMutex = G4MUTEX_INITIALIZER; }
46
47//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
48//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
49
50const G4int G4LivermoreGammaConversion5DModel::maxZ;
51G4double G4LivermoreGammaConversion5DModel::lowEnergyLimit = 2.*CLHEP::electron_mass_c2;
52G4PhysicsFreeVector* G4LivermoreGammaConversion5DModel::data[] = {nullptr};
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) {
102 fParticleChange = GetParticleChangeForGamma();
103 }
104
105 if(IsMaster())
106 {
107 // Initialise element selector
108 InitialiseElementSelectors(particle, cuts);
109 // Access to elements
110 const char* path = G4FindDataDir("G4LEDATA");
111 G4ProductionCutsTable* theCoupleTable =
113 G4int numOfCouples = G4int(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 std::size_t nelm = mat->GetNumberOfElements();
121 for (std::size_t 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 = G4FindDataDir("G4LEDATA");
145 if(!datadir)
146 {
147 G4Exception("G4LivermoreGammaConversion5DModel::ReadData()",
148 "em0006",FatalException,
149 "Environment variable G4LEDATA not defined");
150 return;
151 }
152 }
153 std::ostringstream ost;
154 if(G4EmParameters::Instance()->LivermoreDataDir() == "livermore"){
155 data[Z] = new G4PhysicsFreeVector(true);
156 ost << datadir << "/livermore/pair/pp-cs-" << Z <<".dat";
157 }else{
158 data[Z] = new G4PhysicsFreeVector();
159 ost << datadir << "/epics2017/pair/pp-cs-" << Z <<".dat";
160 }
161
162 std::ifstream fin(ost.str().c_str());
163
164 if( !fin.is_open())
165 {
167 ed << "G4LivermoreGammaConversion5DModel data file <" << ost.str().c_str()
168 << "> is not opened!" << G4endl;
169 G4Exception("G4LivermoreGammaConversion5DModel::ReadData()",
170 "em0003",FatalException,
171 ed,"G4LEDATA version should be G4EMLOW8.0 or later.");
172 return;
173 }
174 else
175 {
176 if(verboseLevel > 1) { G4cout << "File " << ost.str()
177 << " is opened by G4LivermoreGammaConversion5DModel" << G4endl;}
178 data[Z]->Retrieve(fin, true);
179 }
180}
181
182//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
183
186 const G4ParticleDefinition* particle, G4double GammaEnergy, G4double Z,
188{
189 if (verboseLevel > 1)
190 {
191 G4cout << "G4LivermoreGammaConversion5DModel::ComputeCrossSectionPerAtom() Z= "
192 << Z << G4endl;
193 }
194 G4double xs = 0.0;
195 if (GammaEnergy < lowEnergyLimit) { return xs; }
196
197 G4int intZ = std::max(1, std::min(G4lrint(Z), maxZ));
198 G4PhysicsFreeVector* pv = data[intZ];
199 // if element was not initialised
200 // do initialisation safely for MT mode
201 if(!pv)
202 {
203 InitialiseForElement(particle, intZ);
204 pv = data[intZ];
205 if(!pv) { return xs; }
206 }
207 // x-section is taken from the table
208 xs = pv->Value(GammaEnergy);
209 if(verboseLevel > 0)
210 {
211 G4cout << "*** Gamma conversion xs for Z=" << Z << " at energy E(MeV)="
212 << GammaEnergy/MeV << " cs=" << xs/millibarn << " mb" << G4endl;
213 }
214 return xs;
215}
216
217//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
218
220 const G4ParticleDefinition*,
221 G4int Z)
222{
223 G4AutoLock l(&LivermoreGammaConversion5DModelMutex);
224 if(!data[Z]) { ReadData(Z); }
225 l.unlock();
226}
227
228//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
std::vector< const G4Element * > G4ElementVector
const char * G4FindDataDir(const char *)
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
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
const G4int Z[17]
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
static G4EmParameters * Instance()
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:185
size_t GetNumberOfElements() const
Definition: G4Material.hh:181
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:124
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:641
G4bool IsMaster() const
Definition: G4VEmModel.hh:725
void SetCurrentCouple(const G4MaterialCutsCouple *)
Definition: G4VEmModel.hh:468
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:634
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:139
int G4lrint(double ad)
Definition: templates.hh:134