source: trunk/source/processes/electromagnetic/lowenergy/src/G4LowEnergyGammaConversion.cc@ 1014

Last change on this file since 1014 was 1007, checked in by garnier, 17 years ago

update to geant4.9.2

File size: 12.6 KB
RevLine 
[819]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// --------------------------------------------------------------------
27///
28// $Id: G4LowEnergyGammaConversion.cc,v 1.36 2006/06/29 19:40:17 gunter Exp $
[1007]29// GEANT4 tag $Name: geant4-09-02 $
[819]30//
31//
32// --------------------------------------------------------------
33//
34// Author: A. Forti
35// Maria Grazia Pia (Maria.Grazia.Pia@cern.ch)
36//
37// History:
38// --------
39// 02/03/1999 A. Forti 1st implementation
40// 14.03.2000 Veronique Lefebure;
41// Change initialisation of lowestEnergyLimit from 1.22 to 1.022.
42// Note that the hard coded value 1.022 should be used instead of
43// 2*electron_mass_c2 in order to agree with the value of the data bank EPDL97
44// 24.04.01 V.Ivanchenko remove RogueWave
45// 27.07.01 F.Longo correct bug in energy distribution
46// 21.01.03 V.Ivanchenko Cut per region
47// 25.03.03 F.Longo fix in angular distribution of e+/e-
48// 24.04.03 V.Ivanchenko - Cut per region mfpt
49//
50// --------------------------------------------------------------
51
52#include "G4LowEnergyGammaConversion.hh"
53
54#include "Randomize.hh"
55#include "G4ParticleDefinition.hh"
56#include "G4Track.hh"
57#include "G4Step.hh"
58#include "G4ForceCondition.hh"
59#include "G4Gamma.hh"
60#include "G4Electron.hh"
61#include "G4DynamicParticle.hh"
62#include "G4VParticleChange.hh"
63#include "G4ThreeVector.hh"
64#include "G4Positron.hh"
65#include "G4IonisParamElm.hh"
66#include "G4Material.hh"
67#include "G4VCrossSectionHandler.hh"
68#include "G4CrossSectionHandler.hh"
69#include "G4VEMDataSet.hh"
70#include "G4VDataSetAlgorithm.hh"
71#include "G4LogLogInterpolation.hh"
72#include "G4VRangeTest.hh"
73#include "G4RangeTest.hh"
74#include "G4MaterialCutsCouple.hh"
75
76G4LowEnergyGammaConversion::G4LowEnergyGammaConversion(const G4String& processName)
77 : G4VDiscreteProcess(processName),
78 lowEnergyLimit(1.022000*MeV),
79 highEnergyLimit(100*GeV),
80 intrinsicLowEnergyLimit(1.022000*MeV),
81 intrinsicHighEnergyLimit(100*GeV),
82 smallEnergy(2.*MeV)
83
84{
85 if (lowEnergyLimit < intrinsicLowEnergyLimit ||
86 highEnergyLimit > intrinsicHighEnergyLimit)
87 {
88 G4Exception("G4LowEnergyGammaConversion::G4LowEnergyGammaConversion - energy limit outside intrinsic process validity range");
89 }
90
91 // The following pointer is owned by G4DataHandler
92
93 crossSectionHandler = new G4CrossSectionHandler();
94 crossSectionHandler->Initialise(0,1.0220*MeV,100.*GeV,400);
95 meanFreePathTable = 0;
96 rangeTest = new G4RangeTest;
97
98 if (verboseLevel > 0)
99 {
100 G4cout << GetProcessName() << " is created " << G4endl
101 << "Energy range: "
102 << lowEnergyLimit / MeV << " MeV - "
103 << highEnergyLimit / GeV << " GeV"
104 << G4endl;
105 }
106}
107
108G4LowEnergyGammaConversion::~G4LowEnergyGammaConversion()
109{
110 delete meanFreePathTable;
111 delete crossSectionHandler;
112 delete rangeTest;
113}
114
115void G4LowEnergyGammaConversion::BuildPhysicsTable(const G4ParticleDefinition& )
116{
117
118 crossSectionHandler->Clear();
119 G4String crossSectionFile = "pair/pp-cs-";
120 crossSectionHandler->LoadData(crossSectionFile);
121
122 delete meanFreePathTable;
123 meanFreePathTable = crossSectionHandler->BuildMeanFreePathForMaterials();
124}
125
126G4VParticleChange* G4LowEnergyGammaConversion::PostStepDoIt(const G4Track& aTrack,
127 const G4Step& aStep)
128{
129// The energies of the e+ e- secondaries are sampled using the Bethe - Heitler
130// cross sections with Coulomb correction. A modified version of the random
131// number techniques of Butcher & Messel is used (Nuc Phys 20(1960),15).
132
133// Note 1 : Effects due to the breakdown of the Born approximation at low
134// energy are ignored.
135// Note 2 : The differential cross section implicitly takes account of
136// pair creation in both nuclear and atomic electron fields. However triplet
137// prodution is not generated.
138
139 aParticleChange.Initialize(aTrack);
140
141 const G4MaterialCutsCouple* couple = aTrack.GetMaterialCutsCouple();
142
143 const G4DynamicParticle* incidentPhoton = aTrack.GetDynamicParticle();
144 G4double photonEnergy = incidentPhoton->GetKineticEnergy();
145 G4ParticleMomentum photonDirection = incidentPhoton->GetMomentumDirection();
146
147 G4double epsilon ;
148 G4double epsilon0 = electron_mass_c2 / photonEnergy ;
149
150 // Do it fast if photon energy < 2. MeV
151 if (photonEnergy < smallEnergy )
152 {
153 epsilon = epsilon0 + (0.5 - epsilon0) * G4UniformRand();
154 }
155 else
156 {
157 // Select randomly one element in the current material
158 const G4Element* element = crossSectionHandler->SelectRandomElement(couple,photonEnergy);
159
160 if (element == 0)
161 {
162 G4cout << "G4LowEnergyGammaConversion::PostStepDoIt - element = 0" << G4endl;
163 }
164 G4IonisParamElm* ionisation = element->GetIonisation();
165 if (ionisation == 0)
166 {
167 G4cout << "G4LowEnergyGammaConversion::PostStepDoIt - ionisation = 0" << G4endl;
168 }
169
170 // Extract Coulomb factor for this Element
171 G4double fZ = 8. * (ionisation->GetlogZ3());
172 if (photonEnergy > 50. * MeV) fZ += 8. * (element->GetfCoulomb());
173
174 // Limits of the screening variable
175 G4double screenFactor = 136. * epsilon0 / (element->GetIonisation()->GetZ3()) ;
176 G4double screenMax = std::exp ((42.24 - fZ)/8.368) - 0.952 ;
177 G4double screenMin = std::min(4.*screenFactor,screenMax) ;
178
179 // Limits of the energy sampling
180 G4double epsilon1 = 0.5 - 0.5 * std::sqrt(1. - screenMin / screenMax) ;
181 G4double epsilonMin = std::max(epsilon0,epsilon1);
182 G4double epsilonRange = 0.5 - epsilonMin ;
183
184 // Sample the energy rate of the created electron (or positron)
185 G4double screen;
186 G4double gReject ;
187
188 G4double f10 = ScreenFunction1(screenMin) - fZ;
189 G4double f20 = ScreenFunction2(screenMin) - fZ;
190 G4double normF1 = std::max(f10 * epsilonRange * epsilonRange,0.);
191 G4double normF2 = std::max(1.5 * f20,0.);
192
193 do {
194 if (normF1 / (normF1 + normF2) > G4UniformRand() )
195 {
196 epsilon = 0.5 - epsilonRange * std::pow(G4UniformRand(), 0.3333) ;
197 screen = screenFactor / (epsilon * (1. - epsilon));
198 gReject = (ScreenFunction1(screen) - fZ) / f10 ;
199 }
200 else
201 {
202 epsilon = epsilonMin + epsilonRange * G4UniformRand();
203 screen = screenFactor / (epsilon * (1 - epsilon));
204 gReject = (ScreenFunction2(screen) - fZ) / f20 ;
205 }
206 } while ( gReject < G4UniformRand() );
207
208 } // End of epsilon sampling
209
210 // Fix charges randomly
211
212 G4double electronTotEnergy;
213 G4double positronTotEnergy;
214
215 if (CLHEP::RandBit::shootBit())
216 {
217 electronTotEnergy = (1. - epsilon) * photonEnergy;
218 positronTotEnergy = epsilon * photonEnergy;
219 }
220 else
221 {
222 positronTotEnergy = (1. - epsilon) * photonEnergy;
223 electronTotEnergy = epsilon * photonEnergy;
224 }
225
226 // Scattered electron (positron) angles. ( Z - axis along the parent photon)
227 // Universal distribution suggested by L. Urban (Geant3 manual (1993) Phys211),
228 // derived from Tsai distribution (Rev. Mod. Phys. 49, 421 (1977)
229
230 G4double u;
231 const G4double a1 = 0.625;
232 G4double a2 = 3. * a1;
233 // G4double d = 27. ;
234
235 // if (9. / (9. + d) > G4UniformRand())
236 if (0.25 > G4UniformRand())
237 {
238 u = - std::log(G4UniformRand() * G4UniformRand()) / a1 ;
239 }
240 else
241 {
242 u = - std::log(G4UniformRand() * G4UniformRand()) / a2 ;
243 }
244
245 G4double thetaEle = u*electron_mass_c2/electronTotEnergy;
246 G4double thetaPos = u*electron_mass_c2/positronTotEnergy;
247 G4double phi = twopi * G4UniformRand();
248
249 G4double dxEle= std::sin(thetaEle)*std::cos(phi),dyEle= std::sin(thetaEle)*std::sin(phi),dzEle=std::cos(thetaEle);
250 G4double dxPos=-std::sin(thetaPos)*std::cos(phi),dyPos=-std::sin(thetaPos)*std::sin(phi),dzPos=std::cos(thetaPos);
251
252
253 // Kinematics of the created pair:
254 // the electron and positron are assumed to have a symetric angular
255 // distribution with respect to the Z axis along the parent photon
256
257 G4double localEnergyDeposit = 0. ;
258
259 aParticleChange.SetNumberOfSecondaries(2) ;
260 G4double electronKineEnergy = std::max(0.,electronTotEnergy - electron_mass_c2) ;
261
262 // Generate the electron only if with large enough range w.r.t. cuts and safety
263
264 G4double safety = aStep.GetPostStepPoint()->GetSafety();
265
266 if (rangeTest->Escape(G4Electron::Electron(),couple,electronKineEnergy,safety))
267 {
268 G4ThreeVector electronDirection (dxEle, dyEle, dzEle);
269 electronDirection.rotateUz(photonDirection);
270
271 G4DynamicParticle* particle1 = new G4DynamicParticle (G4Electron::Electron(),
272 electronDirection,
273 electronKineEnergy);
274 aParticleChange.AddSecondary(particle1) ;
275 }
276 else
277 {
278 localEnergyDeposit += electronKineEnergy ;
279 }
280
281 // The e+ is always created (even with kinetic energy = 0) for further annihilation
282 G4double positronKineEnergy = std::max(0.,positronTotEnergy - electron_mass_c2) ;
283
284 // Is the local energy deposit correct, if the positron is always created?
285 if (! (rangeTest->Escape(G4Positron::Positron(),couple,positronKineEnergy,safety)))
286 {
287 localEnergyDeposit += positronKineEnergy ;
288 positronKineEnergy = 0. ;
289 }
290
291 G4ThreeVector positronDirection (dxPos, dyPos, dzPos);
292 positronDirection.rotateUz(photonDirection);
293
294 // Create G4DynamicParticle object for the particle2
295 G4DynamicParticle* particle2 = new G4DynamicParticle(G4Positron::Positron(),
296 positronDirection, positronKineEnergy);
297 aParticleChange.AddSecondary(particle2) ;
298
299 aParticleChange.ProposeLocalEnergyDeposit(localEnergyDeposit) ;
300
301 // Kill the incident photon
302 aParticleChange.ProposeMomentumDirection(0.,0.,0.) ;
303 aParticleChange.ProposeEnergy(0.) ;
304 aParticleChange.ProposeTrackStatus(fStopAndKill) ;
305
306 // Reset NbOfInteractionLengthLeft and return aParticleChange
307 return G4VDiscreteProcess::PostStepDoIt(aTrack,aStep);
308}
309
310G4bool G4LowEnergyGammaConversion::IsApplicable(const G4ParticleDefinition& particle)
311{
312 return ( &particle == G4Gamma::Gamma() );
313}
314
315G4double G4LowEnergyGammaConversion::GetMeanFreePath(const G4Track& track,
316 G4double, // previousStepSize
317 G4ForceCondition*)
318{
319 const G4DynamicParticle* photon = track.GetDynamicParticle();
320 G4double energy = photon->GetKineticEnergy();
321 const G4MaterialCutsCouple* couple = track.GetMaterialCutsCouple();
322 size_t materialIndex = couple->GetIndex();
323
324 G4double meanFreePath;
325 if (energy > highEnergyLimit) meanFreePath = meanFreePathTable->FindValue(highEnergyLimit,materialIndex);
326 else if (energy < lowEnergyLimit) meanFreePath = DBL_MAX;
327 else meanFreePath = meanFreePathTable->FindValue(energy,materialIndex);
328 return meanFreePath;
329}
330
331G4double G4LowEnergyGammaConversion::ScreenFunction1(G4double screenVariable)
332{
333 // Compute the value of the screening function 3*phi1 - phi2
334
335 G4double value;
336
337 if (screenVariable > 1.)
338 value = 42.24 - 8.368 * std::log(screenVariable + 0.952);
339 else
340 value = 42.392 - screenVariable * (7.796 - 1.961 * screenVariable);
341
342 return value;
343}
344
345G4double G4LowEnergyGammaConversion::ScreenFunction2(G4double screenVariable)
346{
347 // Compute the value of the screening function 1.5*phi1 - 0.5*phi2
348
349 G4double value;
350
351 if (screenVariable > 1.)
352 value = 42.24 - 8.368 * std::log(screenVariable + 0.952);
353 else
354 value = 41.405 - screenVariable * (5.828 - 0.8945 * screenVariable);
355
356 return value;
357}
Note: See TracBrowser for help on using the repository browser.