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

Last change on this file since 1315 was 1315, checked in by garnier, 14 years ago

update geant4-09-04-beta-cand-01 interfaces-V09-03-09 vis-V09-03-08

File size: 13.3 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// $Id: G4LowEnergyGammaConversion.cc,v 1.39 2009/06/11 15:47:08 mantero Exp $
29// GEANT4 tag $Name: geant4-09-04-beta-cand-01 $
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   G4cout << G4endl;
108   G4cout << "*******************************************************************************" << G4endl;
109   G4cout << "*******************************************************************************" << G4endl;
110   G4cout << "   The class G4LowEnergyGammaConversion is NOT SUPPORTED ANYMORE. " << G4endl;
111   G4cout << "   It will be REMOVED with the next major release of Geant4. " << G4endl;
112   G4cout << "   Please consult: https://twiki.cern.ch/twiki/bin/view/Geant4/LoweProcesses" << G4endl;
113   G4cout << "*******************************************************************************" << G4endl;
114   G4cout << "*******************************************************************************" << G4endl;
115   G4cout << G4endl;
116}
117 
118G4LowEnergyGammaConversion::~G4LowEnergyGammaConversion()
119{
120  delete meanFreePathTable;
121  delete crossSectionHandler;
122  delete rangeTest;
123}
124
125void G4LowEnergyGammaConversion::BuildPhysicsTable(const G4ParticleDefinition& )
126{
127
128  crossSectionHandler->Clear();
129  G4String crossSectionFile = "pair/pp-cs-";
130  crossSectionHandler->LoadData(crossSectionFile);
131
132  delete meanFreePathTable;
133  meanFreePathTable = crossSectionHandler->BuildMeanFreePathForMaterials();
134}
135
136G4VParticleChange* G4LowEnergyGammaConversion::PostStepDoIt(const G4Track& aTrack,
137                                                            const G4Step& aStep)
138{
139// The energies of the e+ e- secondaries are sampled using the Bethe - Heitler
140// cross sections with Coulomb correction. A modified version of the random
141// number techniques of Butcher & Messel is used (Nuc Phys 20(1960),15).
142
143// Note 1 : Effects due to the breakdown of the Born approximation at low
144// energy are ignored.
145// Note 2 : The differential cross section implicitly takes account of
146// pair creation in both nuclear and atomic electron fields. However triplet
147// prodution is not generated.
148
149  aParticleChange.Initialize(aTrack);
150
151  const G4MaterialCutsCouple* couple = aTrack.GetMaterialCutsCouple();
152
153  const G4DynamicParticle* incidentPhoton = aTrack.GetDynamicParticle();
154  G4double photonEnergy = incidentPhoton->GetKineticEnergy();
155  G4ParticleMomentum photonDirection = incidentPhoton->GetMomentumDirection();
156
157  G4double epsilon ;
158  G4double epsilon0 = electron_mass_c2 / photonEnergy ;
159
160  // Do it fast if photon energy < 2. MeV
161  if (photonEnergy < smallEnergy )
162    {
163      epsilon = epsilon0 + (0.5 - epsilon0) * G4UniformRand();
164    }
165  else
166    {
167      // Select randomly one element in the current material
168      const G4Element* element = crossSectionHandler->SelectRandomElement(couple,photonEnergy);
169
170      if (element == 0)
171        {
172          G4cout << "G4LowEnergyGammaConversion::PostStepDoIt - element = 0" << G4endl;
173        }
174      G4IonisParamElm* ionisation = element->GetIonisation();
175       if (ionisation == 0)
176        {
177          G4cout << "G4LowEnergyGammaConversion::PostStepDoIt - ionisation = 0" << G4endl;
178        }
179
180      // Extract Coulomb factor for this Element
181      G4double fZ = 8. * (ionisation->GetlogZ3());
182      if (photonEnergy > 50. * MeV) fZ += 8. * (element->GetfCoulomb());
183
184      // Limits of the screening variable
185      G4double screenFactor = 136. * epsilon0 / (element->GetIonisation()->GetZ3()) ;
186      G4double screenMax = std::exp ((42.24 - fZ)/8.368) - 0.952 ;
187      G4double screenMin = std::min(4.*screenFactor,screenMax) ;
188
189      // Limits of the energy sampling
190      G4double epsilon1 = 0.5 - 0.5 * std::sqrt(1. - screenMin / screenMax) ;
191      G4double epsilonMin = std::max(epsilon0,epsilon1);
192      G4double epsilonRange = 0.5 - epsilonMin ;
193
194      // Sample the energy rate of the created electron (or positron)
195      G4double screen;
196      G4double gReject ;
197
198      G4double f10 = ScreenFunction1(screenMin) - fZ;
199      G4double f20 = ScreenFunction2(screenMin) - fZ;
200      G4double normF1 = std::max(f10 * epsilonRange * epsilonRange,0.);
201      G4double normF2 = std::max(1.5 * f20,0.);
202
203      do {
204        if (normF1 / (normF1 + normF2) > G4UniformRand() )
205          {
206            epsilon = 0.5 - epsilonRange * std::pow(G4UniformRand(), 0.3333) ;
207            screen = screenFactor / (epsilon * (1. - epsilon));
208            gReject = (ScreenFunction1(screen) - fZ) / f10 ;
209          }
210        else
211          {
212            epsilon = epsilonMin + epsilonRange * G4UniformRand();
213            screen = screenFactor / (epsilon * (1 - epsilon));
214            gReject = (ScreenFunction2(screen) - fZ) / f20 ;
215          }
216      } while ( gReject < G4UniformRand() );
217
218    }   //  End of epsilon sampling
219
220  // Fix charges randomly
221
222  G4double electronTotEnergy;
223  G4double positronTotEnergy;
224
225  if (CLHEP::RandBit::shootBit())
226    {
227      electronTotEnergy = (1. - epsilon) * photonEnergy;
228      positronTotEnergy = epsilon * photonEnergy;
229    }
230  else
231    {
232      positronTotEnergy = (1. - epsilon) * photonEnergy;
233      electronTotEnergy = epsilon * photonEnergy;
234    }
235
236  // Scattered electron (positron) angles. ( Z - axis along the parent photon)
237  // Universal distribution suggested by L. Urban (Geant3 manual (1993) Phys211),
238  // derived from Tsai distribution (Rev. Mod. Phys. 49, 421 (1977)
239
240  G4double u;
241  const G4double a1 = 0.625;
242  G4double a2 = 3. * a1;
243  //  G4double d = 27. ;
244
245  //  if (9. / (9. + d) > G4UniformRand())
246  if (0.25 > G4UniformRand())
247    {
248      u = - std::log(G4UniformRand() * G4UniformRand()) / a1 ;
249    }
250  else
251    {
252      u = - std::log(G4UniformRand() * G4UniformRand()) / a2 ;
253    }
254
255  G4double thetaEle = u*electron_mass_c2/electronTotEnergy;
256  G4double thetaPos = u*electron_mass_c2/positronTotEnergy;
257  G4double phi  = twopi * G4UniformRand();
258
259  G4double dxEle= std::sin(thetaEle)*std::cos(phi),dyEle= std::sin(thetaEle)*std::sin(phi),dzEle=std::cos(thetaEle);
260  G4double dxPos=-std::sin(thetaPos)*std::cos(phi),dyPos=-std::sin(thetaPos)*std::sin(phi),dzPos=std::cos(thetaPos);
261 
262 
263  // Kinematics of the created pair:
264  // the electron and positron are assumed to have a symetric angular
265  // distribution with respect to the Z axis along the parent photon
266 
267  G4double localEnergyDeposit = 0. ;
268 
269  aParticleChange.SetNumberOfSecondaries(2) ; 
270  G4double electronKineEnergy = std::max(0.,electronTotEnergy - electron_mass_c2) ;
271 
272  // Generate the electron only if with large enough range w.r.t. cuts and safety
273 
274  G4double safety = aStep.GetPostStepPoint()->GetSafety();
275 
276  if (rangeTest->Escape(G4Electron::Electron(),couple,electronKineEnergy,safety))
277    {
278      G4ThreeVector electronDirection (dxEle, dyEle, dzEle);
279      electronDirection.rotateUz(photonDirection);
280     
281      G4DynamicParticle* particle1 = new G4DynamicParticle (G4Electron::Electron(),
282                                                            electronDirection,
283                                                            electronKineEnergy);
284      aParticleChange.AddSecondary(particle1) ;
285    }
286  else
287    {
288      localEnergyDeposit += electronKineEnergy ;
289    }
290
291  // The e+ is always created (even with kinetic energy = 0) for further annihilation
292  G4double positronKineEnergy = std::max(0.,positronTotEnergy - electron_mass_c2) ;
293
294  // Is the local energy deposit correct, if the positron is always created?
295  if (! (rangeTest->Escape(G4Positron::Positron(),couple,positronKineEnergy,safety)))
296    {
297      localEnergyDeposit += positronKineEnergy ;
298      positronKineEnergy = 0. ;
299    }
300
301  G4ThreeVector positronDirection (dxPos, dyPos, dzPos);
302  positronDirection.rotateUz(photonDirection);   
303 
304  // Create G4DynamicParticle object for the particle2
305  G4DynamicParticle* particle2 = new G4DynamicParticle(G4Positron::Positron(),
306                                                       positronDirection, positronKineEnergy);
307  aParticleChange.AddSecondary(particle2) ; 
308
309  aParticleChange.ProposeLocalEnergyDeposit(localEnergyDeposit) ;
310 
311  // Kill the incident photon
312  aParticleChange.ProposeMomentumDirection(0.,0.,0.) ;
313  aParticleChange.ProposeEnergy(0.) ; 
314  aParticleChange.ProposeTrackStatus(fStopAndKill) ;
315
316  //  Reset NbOfInteractionLengthLeft and return aParticleChange
317  return G4VDiscreteProcess::PostStepDoIt(aTrack,aStep);
318}
319
320G4bool G4LowEnergyGammaConversion::IsApplicable(const G4ParticleDefinition& particle)
321{
322  return ( &particle == G4Gamma::Gamma() ); 
323}
324
325G4double G4LowEnergyGammaConversion::GetMeanFreePath(const G4Track& track, 
326                                                     G4double, // previousStepSize
327                                                     G4ForceCondition*)
328{
329  const G4DynamicParticle* photon = track.GetDynamicParticle();
330  G4double energy = photon->GetKineticEnergy();
331  const G4MaterialCutsCouple* couple = track.GetMaterialCutsCouple();
332  size_t materialIndex = couple->GetIndex();
333
334  G4double meanFreePath;
335  if (energy > highEnergyLimit) meanFreePath = meanFreePathTable->FindValue(highEnergyLimit,materialIndex);
336  else if (energy < lowEnergyLimit) meanFreePath = DBL_MAX;
337  else meanFreePath = meanFreePathTable->FindValue(energy,materialIndex);
338  return meanFreePath;
339}
340
341G4double G4LowEnergyGammaConversion::ScreenFunction1(G4double screenVariable)
342{
343  // Compute the value of the screening function 3*phi1 - phi2
344
345  G4double value;
346 
347  if (screenVariable > 1.)
348    value = 42.24 - 8.368 * std::log(screenVariable + 0.952);
349  else
350    value = 42.392 - screenVariable * (7.796 - 1.961 * screenVariable);
351 
352  return value;
353} 
354
355G4double G4LowEnergyGammaConversion::ScreenFunction2(G4double screenVariable)
356{
357  // Compute the value of the screening function 1.5*phi1 - 0.5*phi2
358 
359  G4double value;
360 
361  if (screenVariable > 1.)
362    value = 42.24 - 8.368 * std::log(screenVariable + 0.952);
363  else
364    value = 41.405 - screenVariable * (5.828 - 0.8945 * screenVariable);
365 
366  return value;
367} 
Note: See TracBrowser for help on using the repository browser.