source: trunk/source/processes/electromagnetic/lowenergy/src/G4LivermoreComptonModel.cc@ 992

Last change on this file since 992 was 991, checked in by garnier, 17 years ago

update

File size: 12.8 KB
Line 
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: G4LivermoreComptonModel.cc,v 1.1 2008/10/30 14:17:46 sincerti Exp $
27// GEANT4 tag $Name: geant4-09-02 $
28//
29
30#include "G4LivermoreComptonModel.hh"
31
32//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
33
34using namespace std;
35
36//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
37
38G4LivermoreComptonModel::G4LivermoreComptonModel(const G4ParticleDefinition*,
39 const G4String& nam)
40:G4VEmModel(nam),isInitialised(false)
41{
42 lowEnergyLimit = 250 * eV; // SI - Could be 10 eV ?
43 highEnergyLimit = 100 * GeV;
44 SetLowEnergyLimit(lowEnergyLimit);
45 SetHighEnergyLimit(highEnergyLimit);
46
47 verboseLevel= 0;
48 // Verbosity scale:
49 // 0 = nothing
50 // 1 = warning for energy non-conservation
51 // 2 = details of energy budget
52 // 3 = calculation of cross sections, file openings, sampling of atoms
53 // 4 = entering in methods
54
55 G4cout << "Livermore Compton model is constructed " << G4endl
56 << "Energy range: "
57 << lowEnergyLimit / keV << " keV - "
58 << highEnergyLimit / GeV << " GeV"
59 << G4endl;
60
61}
62
63//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
64
65G4LivermoreComptonModel::~G4LivermoreComptonModel()
66{
67 delete meanFreePathTable;
68 delete crossSectionHandler;
69 delete scatterFunctionData;
70}
71
72//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
73
74void G4LivermoreComptonModel::Initialise(const G4ParticleDefinition* particle,
75 const G4DataVector& cuts)
76{
77 if (verboseLevel > 3)
78 G4cout << "Calling G4LivermoreComptonModel::Initialise()" << G4endl;
79
80 InitialiseElementSelectors(particle,cuts);
81
82 // Energy limits
83
84 if (LowEnergyLimit() < lowEnergyLimit)
85 {
86 G4cout << "G4LivermoreComptonModel: low energy limit increased from " <<
87 LowEnergyLimit()/eV << " eV to " << lowEnergyLimit << " eV" << G4endl;
88 SetLowEnergyLimit(lowEnergyLimit);
89 }
90
91 if (HighEnergyLimit() > highEnergyLimit)
92 {
93 G4cout << "G4LivermoreComptonModel: high energy limit decreased from " <<
94 HighEnergyLimit()/GeV << " GeV to " << highEnergyLimit << " GeV" << G4endl;
95 SetHighEnergyLimit(highEnergyLimit);
96 }
97
98 // Reading of data files - all materials are read
99
100 crossSectionHandler = new G4CrossSectionHandler;
101 crossSectionHandler->Clear();
102 G4String crossSectionFile = "comp/ce-cs-";
103 crossSectionHandler->LoadData(crossSectionFile);
104
105 G4VDataSetAlgorithm* scatterInterpolation = new G4LogLogInterpolation;
106 G4String scatterFile = "comp/ce-sf-";
107 scatterFunctionData = new G4CompositeEMDataSet(scatterInterpolation, 1., 1.);
108 scatterFunctionData->LoadData(scatterFile);
109
110 // For Doppler broadening
111 shellData.SetOccupancyData();
112 G4String file = "/doppler/shell-doppler";
113 shellData.LoadData(file);
114
115 meanFreePathTable = 0;
116 meanFreePathTable = crossSectionHandler->BuildMeanFreePathForMaterials();
117
118 if (verboseLevel > 2)
119 G4cout << "Loaded cross section files for Livermore Compton model" << G4endl;
120
121 G4cout << "Livermore Compton model is initialized " << G4endl
122 << "Energy range: "
123 << LowEnergyLimit() / keV << " keV - "
124 << HighEnergyLimit() / GeV << " GeV"
125 << G4endl;
126
127 //
128
129 if(isInitialised) return;
130
131 if(pParticleChange)
132 fParticleChange = reinterpret_cast<G4ParticleChangeForGamma*>(pParticleChange);
133 else
134 fParticleChange = new G4ParticleChangeForGamma();
135
136 isInitialised = true;
137}
138
139//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
140
141G4double G4LivermoreComptonModel::ComputeCrossSectionPerAtom(
142 const G4ParticleDefinition*,
143 G4double GammaEnergy,
144 G4double Z, G4double,
145 G4double, G4double)
146{
147 if (verboseLevel > 3)
148 G4cout << "Calling ComputeCrossSectionPerAtom() of G4LivermoreComptonModel" << G4endl;
149
150 G4double cs = crossSectionHandler->FindValue(G4int(Z), GammaEnergy);
151 return cs;
152}
153
154//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
155
156void G4LivermoreComptonModel::SampleSecondaries(std::vector<G4DynamicParticle*>* fvect,
157 const G4MaterialCutsCouple* couple,
158 const G4DynamicParticle* aDynamicGamma,
159 G4double,
160 G4double)
161{
162 // The scattered gamma energy is sampled according to Klein - Nishina formula.
163 // then accepted or rejected depending on the Scattering Function multiplied
164 // by factor from Klein - Nishina formula.
165 // Expression of the angular distribution as Klein Nishina
166 // angular and energy distribution and Scattering fuctions is taken from
167 // D. E. Cullen "A simple model of photon transport" Nucl. Instr. Meth.
168 // Phys. Res. B 101 (1995). Method of sampling with form factors is different
169 // data are interpolated while in the article they are fitted.
170 // Reference to the article is from J. Stepanek New Photon, Positron
171 // and Electron Interaction Data for GEANT in Energy Range from 1 eV to 10
172 // TeV (draft).
173 // The random number techniques of Butcher & Messel are used
174 // (Nucl Phys 20(1960),15).
175
176 if (verboseLevel > 3)
177 G4cout << "Calling SampleSecondaries() of G4LivermoreComptonModel" << G4endl;
178
179 G4double photonEnergy0 = aDynamicGamma->GetKineticEnergy();
180
181 if (photonEnergy0 <= lowEnergyLimit)
182 {
183 fParticleChange->ProposeTrackStatus(fStopAndKill);
184 fParticleChange->SetProposedKineticEnergy(0.);
185 fParticleChange->ProposeLocalEnergyDeposit(photonEnergy0);
186 // SI - IS THE FOLLOWING RETURN NECESSARY ?
187 return ;
188 }
189
190 G4double e0m = photonEnergy0 / electron_mass_c2 ;
191 G4ParticleMomentum photonDirection0 = aDynamicGamma->GetMomentumDirection();
192
193 // Select randomly one element in the current material
194 G4int Z = crossSectionHandler->SelectRandomAtom(couple,photonEnergy0);
195
196 G4double epsilon0 = 1. / (1. + 2. * e0m);
197 G4double epsilon0Sq = epsilon0 * epsilon0;
198 G4double alpha1 = -std::log(epsilon0);
199 G4double alpha2 = 0.5 * (1. - epsilon0Sq);
200
201 G4double wlPhoton = h_Planck*c_light/photonEnergy0;
202
203 // Sample the energy of the scattered photon
204 G4double epsilon;
205 G4double epsilonSq;
206 G4double oneCosT;
207 G4double sinT2;
208 G4double gReject;
209
210 do
211 {
212 if ( alpha1/(alpha1+alpha2) > G4UniformRand())
213 {
214 epsilon = std::exp(-alpha1 * G4UniformRand()); // std::pow(epsilon0,G4UniformRand())
215 epsilonSq = epsilon * epsilon;
216 }
217 else
218 {
219 epsilonSq = epsilon0Sq + (1. - epsilon0Sq) * G4UniformRand();
220 epsilon = std::sqrt(epsilonSq);
221 }
222
223 oneCosT = (1. - epsilon) / ( epsilon * e0m);
224 sinT2 = oneCosT * (2. - oneCosT);
225 G4double x = std::sqrt(oneCosT/2.) / (wlPhoton/cm);
226 G4double scatteringFunction = scatterFunctionData->FindValue(x,Z-1);
227 gReject = (1. - epsilon * sinT2 / (1. + epsilonSq)) * scatteringFunction;
228
229 } while(gReject < G4UniformRand()*Z);
230
231 G4double cosTheta = 1. - oneCosT;
232 G4double sinTheta = std::sqrt (sinT2);
233 G4double phi = twopi * G4UniformRand() ;
234 G4double dirx = sinTheta * std::cos(phi);
235 G4double diry = sinTheta * std::sin(phi);
236 G4double dirz = cosTheta ;
237
238 // Doppler broadening - Method based on:
239 // Y. Namito, S. Ban and H. Hirayama,
240 // "Implementation of the Doppler Broadening of a Compton-Scattered Photon Into the EGS4 Code"
241 // NIM A 349, pp. 489-494, 1994
242
243 // Maximum number of sampling iterations
244 G4int maxDopplerIterations = 1000;
245 G4double bindingE = 0.;
246 G4double photonEoriginal = epsilon * photonEnergy0;
247 G4double photonE = -1.;
248 G4int iteration = 0;
249 G4double eMax = photonEnergy0;
250 do
251 {
252 iteration++;
253 // Select shell based on shell occupancy
254 G4int shell = shellData.SelectRandomShell(Z);
255 bindingE = shellData.BindingEnergy(Z,shell);
256
257 eMax = photonEnergy0 - bindingE;
258
259 // Randomly sample bound electron momentum (memento: the data set is in Atomic Units)
260 G4double pSample = profileData.RandomSelectMomentum(Z,shell);
261 // Rescale from atomic units
262 G4double pDoppler = pSample * fine_structure_const;
263 G4double pDoppler2 = pDoppler * pDoppler;
264 G4double var2 = 1. + oneCosT * e0m;
265 G4double var3 = var2*var2 - pDoppler2;
266 G4double var4 = var2 - pDoppler2 * cosTheta;
267 G4double var = var4*var4 - var3 + pDoppler2 * var3;
268 if (var > 0.)
269 {
270 G4double varSqrt = std::sqrt(var);
271 G4double scale = photonEnergy0 / var3;
272 // Random select either root
273 if (G4UniformRand() < 0.5) photonE = (var4 - varSqrt) * scale;
274 else photonE = (var4 + varSqrt) * scale;
275 }
276 else
277 {
278 photonE = -1.;
279 }
280 } while ( iteration <= maxDopplerIterations &&
281 (photonE < 0. || photonE > eMax || photonE < eMax*G4UniformRand()) );
282
283 // End of recalculation of photon energy with Doppler broadening
284 // Revert to original if maximum number of iterations threshold has been reached
285
286 if (iteration >= maxDopplerIterations)
287 {
288 photonE = photonEoriginal;
289 bindingE = 0.;
290 }
291
292 // Update G4VParticleChange for the scattered photon
293
294 G4ThreeVector photonDirection1(dirx,diry,dirz);
295 photonDirection1.rotateUz(photonDirection0);
296 fParticleChange->ProposeMomentumDirection(photonDirection1) ;
297
298 G4double photonEnergy1 = photonE;
299
300 if (photonEnergy1 > 0.)
301 {
302 fParticleChange->SetProposedKineticEnergy(photonEnergy1) ;
303 }
304 else
305 {
306 fParticleChange->SetProposedKineticEnergy(0.) ;
307 fParticleChange->ProposeTrackStatus(fStopAndKill);
308 }
309
310 // Kinematics of the scattered electron
311 G4double eKineticEnergy = photonEnergy0 - photonEnergy1 - bindingE;
312 G4double eTotalEnergy = eKineticEnergy + electron_mass_c2;
313
314 G4double electronE = photonEnergy0 * (1. - epsilon) + electron_mass_c2;
315 G4double electronP2 = electronE*electronE - electron_mass_c2*electron_mass_c2;
316 G4double sinThetaE = -1.;
317 G4double cosThetaE = 0.;
318 if (electronP2 > 0.)
319 {
320 cosThetaE = (eTotalEnergy + photonEnergy1 )* (1. - epsilon) / std::sqrt(electronP2);
321 sinThetaE = -1. * sqrt(1. - cosThetaE * cosThetaE);
322 }
323
324 G4double eDirX = sinThetaE * std::cos(phi);
325 G4double eDirY = sinThetaE * std::sin(phi);
326 G4double eDirZ = cosThetaE;
327
328 G4ThreeVector eDirection(eDirX,eDirY,eDirZ);
329 eDirection.rotateUz(photonDirection0);
330
331// SI - The range test has been removed wrt original G4LowEnergyCompton class
332
333 fParticleChange->ProposeLocalEnergyDeposit(bindingE);
334
335 G4DynamicParticle* dp = new G4DynamicParticle (G4Electron::Electron(),eDirection,eKineticEnergy) ;
336 fvect->push_back(dp);
337}
338
339//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
340
341G4double G4LivermoreComptonModel::GetMeanFreePath(const G4Track& track,
342 G4double, // previousStepSize
343 G4ForceCondition*)
344{
345 const G4DynamicParticle* photon = track.GetDynamicParticle();
346 G4double energy = photon->GetKineticEnergy();
347 const G4MaterialCutsCouple* couple = track.GetMaterialCutsCouple();
348 size_t materialIndex = couple->GetIndex();
349
350 G4double meanFreePath;
351 if (energy > highEnergyLimit) meanFreePath = meanFreePathTable->FindValue(highEnergyLimit,materialIndex);
352 else if (energy < lowEnergyLimit) meanFreePath = DBL_MAX;
353 else meanFreePath = meanFreePathTable->FindValue(energy,materialIndex);
354 return meanFreePath;
355}
356
Note: See TracBrowser for help on using the repository browser.