Geant4 9.6.0
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
35using namespace std;
36
37//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
38//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
39
41 :G4VEmModel("LivermoreRayleigh"),isInitialised(false),maxZ(100)
42{
43 fParticleChange = 0;
44 lowEnergyLimit = 10 * eV;
45
46 dataCS.resize(maxZ+1,0);
47
49
50 verboseLevel= 0;
51 // Verbosity scale for debugging purposes:
52 // 0 = nothing
53 // 1 = calculation of cross sections, file openings...
54 // 2 = entering in methods
55
56 if(verboseLevel > 0)
57 {
58 G4cout << "G4LivermoreRayleighModel is constructed " << G4endl;
59 }
60
61}
62
63//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
64
66{
67 for(G4int i=0; i<=maxZ; ++i) { delete dataCS[i]; }
68}
69
70//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
71
73 const G4DataVector& cuts)
74{
75 if (verboseLevel > 1)
76 {
77 G4cout << "Calling Initialise() of G4LivermoreRayleighModel." << G4endl
78 << "Energy range: "
79 << LowEnergyLimit() / eV << " eV - "
80 << HighEnergyLimit() / GeV << " GeV"
81 << G4endl;
82 }
83
84 // Initialise element selector
85 InitialiseElementSelectors(particle, cuts);
86
87 // Access to elements
88
89 char* path = getenv("G4LEDATA");
90 G4ProductionCutsTable* theCoupleTable =
92 G4int numOfCouples = theCoupleTable->GetTableSize();
93
94 for(G4int i=0; i<numOfCouples; ++i)
95 {
96 const G4MaterialCutsCouple* couple =
97 theCoupleTable->GetMaterialCutsCouple(i);
98 const G4Material* material = couple->GetMaterial();
99 const G4ElementVector* theElementVector = material->GetElementVector();
100 G4int nelm = material->GetNumberOfElements();
101
102 for (G4int j=0; j<nelm; ++j)
103 {
104 G4int Z = G4lrint((*theElementVector)[j]->GetZ());
105 if(Z < 1) { Z = 1; }
106 else if(Z > maxZ) { Z = maxZ; }
107 if( (!dataCS[Z]) ) { ReadData(Z, path); }
108 }
109 }
110 //
111
112 //fGenerator->PrintGeneratorInformation();
113
114 if(isInitialised) { return; }
115 fParticleChange = GetParticleChangeForGamma();
116 isInitialised = true;
117
118}
119
120//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
121
122void G4LivermoreRayleighModel::ReadData(size_t Z, const char* path)
123{
124 if (verboseLevel > 1)
125 {
126 G4cout << "Calling ReadData() of G4LivermoreRayleighModel"
127 << G4endl;
128 }
129
130 if(dataCS[Z]) { return; }
131
132 const char* datadir = path;
133
134 if(!datadir)
135 {
136 datadir = getenv("G4LEDATA");
137 if(!datadir)
138 {
139 G4Exception("G4LivermoreRayleighModelModel::ReadData()","em0006",
141 "Environment variable G4LEDATA not defined");
142 return;
143 }
144 }
145
146 //
147
148 dataCS[Z] = new G4LPhysicsFreeVector();
149
150 // Activation of spline interpolation
151 //dataCS[Z] ->SetSpline(true);
152
153 std::ostringstream ostCS;
154 ostCS << datadir << "/livermore/rayl/re-cs-" << Z <<".dat";
155 std::ifstream finCS(ostCS.str().c_str());
156
157 if( !finCS .is_open() )
158 {
160 ed << "G4LivermoreRayleighModel data file <" << ostCS.str().c_str()
161 << "> is not opened!" << G4endl;
162 G4Exception("G4LivermoreRayleighModel::ReadData()","em0003",FatalException,
163 ed,"G4LEDATA version should be G4EMLOW6.27 or later.");
164 return;
165 }
166 else
167 {
168 if(verboseLevel > 3) {
169 G4cout << "File " << ostCS.str()
170 << " is opened by G4LivermoreRayleighModel" << G4endl;
171 }
172 dataCS[Z]->Retrieve(finCS, true);
173 }
174}
175
176//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
177
180 G4double GammaEnergy,
183{
184 if (verboseLevel > 1)
185 {
186 G4cout << "Calling ComputeCrossSectionPerAtom() of G4LivermoreRayleighModel"
187 << G4endl;
188 }
189
190 if(GammaEnergy < lowEnergyLimit) { return 0.0; }
191
192 G4double xs = 0.0;
193
194 G4int intZ = G4lrint(Z);
195
196 if(intZ < 1 || intZ > maxZ) { return xs; }
197
198 G4LPhysicsFreeVector* pv = dataCS[intZ];
199
200 // element was not initialised
201 if(!pv)
202 {
203 char* path = getenv("G4LEDATA");
204 ReadData(intZ, path);
205 pv = dataCS[intZ];
206 if(!pv) { return xs; }
207 }
208
209 G4int n = pv->GetVectorLength() - 1;
210 G4double e = GammaEnergy/MeV;
211 if(e >= pv->Energy(n)) {
212 xs = (*pv)[n]/(e*e);
213 } else if(e >= pv->Energy(0)) {
214 xs = pv->Value(e)/(e*e);
215 }
216
217 if(verboseLevel > 0)
218 {
219 G4cout << "****** DEBUG: tcs value for Z=" << Z << " at energy (MeV)=" << e << G4endl;
220 G4cout << " cs (Geant4 internal unit)=" << xs << G4endl;
221 G4cout << " -> first E*E*cs value in CS data file (iu) =" << (*pv)[0] << G4endl;
222 G4cout << " -> last E*E*cs value in CS data file (iu) =" << (*pv)[n] << G4endl;
223 G4cout << "*********************************************************" << G4endl;
224 }
225 return xs;
226}
227
228//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
229
231 std::vector<G4DynamicParticle*>*,
232 const G4MaterialCutsCouple* couple,
233 const G4DynamicParticle* aDynamicGamma,
235{
236 if (verboseLevel > 1) {
237 G4cout << "Calling SampleSecondaries() of G4LivermoreRayleighModel"
238 << G4endl;
239 }
240 G4double photonEnergy0 = aDynamicGamma->GetKineticEnergy();
241
242 // absorption of low-energy gamma
243 if (photonEnergy0 <= lowEnergyLimit)
244 {
245 fParticleChange->ProposeTrackStatus(fStopAndKill);
246 fParticleChange->SetProposedKineticEnergy(0.);
247 fParticleChange->ProposeLocalEnergyDeposit(photonEnergy0);
248 return ;
249 }
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
258 G4ThreeVector photonDirection =
259 GetAngularDistribution()->SampleDirection(aDynamicGamma,
260 photonEnergy0, Z, couple->GetMaterial());
261 fParticleChange->ProposeMomentumDirection(photonDirection);
262}
263
264//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
std::vector< G4Element * > G4ElementVector
@ FatalException
@ fStopAndKill
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4double GetZ() const
Definition: G4Element.hh:131
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0, G4double cut=0, G4double emax=DBL_MAX)
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
const G4Material * GetMaterial() const
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:189
size_t GetNumberOfElements() const
Definition: G4Material.hh:185
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
G4double Value(G4double theEnergy)
size_t GetVectorLength() const
G4double Energy(size_t index) const
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
static G4ProductionCutsTable * GetProductionCutsTable()
virtual G4ThreeVector & SampleDirection(const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, const G4Material *)=0
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:508
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:109
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:529
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:522
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:459
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:515
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:123
void ProposeTrackStatus(G4TrackStatus status)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
int G4lrint(double ad)
Definition: templates.hh:163