source: trunk/source/processes/electromagnetic/lowenergy/src/G4PenelopeGammaConversion.cc@ 1005

Last change on this file since 1005 was 819, checked in by garnier, 17 years ago

import all except CVS

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// --------------------------------------------------------------------
27//
28//
29//
30// --------------------------------------------------------------
31//
32// Author: L.Pandola
33// History:
34// --------
35// 02 Dec 2002 L.Pandola 1st implementation
36// 12 Feb 2003 MG Pia Migration to "cuts per region"
37// 10 Mar 2003 V.Ivanchenko Remove CutPerMaterial warning
38// 13 Mar 2003 L.Pandola Code "cleaned"
39// 25 Mar 2003 L.Pandola Changed the name of the database file to read
40// 24 Apr 2003 V.Ivanchenko Cut per region mfpt
41// 17 Mar 2004 L.Pandola Removed unnecessary calls to std::pow(a,b)
42// --------------------------------------------------------------
43
44#include "G4PenelopeGammaConversion.hh"
45
46#include "Randomize.hh"
47#include "G4ParticleDefinition.hh"
48#include "G4Track.hh"
49#include "G4Step.hh"
50#include "G4ForceCondition.hh"
51#include "G4Gamma.hh"
52#include "G4Electron.hh"
53#include "G4DynamicParticle.hh"
54#include "G4VParticleChange.hh"
55#include "G4ThreeVector.hh"
56#include "G4Positron.hh"
57#include "G4IonisParamElm.hh"
58#include "G4Material.hh"
59#include "G4VCrossSectionHandler.hh"
60#include "G4CrossSectionHandler.hh"
61#include "G4VEMDataSet.hh"
62#include "G4VDataSetAlgorithm.hh"
63#include "G4LogLogInterpolation.hh"
64#include "G4VRangeTest.hh"
65#include "G4RangeTest.hh"
66#include "G4MaterialCutsCouple.hh"
67
68
69G4PenelopeGammaConversion::G4PenelopeGammaConversion(const G4String& processName)
70 : G4VDiscreteProcess(processName),
71 lowEnergyLimit(1.022000*MeV),
72 highEnergyLimit(100*GeV),
73 intrinsicLowEnergyLimit(1.022000*MeV),
74 intrinsicHighEnergyLimit(100*GeV),
75 smallEnergy(1.1*MeV)
76
77{
78 if (lowEnergyLimit < intrinsicLowEnergyLimit ||
79 highEnergyLimit > intrinsicHighEnergyLimit)
80 {
81 G4Exception("G4PenelopeGammaConversion::G4PenelopeGammaConversion - energy limit outside intrinsic process validity range");
82 }
83
84 // The following pointer is owned by G4DataHandler
85 crossSectionHandler = new G4CrossSectionHandler();
86 // Log log interpolation (default)
87 crossSectionHandler->Initialise(0,1.0220*MeV,100.*GeV,400);
88 meanFreePathTable = 0;
89 rangeTest = new G4RangeTest;
90
91 if (verboseLevel > 0)
92 {
93 G4cout << GetProcessName() << " is created " << G4endl
94 << "Energy range: "
95 << lowEnergyLimit / MeV << " MeV - "
96 << highEnergyLimit / GeV << " GeV"
97 << G4endl;
98 }
99}
100
101G4PenelopeGammaConversion::~G4PenelopeGammaConversion()
102{
103 delete meanFreePathTable;
104 delete crossSectionHandler;
105 delete rangeTest;
106}
107
108void G4PenelopeGammaConversion::BuildPhysicsTable(const G4ParticleDefinition& )
109{
110
111 crossSectionHandler->Clear();
112 G4String crossSectionFile = "penelope/pp-cs-pen-";
113 crossSectionHandler->LoadData(crossSectionFile);
114 delete meanFreePathTable;
115 meanFreePathTable = crossSectionHandler->BuildMeanFreePathForMaterials();
116}
117
118G4VParticleChange* G4PenelopeGammaConversion::PostStepDoIt(const G4Track& aTrack,
119 const G4Step& aStep)
120{
121 aParticleChange.Initialize(aTrack);
122
123 const G4MaterialCutsCouple* couple = aTrack.GetMaterialCutsCouple();
124
125 const G4DynamicParticle* incidentPhoton = aTrack.GetDynamicParticle();
126 G4double photonEnergy = incidentPhoton->GetKineticEnergy();
127 G4ParticleMomentum photonDirection = incidentPhoton->GetMomentumDirection();
128
129 G4double eps ;
130 G4double eki = electron_mass_c2 / photonEnergy ;
131
132 // Do it fast if photon energy < 1.1 MeV
133 if (photonEnergy < smallEnergy )
134 {
135 eps = eki + (1-2*eki) * G4UniformRand();
136 }
137 else
138 {
139 // Select randomly one element in the current material
140 const G4Element* element = crossSectionHandler->SelectRandomElement(couple,photonEnergy);
141
142 if (element == 0)
143 {
144 G4cout << "G4PenelopeGammaConversion::PostStepDoIt - element = 0" << G4endl;
145 }
146 G4IonisParamElm* ionisation = element->GetIonisation();
147 if (ionisation == 0)
148 {
149 G4cout << "G4PenelopeGammaConversion::PostStepDoIt - ionisation = 0" << G4endl;
150 }
151
152 //Low energy and Coulomb corrections
153 G4double Z=ionisation->GetZ();
154 G4double ZAlpha = Z*fine_structure_const;
155 G4double ScreenRadius = GetScreeningRadius(Z);
156 G4double funct1=0,g0=0;
157 G4double g1min=0,g2min=0;
158 funct1 = 4.0*std::log(ScreenRadius);
159 g0 = funct1-4*CoulombCorrection(ZAlpha)+LowEnergyCorrection(ZAlpha,eki);
160 G4double bmin = 2*eki*ScreenRadius;
161 g1min=g0+ScreenFunction(bmin,1);
162 g2min=g0+ScreenFunction(bmin,2);
163 G4double xr,a1,p1;
164 xr=0.5-eki;
165 a1=(2.0/3.0)*g1min*xr*xr;
166 p1=a1/(a1+g2min);
167
168 //Random sampling of eps
169 G4double rand1,rand2,rand3,b;
170 G4double g1;
171
172 do{
173 rand1 = G4UniformRand();
174 if (rand1 < p1) {
175 rand2 = 2.0*G4UniformRand()-1.0;
176 if (rand2 < 0) {
177 eps = 0.5 - xr*std::pow(std::abs(rand2),(1./3.));
178 }
179 else
180 {
181 eps = 0.5 + xr*std::pow(rand2,(1./3.));
182 }
183 b = (eki*ScreenRadius)/(2*eps*(1.0-eps));
184 g1 = g0+ScreenFunction(b,1);
185 if (g1 < 0) g1=0;
186 rand3 = G4UniformRand()*g1min;
187 }
188 else
189 {
190 eps = eki+2.0*xr*G4UniformRand();
191 b = (eki*ScreenRadius)/(2*eps*(1.0-eps));
192 g1 = g0+ScreenFunction(b,2);
193 if (g1 < 0) g1=0;
194 rand3 = G4UniformRand()*g2min;
195 }
196 } while (rand3>g1);
197 } //End of eps sampling
198
199 G4double electronTotEnergy;
200 G4double positronTotEnergy;
201
202 electronTotEnergy = eps*photonEnergy;
203 positronTotEnergy = (1.0-eps)*photonEnergy;
204
205 // Scattered electron (positron) angles. ( Z - axis along the parent photon)
206
207 //electron kinematics
208 G4double costheta_el,costheta_po;
209 G4double phi_el,phi_po;
210 G4double electronKineEnergy = std::max(0.,electronTotEnergy - electron_mass_c2) ;
211 costheta_el = G4UniformRand()*2.0-1.0;
212 G4double kk = std::sqrt(electronKineEnergy*(electronKineEnergy+2.*electron_mass_c2));
213 costheta_el = (costheta_el*electronTotEnergy+kk)/(electronTotEnergy+costheta_el*kk);
214 phi_el = twopi * G4UniformRand() ;
215 G4double dirX_el = std::sqrt(1.-costheta_el*costheta_el) * std::cos(phi_el);
216 G4double dirY_el = std::sqrt(1.-costheta_el*costheta_el) * std::sin(phi_el);
217 G4double dirZ_el = costheta_el;
218
219 //positron kinematics
220 G4double positronKineEnergy = std::max(0.,positronTotEnergy - electron_mass_c2) ;
221 costheta_po = G4UniformRand()*2.0-1.0;
222 kk = std::sqrt(positronKineEnergy*(positronKineEnergy+2.*electron_mass_c2));
223 costheta_po = (costheta_po*positronTotEnergy+kk)/(positronTotEnergy+costheta_po*kk);
224 phi_po = twopi * G4UniformRand() ;
225 G4double dirX_po = std::sqrt(1.-costheta_po*costheta_po) * std::cos(phi_po);
226 G4double dirY_po = std::sqrt(1.-costheta_po*costheta_po) * std::sin(phi_po);
227 G4double dirZ_po = costheta_po;
228
229// Kinematics of the created pair:
230// the electron and positron are assumed to have a symetric angular
231// distribution with respect to the Z axis along the parent photon
232
233 G4double localEnergyDeposit = 0. ;
234
235 aParticleChange.SetNumberOfSecondaries(2) ;
236
237
238 // Generate the electron only if with large enough range w.r.t. cuts and safety
239
240 G4double safety = aStep.GetPostStepPoint()->GetSafety();
241
242 if (rangeTest->Escape(G4Electron::Electron(),couple,electronKineEnergy,safety))
243 {
244 G4ThreeVector electronDirection ( dirX_el, dirY_el, dirZ_el);
245 electronDirection.rotateUz(photonDirection);
246 G4DynamicParticle* particle1 = new G4DynamicParticle (G4Electron::Electron(),
247 electronDirection,
248 electronKineEnergy);
249 aParticleChange.AddSecondary(particle1) ;
250 }
251 else
252 {
253 localEnergyDeposit += electronKineEnergy ;
254 }
255
256
257 if (! (rangeTest->Escape(G4Positron::Positron(),couple,positronKineEnergy,safety)))
258 {
259 localEnergyDeposit += positronKineEnergy ;
260 positronKineEnergy = 0. ;
261 }
262 G4ThreeVector positronDirection(dirX_po,dirY_po,dirZ_po);
263 positronDirection.rotateUz(photonDirection);
264
265 // Create G4DynamicParticle object for the particle2
266 G4DynamicParticle* particle2 = new G4DynamicParticle(G4Positron::Positron(),
267 positronDirection, positronKineEnergy);
268 aParticleChange.AddSecondary(particle2) ;
269
270 aParticleChange.ProposeLocalEnergyDeposit(localEnergyDeposit) ;
271
272 // Kill the incident photon
273 aParticleChange.ProposeMomentumDirection(0.,0.,0.) ;
274 aParticleChange.ProposeEnergy(0.) ;
275 aParticleChange.ProposeTrackStatus(fStopAndKill) ;
276
277 // Reset NbOfInteractionLengthLeft and return aParticleChange
278 return G4VDiscreteProcess::PostStepDoIt(aTrack,aStep);
279}
280
281G4bool G4PenelopeGammaConversion::IsApplicable(const G4ParticleDefinition& particle)
282{
283 return ( &particle == G4Gamma::Gamma() );
284}
285
286G4double G4PenelopeGammaConversion::GetMeanFreePath(const G4Track& track,
287 G4double, // previousStepSize
288 G4ForceCondition*)
289{
290 const G4DynamicParticle* photon = track.GetDynamicParticle();
291 G4double energy = photon->GetKineticEnergy();
292 const G4MaterialCutsCouple* couple = track.GetMaterialCutsCouple();
293 size_t materialIndex = couple->GetIndex();
294
295 G4double meanFreePath;
296 if (energy > highEnergyLimit) meanFreePath = meanFreePathTable->FindValue(highEnergyLimit,materialIndex);
297 else if (energy < lowEnergyLimit) meanFreePath = DBL_MAX;
298 else meanFreePath = meanFreePathTable->FindValue(energy,materialIndex);
299 return meanFreePath;
300}
301
302G4double G4PenelopeGammaConversion::ScreenFunction(G4double b,G4int icase)
303{
304 G4double bsquare=b*b;
305 G4double a0,f1,f2,g1,g2;
306 f1=2.0-2*std::log(1+bsquare);
307 f2=f1-(2.0/3.0);
308 if (b < 1.0e-10)
309 {
310 f1=f1-twopi*b;
311 }
312 else
313 {
314 a0 = 4*b*std::atan(1.0/b);
315 f1 = f1 - a0;
316 f2 = f2+2*bsquare*(4.0-a0-3*std::log((1+bsquare)/bsquare));
317 }
318 g1=0.5*(3*f1-f2);
319 g2=0.25*(3*f1+f2);
320 if (icase==1) {
321 return g1;
322 }
323 else
324 {
325 return g2;
326 }
327}
328
329G4double G4PenelopeGammaConversion::CoulombCorrection(G4double a)
330{
331 G4double fc=0;
332 G4double b[7] = {0.202059,-0.03693,0.00835,-0.00201,0.00049,-0.00012,0.00003};
333 G4double aSquared = a*a;
334 G4double aFourth = aSquared*aSquared;
335 G4double aEighth = aFourth*aFourth;
336
337 fc = ((1.0/(1.0+a*a))+b[0]+b[1]*aSquared+b[2]*aFourth+b[3]*(aSquared*aFourth)+
338 b[4]*aEighth+b[5]*(aEighth*aSquared)+b[6]*(aEighth*aFourth));
339 fc=aSquared*fc;
340 return fc;
341}
342
343G4double G4PenelopeGammaConversion::LowEnergyCorrection(G4double a,G4double eki)
344{
345 G4double f0=0,t=0;
346 G4double b[12] = {-1.744,-12.10,11.18,8.523,73.26,-41.41,-13.52,-121.1,94.41,8.946,62.05,-63.41};
347 t=std::sqrt(2.0*eki);
348 G4double tSq = t*t;
349 f0=(b[0]+b[1]*a+b[2]*a*a)*t+(b[3]+b[4]*a+b[5]*a*a)*(tSq)+(b[6]+b[7]*a+b[8]*a*a)*(tSq*t)+
350 (b[9]+b[10]*a+b[11]*a*a)*(tSq*tSq);
351 return f0;
352}
353
354G4double G4PenelopeGammaConversion::GetScreeningRadius(G4double Z)
355{
356 char* path = getenv("G4LEDATA");
357 if (!path)
358 {
359 G4String excep = "G4PenelopeGammaConversion - G4LEDATA environment variable not set!";
360 G4Exception(excep);
361 }
362 G4String pathString(path);
363 G4String pathFile = pathString + "/penelope/pp-pen.dat";
364 std::ifstream file(pathFile);
365 std::filebuf* lsdp = file.rdbuf();
366
367 if (!(lsdp->is_open()))
368 {
369 G4String excep = "G4PenelopeGammaConversion - data file " + pathFile + "not found!";
370 G4Exception(excep);
371 }
372 G4int k;
373 G4double a1,a2;
374 while(!file.eof()) {
375 file >> k >> a1 >> a2;
376 if ((G4double) k == Z)
377 {
378 return a1;
379 }
380 }
381 G4String excep = "G4PenelopeGammaConversion - Screening Radius for not found in the data file";
382 G4Exception(excep);
383 return 0;
384}
Note: See TracBrowser for help on using the repository browser.