source: trunk/source/processes/electromagnetic/lowenergy/src/G4LivermoreGammaConversionModel.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: G4LivermoreGammaConversionModel.cc,v 1.1 2008/10/30 14:16:35 sincerti Exp $
27// GEANT4 tag $Name: geant4-09-02 $
28//
29
30#include "G4LivermoreGammaConversionModel.hh"
31
32//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
33
34using namespace std;
35
36//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
37
38G4LivermoreGammaConversionModel::G4LivermoreGammaConversionModel(const G4ParticleDefinition*,
39 const G4String& nam)
40:G4VEmModel(nam),smallEnergy(2.*MeV),isInitialised(false)
41{
42 lowEnergyLimit = 1.022000 * MeV;
43 highEnergyLimit = 100 * GeV;
44
45 G4cout << "Livermore Gamma conversion is constructed " << G4endl
46 << "Energy range: "
47 << lowEnergyLimit / keV << " keV - "
48 << highEnergyLimit / GeV << " GeV"
49 << G4endl;
50
51 verboseLevel= 0;
52 // Verbosity scale:
53 // 0 = nothing
54 // 1 = warning for energy non-conservation
55 // 2 = details of energy budget
56 // 3 = calculation of cross sections, file openings, sampling of atoms
57 // 4 = entering in methods
58
59}
60
61//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
62
63G4LivermoreGammaConversionModel::~G4LivermoreGammaConversionModel()
64{
65 delete meanFreePathTable;
66 delete crossSectionHandler;
67}
68
69//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
70
71void G4LivermoreGammaConversionModel::Initialise(const G4ParticleDefinition* particle,
72 const G4DataVector& cuts)
73{
74 if (verboseLevel > 3)
75 G4cout << "Calling G4LivermoreGammaConversionModel::Initialise()" << G4endl;
76
77 InitialiseElementSelectors(particle,cuts);
78
79 // Energy limits
80
81 if (LowEnergyLimit() < lowEnergyLimit)
82 {
83 G4cout << "G4LivermoreGammaConversionModel: low energy limit increased from " <<
84 LowEnergyLimit()/eV << " eV to " << lowEnergyLimit << " eV" << G4endl;
85 SetLowEnergyLimit(lowEnergyLimit);
86 }
87
88 if (HighEnergyLimit() > highEnergyLimit)
89 {
90 G4cout << "G4LivermoreGammaConversionModel: high energy limit decreased from " <<
91 HighEnergyLimit()/GeV << " GeV to " << highEnergyLimit << " GeV" << G4endl;
92 SetHighEnergyLimit(highEnergyLimit);
93 }
94
95 // Read data tables for all materials
96
97 crossSectionHandler = new G4CrossSectionHandler();
98 crossSectionHandler->Initialise(0,1.0220*MeV,100.*GeV,400);
99 G4String crossSectionFile = "pair/pp-cs-";
100 crossSectionHandler->LoadData(crossSectionFile);
101
102 meanFreePathTable = 0;
103 meanFreePathTable = crossSectionHandler->BuildMeanFreePathForMaterials();
104
105 //
106
107 if (verboseLevel > 2)
108 G4cout << "Loaded cross section files for PenelopeGammaConversion" << G4endl;
109
110 G4cout << "Livermore Gamma Conversion model is initialized " << G4endl
111 << "Energy range: "
112 << LowEnergyLimit() / MeV << " MeV - "
113 << HighEnergyLimit() / GeV << " GeV"
114 << G4endl;
115
116 if(isInitialised) return;
117
118 if(pParticleChange)
119 fParticleChange = reinterpret_cast<G4ParticleChangeForGamma*>(pParticleChange);
120 else
121 fParticleChange = new G4ParticleChangeForGamma();
122 isInitialised = true;}
123
124//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
125
126G4double G4LivermoreGammaConversionModel::ComputeCrossSectionPerAtom(
127 const G4ParticleDefinition*,
128 G4double GammaEnergy,
129 G4double Z, G4double,
130 G4double, G4double)
131{
132 if (verboseLevel > 3)
133 G4cout << "Calling ComputeCrossSectionPerAtom() of G4LivermoreGammaConversionModel" << G4endl;
134
135 G4double cs = crossSectionHandler->FindValue(G4int(Z), GammaEnergy);
136 return cs;
137}
138
139//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
140
141void G4LivermoreGammaConversionModel::SampleSecondaries(std::vector<G4DynamicParticle*>* fvect,
142 const G4MaterialCutsCouple* couple,
143 const G4DynamicParticle* aDynamicGamma,
144 G4double,
145 G4double)
146{
147
148// The energies of the e+ e- secondaries are sampled using the Bethe - Heitler
149// cross sections with Coulomb correction. A modified version of the random
150// number techniques of Butcher & Messel is used (Nuc Phys 20(1960),15).
151
152// Note 1 : Effects due to the breakdown of the Born approximation at low
153// energy are ignored.
154// Note 2 : The differential cross section implicitly takes account of
155// pair creation in both nuclear and atomic electron fields. However triplet
156// prodution is not generated.
157
158 if (verboseLevel > 3)
159 G4cout << "Calling SampleSecondaries() of G4LivermoreGammaConversionModel" << G4endl;
160
161 G4double photonEnergy = aDynamicGamma->GetKineticEnergy();
162 G4ParticleMomentum photonDirection = aDynamicGamma->GetMomentumDirection();
163
164 G4double epsilon ;
165 G4double epsilon0 = electron_mass_c2 / photonEnergy ;
166
167 // Do it fast if photon energy < 2. MeV
168 if (photonEnergy < smallEnergy )
169 {
170 epsilon = epsilon0 + (0.5 - epsilon0) * G4UniformRand();
171 }
172 else
173 {
174 // Select randomly one element in the current material
175 const G4Element* element = crossSectionHandler->SelectRandomElement(couple,photonEnergy);
176
177 if (element == 0)
178 {
179 G4cout << "G4LivermoreGammaConversionModel::SampleSecondaries - element = 0" << G4endl;
180 }
181 G4IonisParamElm* ionisation = element->GetIonisation();
182 if (ionisation == 0)
183 {
184 G4cout << "G4LivermoreGammaConversionModel::SampleSecondaries - ionisation = 0" << G4endl;
185 }
186
187 // Extract Coulomb factor for this Element
188 G4double fZ = 8. * (ionisation->GetlogZ3());
189 if (photonEnergy > 50. * MeV) fZ += 8. * (element->GetfCoulomb());
190
191 // Limits of the screening variable
192 G4double screenFactor = 136. * epsilon0 / (element->GetIonisation()->GetZ3()) ;
193 G4double screenMax = std::exp ((42.24 - fZ)/8.368) - 0.952 ;
194 G4double screenMin = std::min(4.*screenFactor,screenMax) ;
195
196 // Limits of the energy sampling
197 G4double epsilon1 = 0.5 - 0.5 * std::sqrt(1. - screenMin / screenMax) ;
198 G4double epsilonMin = std::max(epsilon0,epsilon1);
199 G4double epsilonRange = 0.5 - epsilonMin ;
200
201 // Sample the energy rate of the created electron (or positron)
202 G4double screen;
203 G4double gReject ;
204
205 G4double f10 = ScreenFunction1(screenMin) - fZ;
206 G4double f20 = ScreenFunction2(screenMin) - fZ;
207 G4double normF1 = std::max(f10 * epsilonRange * epsilonRange,0.);
208 G4double normF2 = std::max(1.5 * f20,0.);
209
210 do {
211 if (normF1 / (normF1 + normF2) > G4UniformRand() )
212 {
213 epsilon = 0.5 - epsilonRange * std::pow(G4UniformRand(), 0.3333) ;
214 screen = screenFactor / (epsilon * (1. - epsilon));
215 gReject = (ScreenFunction1(screen) - fZ) / f10 ;
216 }
217 else
218 {
219 epsilon = epsilonMin + epsilonRange * G4UniformRand();
220 screen = screenFactor / (epsilon * (1 - epsilon));
221 gReject = (ScreenFunction2(screen) - fZ) / f20 ;
222 }
223 } while ( gReject < G4UniformRand() );
224
225 } // End of epsilon sampling
226
227 // Fix charges randomly
228
229 G4double electronTotEnergy;
230 G4double positronTotEnergy;
231
232 if (CLHEP::RandBit::shootBit())
233 {
234 electronTotEnergy = (1. - epsilon) * photonEnergy;
235 positronTotEnergy = epsilon * photonEnergy;
236 }
237 else
238 {
239 positronTotEnergy = (1. - epsilon) * photonEnergy;
240 electronTotEnergy = epsilon * photonEnergy;
241 }
242
243 // Scattered electron (positron) angles. ( Z - axis along the parent photon)
244 // Universal distribution suggested by L. Urban (Geant3 manual (1993) Phys211),
245 // derived from Tsai distribution (Rev. Mod. Phys. 49, 421 (1977)
246
247 G4double u;
248 const G4double a1 = 0.625;
249 G4double a2 = 3. * a1;
250 // G4double d = 27. ;
251
252 // if (9. / (9. + d) > G4UniformRand())
253 if (0.25 > G4UniformRand())
254 {
255 u = - std::log(G4UniformRand() * G4UniformRand()) / a1 ;
256 }
257 else
258 {
259 u = - std::log(G4UniformRand() * G4UniformRand()) / a2 ;
260 }
261
262 G4double thetaEle = u*electron_mass_c2/electronTotEnergy;
263 G4double thetaPos = u*electron_mass_c2/positronTotEnergy;
264 G4double phi = twopi * G4UniformRand();
265
266 G4double dxEle= std::sin(thetaEle)*std::cos(phi),dyEle= std::sin(thetaEle)*std::sin(phi),dzEle=std::cos(thetaEle);
267 G4double dxPos=-std::sin(thetaPos)*std::cos(phi),dyPos=-std::sin(thetaPos)*std::sin(phi),dzPos=std::cos(thetaPos);
268
269
270 // Kinematics of the created pair:
271 // the electron and positron are assumed to have a symetric angular
272 // distribution with respect to the Z axis along the parent photon
273
274// aParticleChange.SetNumberOfSecondaries(2) ;
275 G4double electronKineEnergy = std::max(0.,electronTotEnergy - electron_mass_c2) ;
276
277// SI - The range test has been removed wrt original G4LowEnergyGammaconversion class
278
279 G4ThreeVector electronDirection (dxEle, dyEle, dzEle);
280 electronDirection.rotateUz(photonDirection);
281
282 G4DynamicParticle* particle1 = new G4DynamicParticle (G4Electron::Electron(),
283 electronDirection,
284 electronKineEnergy);
285
286 // The e+ is always created (even with kinetic energy = 0) for further annihilation
287 G4double positronKineEnergy = std::max(0.,positronTotEnergy - electron_mass_c2) ;
288
289// SI - The range test has been removed wrt original G4LowEnergyGammaconversion class
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 // Fill output vector
298 fvect->push_back(particle1);
299 fvect->push_back(particle2);
300
301 // kill incident photon
302 fParticleChange->SetProposedKineticEnergy(0.);
303 fParticleChange->ProposeTrackStatus(fStopAndKill);
304
305}
306
307//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
308
309G4double G4LivermoreGammaConversionModel::ScreenFunction1(G4double screenVariable)
310{
311 // Compute the value of the screening function 3*phi1 - phi2
312
313 G4double value;
314
315 if (screenVariable > 1.)
316 value = 42.24 - 8.368 * std::log(screenVariable + 0.952);
317 else
318 value = 42.392 - screenVariable * (7.796 - 1.961 * screenVariable);
319
320 return value;
321}
322
323//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
324
325G4double G4LivermoreGammaConversionModel::ScreenFunction2(G4double screenVariable)
326{
327 // Compute the value of the screening function 1.5*phi1 - 0.5*phi2
328
329 G4double value;
330
331 if (screenVariable > 1.)
332 value = 42.24 - 8.368 * std::log(screenVariable + 0.952);
333 else
334 value = 41.405 - screenVariable * (5.828 - 0.8945 * screenVariable);
335
336 return value;
337}
338
339//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
340
341G4double G4LivermoreGammaConversionModel::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.