source: trunk/source/processes/electromagnetic/lowenergy/src/G4LivermoreGammaConversionModel.cc @ 989

Last change on this file since 989 was 968, checked in by garnier, 15 years ago

fichier ajoutes

File size: 13.0 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.2 2009/01/21 10:58:13 sincerti Exp $
27// GEANT4 tag $Name: geant4-09-02-ref-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),crossSectionHandler(0),meanFreePathTable(0)
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  if (meanFreePathTable) delete meanFreePathTable;
66  if (crossSectionHandler) 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  if (crossSectionHandler)
78  {
79    crossSectionHandler->Clear();
80    delete crossSectionHandler;
81  }
82
83  // Energy limits
84 
85  if (LowEnergyLimit() < lowEnergyLimit)
86  {
87    G4cout << "G4LivermoreGammaConversionModel: low energy limit increased from " << 
88        LowEnergyLimit()/eV << " eV to " << lowEnergyLimit << " eV" << G4endl;
89    SetLowEnergyLimit(lowEnergyLimit);
90  }
91
92  if (HighEnergyLimit() > highEnergyLimit)
93  {
94    G4cout << "G4LivermoreGammaConversionModel: high energy limit decreased from " << 
95        HighEnergyLimit()/GeV << " GeV to " << highEnergyLimit << " GeV" << G4endl;
96    SetHighEnergyLimit(highEnergyLimit);
97  }
98
99  // Read data tables for all materials
100 
101  crossSectionHandler = new G4CrossSectionHandler();
102  crossSectionHandler->Initialise(0,1.0220*MeV,100.*GeV,400);
103  G4String crossSectionFile = "pair/pp-cs-";
104  crossSectionHandler->LoadData(crossSectionFile);
105
106  meanFreePathTable = 0;
107  meanFreePathTable = crossSectionHandler->BuildMeanFreePathForMaterials();
108
109  //
110 
111  if (verboseLevel > 2) 
112    G4cout << "Loaded cross section files for PenelopeGammaConversion" << G4endl;
113
114  InitialiseElementSelectors(particle,cuts);
115
116  G4cout << "Livermore Gamma Conversion model is initialized " << G4endl
117         << "Energy range: "
118         << LowEnergyLimit() / MeV << " MeV - "
119         << HighEnergyLimit() / GeV << " GeV"
120         << G4endl;
121
122  if(isInitialised) return;
123
124  if(pParticleChange)
125    fParticleChange = reinterpret_cast<G4ParticleChangeForGamma*>(pParticleChange);
126  else
127    fParticleChange = new G4ParticleChangeForGamma();
128  isInitialised = true;}
129
130//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
131
132G4double G4LivermoreGammaConversionModel::ComputeCrossSectionPerAtom(
133                                       const G4ParticleDefinition*,
134                                             G4double GammaEnergy,
135                                             G4double Z, G4double,
136                                             G4double, G4double)
137{
138  if (verboseLevel > 3)
139    G4cout << "Calling ComputeCrossSectionPerAtom() of G4LivermoreGammaConversionModel" << G4endl;
140
141  G4double cs = crossSectionHandler->FindValue(G4int(Z), GammaEnergy);
142  return cs;
143}
144
145//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
146
147void G4LivermoreGammaConversionModel::SampleSecondaries(std::vector<G4DynamicParticle*>* fvect,
148                                              const G4MaterialCutsCouple* couple,
149                                              const G4DynamicParticle* aDynamicGamma,
150                                              G4double,
151                                              G4double)
152{
153
154// The energies of the e+ e- secondaries are sampled using the Bethe - Heitler
155// cross sections with Coulomb correction. A modified version of the random
156// number techniques of Butcher & Messel is used (Nuc Phys 20(1960),15).
157
158// Note 1 : Effects due to the breakdown of the Born approximation at low
159// energy are ignored.
160// Note 2 : The differential cross section implicitly takes account of
161// pair creation in both nuclear and atomic electron fields. However triplet
162// prodution is not generated.
163
164  if (verboseLevel > 3)
165    G4cout << "Calling SampleSecondaries() of G4LivermoreGammaConversionModel" << G4endl;
166
167  G4double photonEnergy = aDynamicGamma->GetKineticEnergy();
168  G4ParticleMomentum photonDirection = aDynamicGamma->GetMomentumDirection();
169
170  G4double epsilon ;
171  G4double epsilon0 = electron_mass_c2 / photonEnergy ;
172
173  // Do it fast if photon energy < 2. MeV
174  if (photonEnergy < smallEnergy )
175    {
176      epsilon = epsilon0 + (0.5 - epsilon0) * G4UniformRand();
177    }
178  else
179    {
180      // Select randomly one element in the current material
181      const G4Element* element = crossSectionHandler->SelectRandomElement(couple,photonEnergy);
182
183      if (element == 0)
184        {
185          G4cout << "G4LivermoreGammaConversionModel::SampleSecondaries - element = 0" << G4endl;
186        }
187      G4IonisParamElm* ionisation = element->GetIonisation();
188       if (ionisation == 0)
189        {
190          G4cout << "G4LivermoreGammaConversionModel::SampleSecondaries - ionisation = 0" << G4endl;
191        }
192
193      // Extract Coulomb factor for this Element
194      G4double fZ = 8. * (ionisation->GetlogZ3());
195      if (photonEnergy > 50. * MeV) fZ += 8. * (element->GetfCoulomb());
196
197      // Limits of the screening variable
198      G4double screenFactor = 136. * epsilon0 / (element->GetIonisation()->GetZ3()) ;
199      G4double screenMax = std::exp ((42.24 - fZ)/8.368) - 0.952 ;
200      G4double screenMin = std::min(4.*screenFactor,screenMax) ;
201
202      // Limits of the energy sampling
203      G4double epsilon1 = 0.5 - 0.5 * std::sqrt(1. - screenMin / screenMax) ;
204      G4double epsilonMin = std::max(epsilon0,epsilon1);
205      G4double epsilonRange = 0.5 - epsilonMin ;
206
207      // Sample the energy rate of the created electron (or positron)
208      G4double screen;
209      G4double gReject ;
210
211      G4double f10 = ScreenFunction1(screenMin) - fZ;
212      G4double f20 = ScreenFunction2(screenMin) - fZ;
213      G4double normF1 = std::max(f10 * epsilonRange * epsilonRange,0.);
214      G4double normF2 = std::max(1.5 * f20,0.);
215
216      do {
217        if (normF1 / (normF1 + normF2) > G4UniformRand() )
218          {
219            epsilon = 0.5 - epsilonRange * std::pow(G4UniformRand(), 0.3333) ;
220            screen = screenFactor / (epsilon * (1. - epsilon));
221            gReject = (ScreenFunction1(screen) - fZ) / f10 ;
222          }
223        else
224          {
225            epsilon = epsilonMin + epsilonRange * G4UniformRand();
226            screen = screenFactor / (epsilon * (1 - epsilon));
227            gReject = (ScreenFunction2(screen) - fZ) / f20 ;
228          }
229      } while ( gReject < G4UniformRand() );
230
231    }   //  End of epsilon sampling
232
233  // Fix charges randomly
234
235  G4double electronTotEnergy;
236  G4double positronTotEnergy;
237
238  if (CLHEP::RandBit::shootBit())
239    {
240      electronTotEnergy = (1. - epsilon) * photonEnergy;
241      positronTotEnergy = epsilon * photonEnergy;
242    }
243  else
244    {
245      positronTotEnergy = (1. - epsilon) * photonEnergy;
246      electronTotEnergy = epsilon * photonEnergy;
247    }
248
249  // Scattered electron (positron) angles. ( Z - axis along the parent photon)
250  // Universal distribution suggested by L. Urban (Geant3 manual (1993) Phys211),
251  // derived from Tsai distribution (Rev. Mod. Phys. 49, 421 (1977)
252
253  G4double u;
254  const G4double a1 = 0.625;
255  G4double a2 = 3. * a1;
256  //  G4double d = 27. ;
257
258  //  if (9. / (9. + d) > G4UniformRand())
259  if (0.25 > G4UniformRand())
260    {
261      u = - std::log(G4UniformRand() * G4UniformRand()) / a1 ;
262    }
263  else
264    {
265      u = - std::log(G4UniformRand() * G4UniformRand()) / a2 ;
266    }
267
268  G4double thetaEle = u*electron_mass_c2/electronTotEnergy;
269  G4double thetaPos = u*electron_mass_c2/positronTotEnergy;
270  G4double phi  = twopi * G4UniformRand();
271
272  G4double dxEle= std::sin(thetaEle)*std::cos(phi),dyEle= std::sin(thetaEle)*std::sin(phi),dzEle=std::cos(thetaEle);
273  G4double dxPos=-std::sin(thetaPos)*std::cos(phi),dyPos=-std::sin(thetaPos)*std::sin(phi),dzPos=std::cos(thetaPos);
274 
275 
276  // Kinematics of the created pair:
277  // the electron and positron are assumed to have a symetric angular
278  // distribution with respect to the Z axis along the parent photon
279 
280//  aParticleChange.SetNumberOfSecondaries(2) ;
281  G4double electronKineEnergy = std::max(0.,electronTotEnergy - electron_mass_c2) ;
282 
283// SI - The range test has been removed wrt original G4LowEnergyGammaconversion class
284
285  G4ThreeVector electronDirection (dxEle, dyEle, dzEle);
286  electronDirection.rotateUz(photonDirection);
287     
288  G4DynamicParticle* particle1 = new G4DynamicParticle (G4Electron::Electron(),
289                                                            electronDirection,
290                                                            electronKineEnergy);
291
292  // The e+ is always created (even with kinetic energy = 0) for further annihilation
293  G4double positronKineEnergy = std::max(0.,positronTotEnergy - electron_mass_c2) ;
294
295// SI - The range test has been removed wrt original G4LowEnergyGammaconversion class
296
297  G4ThreeVector positronDirection (dxPos, dyPos, dzPos);
298  positronDirection.rotateUz(photonDirection);   
299 
300  // Create G4DynamicParticle object for the particle2
301  G4DynamicParticle* particle2 = new G4DynamicParticle(G4Positron::Positron(),
302                                                       positronDirection, positronKineEnergy);
303  // Fill output vector
304  fvect->push_back(particle1);
305  fvect->push_back(particle2);
306
307  // kill incident photon
308  fParticleChange->SetProposedKineticEnergy(0.);
309  fParticleChange->ProposeTrackStatus(fStopAndKill);   
310
311}
312
313//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
314
315G4double G4LivermoreGammaConversionModel::ScreenFunction1(G4double screenVariable)
316{
317  // Compute the value of the screening function 3*phi1 - phi2
318
319  G4double value;
320 
321  if (screenVariable > 1.)
322    value = 42.24 - 8.368 * std::log(screenVariable + 0.952);
323  else
324    value = 42.392 - screenVariable * (7.796 - 1.961 * screenVariable);
325 
326  return value;
327} 
328
329//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
330
331G4double G4LivermoreGammaConversionModel::ScreenFunction2(G4double screenVariable)
332{
333  // Compute the value of the screening function 1.5*phi1 - 0.5*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 = 41.405 - screenVariable * (5.828 - 0.8945 * screenVariable);
341 
342  return value;
343} 
344
345//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
346
347G4double G4LivermoreGammaConversionModel::GetMeanFreePath(const G4Track& track, 
348                                                     G4double, // previousStepSize
349                                                     G4ForceCondition*)
350{
351  const G4DynamicParticle* photon = track.GetDynamicParticle();
352  G4double energy = photon->GetKineticEnergy();
353  const G4MaterialCutsCouple* couple = track.GetMaterialCutsCouple();
354  size_t materialIndex = couple->GetIndex();
355
356  G4double meanFreePath;
357  if (energy > highEnergyLimit) meanFreePath = meanFreePathTable->FindValue(highEnergyLimit,materialIndex);
358  else if (energy < lowEnergyLimit) meanFreePath = DBL_MAX;
359  else meanFreePath = meanFreePathTable->FindValue(energy,materialIndex);
360  return meanFreePath;
361}
362
Note: See TracBrowser for help on using the repository browser.