source: trunk/source/processes/electromagnetic/standard/include/G4PairProductionRelModel.hh

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

update ti head

File size: 8.1 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: G4PairProductionRelModel.hh,v 1.9 2010/10/26 10:35:22 vnivanch Exp $
27// GEANT4 tag $Name: emstand-V09-03-24 $
28//
29// -------------------------------------------------------------------
30//
31// GEANT4 Class header file
32//
33//
34// File name:     G4PairProductionRelModel
35//
36// Author:        Andreas Schaelicke
37//
38// Creation date: 02.04.2009
39//
40// Modifications:
41//
42// Class Description:
43//
44// Implementation of gamma convertion to e+e- in the field of a nucleus
45// relativistic approximation
46//
47
48// -------------------------------------------------------------------
49//
50
51#ifndef G4PairProductionRelModel_h
52#define G4PairProductionRelModel_h 1
53
54#include "G4VEmModel.hh"
55#include "G4PhysicsTable.hh"
56#include "G4NistManager.hh"
57
58class G4ParticleChangeForGamma;
59
60class G4PairProductionRelModel : public G4VEmModel
61{
62
63public:
64
65  G4PairProductionRelModel(const G4ParticleDefinition* p = 0, 
66                           const G4String& nam = "BetheHeitlerLPM");
67 
68  virtual ~G4PairProductionRelModel();
69
70  virtual void Initialise(const G4ParticleDefinition*, const G4DataVector&);
71
72  virtual G4double ComputeCrossSectionPerAtom(
73                                const G4ParticleDefinition*,
74                                      G4double kinEnergy, 
75                                      G4double Z, 
76                                      G4double A=0., 
77                                      G4double cut=0.,
78                                      G4double emax=DBL_MAX);
79
80  virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
81                                 const G4MaterialCutsCouple*,
82                                 const G4DynamicParticle*,
83                                 G4double tmin,
84                                 G4double maxEnergy);
85
86  virtual void SetupForMaterial(const G4ParticleDefinition*,
87                                const G4Material*,G4double);
88
89  // * fast inline functions *
90  inline void SetCurrentElement(G4double /*Z*/);
91
92  // set / get methods
93  inline void SetLPMconstant(G4double val);
94  inline G4double LPMconstant() const;
95
96  inline void SetLPMflag(G4bool);
97  inline G4bool LPMflag() const;
98
99protected:
100  // screening functions
101  inline G4double Phi1(G4double delta) const;
102  inline G4double Phi2(G4double delta) const;
103  inline G4double DeltaMax() const;
104  inline G4double DeltaMin(G4double) const;
105  // lpm functions
106  void  CalcLPMFunctions(G4double k, G4double eplus);
107
108  // obsolete
109  G4double ScreenFunction1(G4double ScreenVariable);
110  G4double ScreenFunction2(G4double ScreenVariable);
111
112  G4double ComputeXSectionPerAtom(G4double totalEnergy, G4double Z);
113
114  G4double ComputeDXSectionPerAtom(G4double eplusEnergy, G4double totalEnergy, G4double Z);
115  G4double ComputeRelDXSectionPerAtom(G4double eplusEnergy, G4double totalEnergy, G4double Z);
116
117  // hide assignment operator
118  G4PairProductionRelModel & operator=(const G4PairProductionRelModel &right);
119  G4PairProductionRelModel(const  G4PairProductionRelModel&);
120
121  G4NistManager*              nist;
122
123  G4ParticleDefinition*     theGamma;
124  G4ParticleDefinition*     theElectron;
125  G4ParticleDefinition*     thePositron;
126  G4ParticleChangeForGamma* fParticleChange;
127
128  G4double fLPMconstant;
129  G4bool   fLPMflag;
130
131  // cash
132  G4double z13, z23, lnZ;
133  G4double Fel, Finel, fCoulomb; 
134  G4double currentZ;
135
136  // LPM effect
137  G4double lpmEnergy;
138  G4double xiLPM, phiLPM, gLPM;
139
140  // consts
141  G4bool   use_completescreening;
142
143  static const G4double xgi[8], wgi[8];
144  static const G4double Fel_light[5];
145  static const G4double Finel_light[5];
146  static const G4double facFel;
147  static const G4double facFinel;
148
149  static const G4double preS1, logTwo;
150
151};
152
153
154//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
155
156inline 
157void G4PairProductionRelModel::SetLPMconstant(G4double val) 
158{
159  fLPMconstant = val;
160}
161
162//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
163
164inline 
165G4double G4PairProductionRelModel::LPMconstant() const 
166{
167  return fLPMconstant;
168}
169
170//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
171
172inline 
173void G4PairProductionRelModel::SetLPMflag(G4bool val) 
174{
175  fLPMflag = val;
176}
177
178//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
179
180inline 
181G4bool G4PairProductionRelModel::LPMflag() const 
182{
183  return fLPMflag;
184}
185
186//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
187
188inline void G4PairProductionRelModel::SetCurrentElement(G4double Z)
189{
190  if(Z != currentZ) {
191    currentZ = Z;
192
193    G4int iz = G4int(Z);
194    z13 = nist->GetZ13(iz);
195    z23 = z13*z13;
196    lnZ = nist->GetLOGZ(iz);
197
198    if (iz <= 4) {
199      Fel = Fel_light[iz]; 
200      Finel = Finel_light[iz] ; 
201    }
202    else {
203      Fel = facFel - lnZ/3. ;
204      Finel = facFinel - 2.*lnZ/3. ;
205    }
206    fCoulomb=GetCurrentElement()->GetfCoulomb();
207  }
208}
209
210//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
211
212inline G4double G4PairProductionRelModel::Phi1(G4double delta) const
213{
214   G4double screenVal;
215
216   if (delta > 1.)
217     screenVal = 21.12 - 4.184*std::log(delta+0.952);
218   else
219     screenVal = 20.868 - delta*(3.242 - 0.625*delta);
220
221   return screenVal;
222}
223
224//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
225
226inline G4double G4PairProductionRelModel::Phi2(G4double delta) const
227{
228   G4double screenVal;
229
230   if (delta > 1.)
231     screenVal = 21.12 - 4.184*std::log(delta+0.952);
232   else
233     screenVal = 20.209 - delta*(1.930 + 0.086*delta);
234
235   return screenVal;
236}
237
238inline G4double G4PairProductionRelModel::ScreenFunction1(G4double ScreenVariable)
239
240// compute the value of the screening function 3*PHI1 - PHI2
241
242{
243   G4double screenVal;
244
245   if (ScreenVariable > 1.)
246     screenVal = 42.24 - 8.368*std::log(ScreenVariable+0.952);
247   else
248     screenVal = 42.392 - ScreenVariable*(7.796 - 1.961*ScreenVariable);
249
250   return screenVal;
251}
252
253//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
254
255inline G4double G4PairProductionRelModel::ScreenFunction2(G4double ScreenVariable)
256
257// compute the value of the screening function 1.5*PHI1 + 0.5*PHI2
258
259{
260   G4double screenVal;
261
262   if (ScreenVariable > 1.)
263     screenVal = 42.24 - 8.368*std::log(ScreenVariable+0.952);
264   else
265     screenVal = 41.405 - ScreenVariable*(5.828 - 0.8945*ScreenVariable);
266
267   return screenVal;
268}
269
270//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
271
272
273inline G4double G4PairProductionRelModel::DeltaMax() const
274{
275  // k > 50 MeV
276  G4double FZ = 8.*(lnZ/3. + fCoulomb);
277  return std::exp( (42.24-FZ)/8.368 ) + 0.952;
278}
279
280inline G4double G4PairProductionRelModel::DeltaMin(G4double k) const
281{
282  return 4.*136./z13*(electron_mass_c2/k);
283}
284
285
286
287#endif
Note: See TracBrowser for help on using the repository browser.