Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4LivermoreRayleighModel.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: Sebastien Incerti
27// 31 March 2012
28// on base of G4LivermoreRayleighModel
29//
30
32#include "G4SystemOfUnits.hh"
34#include "G4EmParameters.hh"
35#include "G4AutoLock.hh"
36
37//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
38
39using namespace std;
40namespace { G4Mutex LivermoreRayleighModelMutex = G4MUTEX_INITIALIZER; }
41
42//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
43
44G4PhysicsFreeVector* G4LivermoreRayleighModel::dataCS[] = {nullptr};
45
47 :G4VEmModel("LivermoreRayleigh"),maxZ(100),isInitialised(false)
48{
49 fParticleChange = nullptr;
50 lowEnergyLimit = 10 * CLHEP::eV;
51
53
54 verboseLevel= 0;
55 // Verbosity scale for debugging purposes:
56 // 0 = nothing
57 // 1 = calculation of cross sections, file openings...
58 // 2 = entering in methods
59
60 if(verboseLevel > 0)
61 {
62 G4cout << "G4LivermoreRayleighModel is constructed " << G4endl;
63 }
64
65}
66
67//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
68
70{
71 if(IsMaster())
72 {
73 for(G4int i = 0; i <= maxZ; ++i)
74 {
75 if(dataCS[i])
76 {
77 delete dataCS[i];
78 dataCS[i] = nullptr;
79 }
80 }
81 }
82}
83
84//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
85
87 const G4DataVector& cuts)
88{
89 if (verboseLevel > 1)
90 {
91 G4cout << "Calling Initialise() of G4LivermoreRayleighModel." << G4endl
92 << "Energy range: "
93 << LowEnergyLimit() / eV << " eV - "
94 << HighEnergyLimit() / GeV << " GeV"
95 << G4endl;
96 }
97
98 if(IsMaster()) {
99 // Initialise element selector
100 InitialiseElementSelectors(particle, cuts);
101
102 // Access to elements
103 const char* path = G4FindDataDir("G4LEDATA");
104 const G4ElementTable* elemTable = G4Element::GetElementTable();
105 std::size_t numElems = (*elemTable).size();
106 for(std::size_t ie = 0; ie < numElems; ++ie)
107 {
108 const G4Element* elem = (*elemTable)[ie];
109 const G4int Z = std::min(maxZ, elem->GetZasInt());
110 if(dataCS[Z] == nullptr)
111 {
112 ReadData(Z, path);
113 }
114 }
115 }
116 if(isInitialised) { return; }
117 fParticleChange = GetParticleChangeForGamma();
118 isInitialised = true;
119}
120
121//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
122
124 G4VEmModel* masterModel)
125{
127}
128
129//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
130
131void G4LivermoreRayleighModel::ReadData(std::size_t Z, const char* path)
132{
133 if (verboseLevel > 1)
134 {
135 G4cout << "Calling ReadData() of G4LivermoreRayleighModel"
136 << G4endl;
137 }
138
139 if(nullptr != dataCS[Z]) { return; }
140
141 const char* datadir = path;
142
143 if(datadir == nullptr)
144 {
145 datadir = G4FindDataDir("G4LEDATA");
146 if(datadir == nullptr)
147 {
148 G4Exception("G4LivermoreRayleighModelModel::ReadData()","em0006",
150 "Environment variable G4LEDATA not defined");
151 return;
152 }
153 }
154 dataCS[Z] = new G4PhysicsFreeVector();
155
156 std::ostringstream ostCS;
157 if(G4EmParameters::Instance()->LivermoreDataDir() == "livermore"){
158 ostCS << datadir << "/livermore/rayl/re-cs-" << Z <<".dat";
159 }else{
160 ostCS << datadir << "/epics2017/rayl/re-cs-" << Z <<".dat";
161 }
162
163 std::ifstream finCS(ostCS.str().c_str());
164
165 if( !finCS .is_open() )
166 {
168 ed << "G4LivermoreRayleighModel data file <" << ostCS.str().c_str()
169 << "> is not opened!" << G4endl;
170 G4Exception("G4LivermoreRayleighModel::ReadData()","em0003",FatalException,
171 ed,"G4LEDATA version should be G4EMLOW8.0 or later.");
172 return;
173 }
174 else
175 {
176 if(verboseLevel > 3) {
177 G4cout << "File " << ostCS.str()
178 << " is opened by G4LivermoreRayleighModel" << G4endl;
179 }
180 dataCS[Z]->Retrieve(finCS, true);
181 }
182}
183
184//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
185
188 G4double GammaEnergy,
191{
192 if (verboseLevel > 1)
193 {
194 G4cout << "G4LivermoreRayleighModel::ComputeCrossSectionPerAtom()"
195 << G4endl;
196 }
197
198 if(GammaEnergy < lowEnergyLimit) { return 0.0; }
199
200 G4double xs = 0.0;
201 G4int intZ = G4lrint(Z);
202 if(intZ < 1 || intZ > maxZ) { return xs; }
203
204 G4PhysicsFreeVector* pv = dataCS[intZ];
205
206 // if element was not initialised
207 // do initialisation safely for MT mode
208 if(nullptr == pv) {
209 InitialiseForElement(0, intZ);
210 pv = dataCS[intZ];
211 if(nullptr == pv) { return xs; }
212 }
213
214 G4int n = G4int(pv->GetVectorLength() - 1);
215 G4double e = GammaEnergy/MeV;
216 if(e >= pv->Energy(n)) {
217 xs = (*pv)[n]/(e*e);
218 } else if(e >= pv->Energy(0)) {
219 xs = pv->Value(e)/(e*e);
220 }
221
222 if(verboseLevel > 1)
223 {
224 G4cout << "****** DEBUG: tcs value for Z=" << Z << " at energy (MeV)="
225 << e << G4endl;
226 G4cout << " cs (Geant4 internal unit)=" << xs << G4endl;
227 G4cout << " -> first E*E*cs value in CS data file (iu) =" << (*pv)[0]
228 << G4endl;
229 G4cout << " -> last E*E*cs value in CS data file (iu) =" << (*pv)[n]
230 << G4endl;
231 G4cout << "*********************************************************"
232 << G4endl;
233 }
234 return xs;
235}
236
237//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
238
240 std::vector<G4DynamicParticle*>*,
241 const G4MaterialCutsCouple* couple,
242 const G4DynamicParticle* aDynamicGamma,
244{
245 if (verboseLevel > 1) {
246 G4cout << "Calling SampleSecondaries() of G4LivermoreRayleighModel"
247 << G4endl;
248 }
249 G4double photonEnergy0 = aDynamicGamma->GetKineticEnergy();
250
251 // Select randomly one element in the current material
252 const G4ParticleDefinition* particle = aDynamicGamma->GetDefinition();
253 const G4Element* elm = SelectRandomAtom(couple,particle,photonEnergy0);
254 G4int Z = G4lrint(elm->GetZ());
255
256 // Sample the angle of the scattered photon
257 G4ThreeVector photonDirection =
258 GetAngularDistribution()->SampleDirection(aDynamicGamma,
259 photonEnergy0,
260 Z, couple->GetMaterial());
261 fParticleChange->ProposeMomentumDirection(photonDirection);
262}
263
264//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
265
266void
268 G4int Z)
269{
270 G4AutoLock l(&LivermoreRayleighModelMutex);
271 if(nullptr == dataCS[Z]) { ReadData(Z); }
272 l.unlock();
273}
274
275//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
std::vector< G4Element * > G4ElementTable
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 elem(i, j)
#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
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
static G4ElementTable * GetElementTable()
Definition: G4Element.cc:403
G4double GetZ() const
Definition: G4Element.hh:131
static G4EmParameters * Instance()
G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0, G4double cut=0, G4double emax=DBL_MAX) override
void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
void InitialiseForElement(const G4ParticleDefinition *, G4int Z) override
void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel) override
const G4Material * GetMaterial() const
void ProposeMomentumDirection(const G4ThreeVector &Pfinal)
G4double Energy(const std::size_t index) const
G4bool Retrieve(std::ifstream &fIn, G4bool ascii=false)
G4double Value(const G4double energy, std::size_t &lastidx) const
std::size_t GetVectorLength() const
virtual G4ThreeVector & SampleDirection(const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, const G4Material *)=0
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
Definition: G4VEmModel.hh:831
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:600
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:124
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:641
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:823
G4bool IsMaster() const
Definition: G4VEmModel.hh:725
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:634
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:561
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:607
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:139
int G4lrint(double ad)
Definition: templates.hh:134