Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
G4LivermoreComptonModel.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// $Id$
27//
28//
29// Author: Sebastien Incerti
30// 30 October 2008
31// on base of G4LowEnergyCompton developed by A.Forti and M.G.Pia
32//
33// History:
34// --------
35// 18 Apr 2009 V Ivanchenko Cleanup initialisation and generation of secondaries:
36// - apply internal high-energy limit only in constructor
37// - do not apply low-energy limit (default is 0)
38// - remove GetMeanFreePath method and table
39// - added protection against numerical problem in energy sampling
40// - use G4ElementSelector
41// 26 Dec 2010 V Ivanchenko Load data tables only once to avoid memory leak
42// 30 May 2011 V Ivanchenko Migration to model design for deexcitation
43
46#include "G4SystemOfUnits.hh"
47#include "G4Electron.hh"
49#include "G4LossTableManager.hh"
51#include "G4AtomicShell.hh"
55#include "G4Gamma.hh"
56
57//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
58
59using namespace std;
60
61//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
62
64 const G4String& nam)
65 :G4VEmModel(nam),fParticleChange(0),isInitialised(false),
66 scatterFunctionData(0),
67 crossSectionHandler(0),fAtomDeexcitation(0)
68{
69 lowEnergyLimit = 250 * eV;
70 highEnergyLimit = 100 * GeV;
71
72 verboseLevel=0 ;
73 // Verbosity scale:
74 // 0 = nothing
75 // 1 = warning for energy non-conservation
76 // 2 = details of energy budget
77 // 3 = calculation of cross sections, file openings, sampling of atoms
78 // 4 = entering in methods
79
80 if( verboseLevel>0 ) {
81 G4cout << "Livermore Compton model is constructed " << G4endl
82 << "Energy range: "
83 << lowEnergyLimit / eV << " eV - "
84 << highEnergyLimit / GeV << " GeV"
85 << G4endl;
86 }
87
88 //Mark this model as "applicable" for atomic deexcitation
90
91}
92
93//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
94
96{
97 delete crossSectionHandler;
98 delete scatterFunctionData;
99}
100
101//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
102
104 const G4DataVector& cuts)
105{
106 if (verboseLevel > 2) {
107 G4cout << "Calling G4LivermoreComptonModel::Initialise()" << G4endl;
108 }
109
110 if (crossSectionHandler)
111 {
112 crossSectionHandler->Clear();
113 delete crossSectionHandler;
114 }
115 delete scatterFunctionData;
116
117 // Reading of data files - all materials are read
118 crossSectionHandler = new G4CrossSectionHandler;
119 G4String crossSectionFile = "comp/ce-cs-";
120 crossSectionHandler->LoadData(crossSectionFile);
121
122 G4VDataSetAlgorithm* scatterInterpolation = new G4LogLogInterpolation;
123 G4String scatterFile = "comp/ce-sf-";
124 scatterFunctionData = new G4CompositeEMDataSet(scatterInterpolation, 1., 1.);
125 scatterFunctionData->LoadData(scatterFile);
126
127 // For Doppler broadening
128 shellData.SetOccupancyData();
129 G4String file = "/doppler/shell-doppler";
130 shellData.LoadData(file);
131
132 InitialiseElementSelectors(particle,cuts);
133
134 if (verboseLevel > 2) {
135 G4cout << "Loaded cross section files for Livermore Compton model" << G4endl;
136 }
137
138 if(isInitialised) { return; }
139 isInitialised = true;
140
142
143 fAtomDeexcitation = G4LossTableManager::Instance()->AtomDeexcitation();
144
145 if( verboseLevel>0 ) {
146 G4cout << "Livermore Compton model is initialized " << G4endl
147 << "Energy range: "
148 << LowEnergyLimit() / eV << " eV - "
149 << HighEnergyLimit() / GeV << " GeV"
150 << G4endl;
151 }
152}
153
154//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
155
158 G4double GammaEnergy,
161{
162 if (verboseLevel > 3) {
163 G4cout << "Calling ComputeCrossSectionPerAtom() of G4LivermoreComptonModel" << G4endl;
164 }
165 if (GammaEnergy < lowEnergyLimit || GammaEnergy > highEnergyLimit) { return 0.0; }
166
167 G4double cs = crossSectionHandler->FindValue(G4int(Z), GammaEnergy);
168 return cs;
169}
170
171//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
172
173void G4LivermoreComptonModel::SampleSecondaries(std::vector<G4DynamicParticle*>* fvect,
174 const G4MaterialCutsCouple* couple,
175 const G4DynamicParticle* aDynamicGamma,
177{
178
179 // The scattered gamma energy is sampled according to Klein - Nishina formula.
180 // then accepted or rejected depending on the Scattering Function multiplied
181 // by factor from Klein - Nishina formula.
182 // Expression of the angular distribution as Klein Nishina
183 // angular and energy distribution and Scattering fuctions is taken from
184 // D. E. Cullen "A simple model of photon transport" Nucl. Instr. Meth.
185 // Phys. Res. B 101 (1995). Method of sampling with form factors is different
186 // data are interpolated while in the article they are fitted.
187 // Reference to the article is from J. Stepanek New Photon, Positron
188 // and Electron Interaction Data for GEANT in Energy Range from 1 eV to 10
189 // TeV (draft).
190 // The random number techniques of Butcher & Messel are used
191 // (Nucl Phys 20(1960),15).
192
193 G4double photonEnergy0 = aDynamicGamma->GetKineticEnergy();
194
195 if (verboseLevel > 3) {
196 G4cout << "G4LivermoreComptonModel::SampleSecondaries() E(MeV)= "
197 << photonEnergy0/MeV << " in " << couple->GetMaterial()->GetName()
198 << G4endl;
199 }
200
201 // low-energy gamma is absorpted by this process
202 if (photonEnergy0 <= lowEnergyLimit)
203 {
207 return ;
208 }
209
210 G4double e0m = photonEnergy0 / electron_mass_c2 ;
211 G4ParticleMomentum photonDirection0 = aDynamicGamma->GetMomentumDirection();
212
213 // Select randomly one element in the current material
214 const G4ParticleDefinition* particle = aDynamicGamma->GetDefinition();
215 const G4Element* elm = SelectRandomAtom(couple,particle,photonEnergy0);
216 G4int Z = (G4int)elm->GetZ();
217
218 G4double epsilon0Local = 1. / (1. + 2. * e0m);
219 G4double epsilon0Sq = epsilon0Local * epsilon0Local;
220 G4double alpha1 = -std::log(epsilon0Local);
221 G4double alpha2 = 0.5 * (1. - epsilon0Sq);
222
223 G4double wlPhoton = h_Planck*c_light/photonEnergy0;
224
225 // Sample the energy of the scattered photon
226 G4double epsilon;
227 G4double epsilonSq;
228 G4double oneCosT;
229 G4double sinT2;
230 G4double gReject;
231
232 do
233 {
234 if ( alpha1/(alpha1+alpha2) > G4UniformRand())
235 {
236 // std::pow(epsilon0Local,G4UniformRand())
237 epsilon = std::exp(-alpha1 * G4UniformRand());
238 epsilonSq = epsilon * epsilon;
239 }
240 else
241 {
242 epsilonSq = epsilon0Sq + (1. - epsilon0Sq) * G4UniformRand();
243 epsilon = std::sqrt(epsilonSq);
244 }
245
246 oneCosT = (1. - epsilon) / ( epsilon * e0m);
247 sinT2 = oneCosT * (2. - oneCosT);
248 G4double x = std::sqrt(oneCosT/2.) / (wlPhoton/cm);
249 G4double scatteringFunction = scatterFunctionData->FindValue(x,Z-1);
250 gReject = (1. - epsilon * sinT2 / (1. + epsilonSq)) * scatteringFunction;
251
252 } while(gReject < G4UniformRand()*Z);
253
254 G4double cosTheta = 1. - oneCosT;
255 G4double sinTheta = std::sqrt (sinT2);
256 G4double phi = twopi * G4UniformRand() ;
257 G4double dirx = sinTheta * std::cos(phi);
258 G4double diry = sinTheta * std::sin(phi);
259 G4double dirz = cosTheta ;
260
261 // Doppler broadening - Method based on:
262 // Y. Namito, S. Ban and H. Hirayama,
263 // "Implementation of the Doppler Broadening of a Compton-Scattered Photon
264 // into the EGS4 Code", NIM A 349, pp. 489-494, 1994
265
266 // Maximum number of sampling iterations
267 G4int maxDopplerIterations = 1000;
268 G4double bindingE = 0.;
269 G4double photonEoriginal = epsilon * photonEnergy0;
270 G4double photonE = -1.;
271 G4int iteration = 0;
272 G4double eMax = photonEnergy0;
273
274 G4int shellIdx = 0;
275
276 do
277 {
278 ++iteration;
279 // Select shell based on shell occupancy
280 shellIdx = shellData.SelectRandomShell(Z);
281 bindingE = shellData.BindingEnergy(Z,shellIdx);
282
283 eMax = photonEnergy0 - bindingE;
284
285 // Randomly sample bound electron momentum
286 // (memento: the data set is in Atomic Units)
287 G4double pSample = profileData.RandomSelectMomentum(Z,shellIdx);
288 // Rescale from atomic units
289 G4double pDoppler = pSample * fine_structure_const;
290 G4double pDoppler2 = pDoppler * pDoppler;
291 G4double var2 = 1. + oneCosT * e0m;
292 G4double var3 = var2*var2 - pDoppler2;
293 G4double var4 = var2 - pDoppler2 * cosTheta;
294 G4double var = var4*var4 - var3 + pDoppler2 * var3;
295 if (var > 0.)
296 {
297 G4double varSqrt = std::sqrt(var);
298 G4double scale = photonEnergy0 / var3;
299 // Random select either root
300 if (G4UniformRand() < 0.5) { photonE = (var4 - varSqrt) * scale; }
301 else { photonE = (var4 + varSqrt) * scale; }
302 }
303 else
304 {
305 photonE = -1.;
306 }
307 } while ( iteration <= maxDopplerIterations &&
308
309// JB : corrected the following condition
310// (photonE < 0. || photonE > eMax || photonE < eMax*G4UniformRand()) );
311
312 (photonE > eMax ) );
313
314 // End of recalculation of photon energy with Doppler broadening
315 // Revert to original if maximum number of iterations threshold has been reached
316
317 if (iteration >= maxDopplerIterations)
318 {
319 photonE = photonEoriginal;
320 bindingE = 0.;
321 }
322
323 // Update G4VParticleChange for the scattered photon
324
325 G4ThreeVector photonDirection1(dirx,diry,dirz);
326 photonDirection1.rotateUz(photonDirection0);
327 fParticleChange->ProposeMomentumDirection(photonDirection1) ;
328
329 G4double photonEnergy1 = photonE;
330
331 if (photonEnergy1 > 0.)
332 {
334 }
335 else
336 {
337 photonEnergy1 = 0.;
340 }
341
342 // Kinematics of the scattered electron
343 G4double eKineticEnergy = photonEnergy0 - photonEnergy1 - bindingE;
344
345 // protection against negative final energy: no e- is created
346 if(eKineticEnergy < 0.0) {
347 bindingE = photonEnergy0 - photonEnergy1;
348
349 } else {
350 G4double eTotalEnergy = eKineticEnergy + electron_mass_c2;
351
352 G4double electronE = photonEnergy0 * (1. - epsilon) + electron_mass_c2;
353 G4double electronP2 = electronE*electronE - electron_mass_c2*electron_mass_c2;
354 G4double sinThetaE = -1.;
355 G4double cosThetaE = 0.;
356 if (electronP2 > 0.)
357 {
358 cosThetaE = (eTotalEnergy + photonEnergy1 )* (1. - epsilon) / std::sqrt(electronP2);
359 sinThetaE = -sqrt((1. - cosThetaE) * (1. + cosThetaE));
360 }
361
362 G4double eDirX = sinThetaE * std::cos(phi);
363 G4double eDirY = sinThetaE * std::sin(phi);
364 G4double eDirZ = cosThetaE;
365
366 G4ThreeVector eDirection(eDirX,eDirY,eDirZ);
367 eDirection.rotateUz(photonDirection0);
369 eDirection,eKineticEnergy) ;
370 fvect->push_back(dp);
371 }
372
373 // sample deexcitation
374 //
375 if(fAtomDeexcitation && iteration < maxDopplerIterations) {
376 G4int index = couple->GetIndex();
377 if(fAtomDeexcitation->CheckDeexcitationActiveRegion(index)) {
378 size_t nbefore = fvect->size();
380 const G4AtomicShell* shell = fAtomDeexcitation->GetAtomicShell(Z, as);
381 fAtomDeexcitation->GenerateParticles(fvect, shell, Z, index);
382 size_t nafter = fvect->size();
383 if(nafter > nbefore) {
384 for (size_t i=nbefore; i<nafter; ++i) {
385 bindingE -= ((*fvect)[i])->GetKineticEnergy();
386 }
387 }
388 }
389 }
390 if(bindingE < 0.0) { bindingE = 0.0; }
392}
393
G4AtomicShellEnumerator
@ fStopAndKill
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:53
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:72
G4double RandomSelectMomentum(G4int Z, G4int shellIndex) const
const G4ThreeVector & GetMomentumDirection() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
static G4Electron * Electron()
Definition: G4Electron.cc:94
G4double GetZ() const
Definition: G4Element.hh:131
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
G4LivermoreComptonModel(const G4ParticleDefinition *p=0, const G4String &nam="LivermoreCompton")
G4ParticleChangeForGamma * fParticleChange
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0, G4double cut=0, G4double emax=DBL_MAX)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
static G4LossTableManager * Instance()
G4VAtomDeexcitation * AtomDeexcitation()
const G4Material * GetMaterial() const
const G4String & GetName() const
Definition: G4Material.hh:177
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
void SetOccupancyData()
Definition: G4ShellData.hh:70
G4double BindingEnergy(G4int Z, G4int shellIndex) const
Definition: G4ShellData.cc:166
void LoadData(const G4String &fileName)
Definition: G4ShellData.cc:234
G4int SelectRandomShell(G4int Z) const
Definition: G4ShellData.cc:363
G4bool CheckDeexcitationActiveRegion(G4int coupleIndex)
virtual const G4AtomicShell * GetAtomicShell(G4int Z, G4AtomicShellEnumerator shell)=0
void GenerateParticles(std::vector< G4DynamicParticle * > *secVect, const G4AtomicShell *, G4int Z, G4int coupleIndex)
G4double FindValue(G4int Z, G4double e) const
void LoadData(const G4String &dataFile)
virtual G4double FindValue(G4double x, G4int componentId=0) const =0
virtual G4bool LoadData(const G4String &fileName)=0
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 SetDeexcitationFlag(G4bool val)
Definition: G4VEmModel.hh:641
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:123
void ProposeTrackStatus(G4TrackStatus status)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)