Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4BetheHeitlerModel.hh
Go to the documentation of this file.
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: G4BetheHeitlerModel.hh 100399 2016-10-20 07:38:12Z gcosmo $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class header file
31 //
32 //
33 // File name: G4BetheHeitlerModel
34 //
35 // Author: Vladimir Ivanchenko on base of Michel Maire code
36 //
37 // Creation date: 19.04.2005
38 //
39 // Modifications:
40 // 02-02-06 Remove InitialiseCrossSectionPerAtom();
41 //
42 // Class Description:
43 //
44 // Implementation of gamma convertion to e+e- in the field of a nucleus
45 //
46 
47 // -------------------------------------------------------------------
48 //
49 
50 #ifndef G4BetheHeitlerModel_h
51 #define G4BetheHeitlerModel_h 1
52 
53 #include "G4VEmModel.hh"
54 #include "G4PhysicsTable.hh"
55 #include "G4Log.hh"
56 
58 class G4Pow;
59 
61 {
62 
63 public:
64 
65  explicit G4BetheHeitlerModel(const G4ParticleDefinition* p = 0,
66  const G4String& nam = "BetheHeitler");
67 
68  virtual ~G4BetheHeitlerModel();
69 
70  virtual void Initialise(const G4ParticleDefinition*,
71  const G4DataVector&) override;
72 
73  virtual void InitialiseLocal(const G4ParticleDefinition*,
74  G4VEmModel* masterModel) override;
75 
77  const G4ParticleDefinition*,
78  G4double kinEnergy,
79  G4double Z,
80  G4double A=0.,
81  G4double cut=0.,
82  G4double emax=DBL_MAX) override;
83 
84  virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
85  const G4MaterialCutsCouple*,
86  const G4DynamicParticle*,
87  G4double tmin,
88  G4double maxEnergy) override;
89 
90 private:
91 
92  G4double ScreenFunction1(G4double ScreenVariable);
93 
94  G4double ScreenFunction2(G4double ScreenVariable);
95 
96  // hide assignment operator
97  G4BetheHeitlerModel & operator=(const G4BetheHeitlerModel &right) = delete;
99 
100  G4Pow* g4calc;
101  G4ParticleDefinition* theGamma;
102  G4ParticleDefinition* theElectron;
103  G4ParticleDefinition* thePositron;
104  G4ParticleChangeForGamma* fParticleChange;
105 };
106 
107 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
108 
109 inline G4double G4BetheHeitlerModel::ScreenFunction1(G4double ScreenVariable)
110 
111 // compute the value of the screening function 3*PHI1 - PHI2
112 
113 {
114  G4double screenVal;
115 
116  if (ScreenVariable > 1.)
117  screenVal = 42.24 - 8.368*G4Log(ScreenVariable+0.952);
118  else
119  screenVal = 42.392 - ScreenVariable*(7.796 - 1.961*ScreenVariable);
120 
121  return screenVal;
122 }
123 
124 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
125 
126 inline G4double G4BetheHeitlerModel::ScreenFunction2(G4double ScreenVariable)
127 
128 // compute the value of the screening function 1.5*PHI1 - 0.5*PHI2
129 
130 {
131  G4double screenVal;
132 
133  if (ScreenVariable > 1.)
134  screenVal = 42.24 - 8.368*G4Log(ScreenVariable+0.952);
135  else
136  screenVal = 41.405 - ScreenVariable*(5.828 - 0.8945*ScreenVariable);
137 
138  return screenVal;
139 }
140 
141 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
142 
143 #endif
Definition: G4Pow.hh:56
const char * p
Definition: xmltok.h:285
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cut=0., G4double emax=DBL_MAX) override
virtual void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel) override
double A(double temperature)
G4BetheHeitlerModel(const G4ParticleDefinition *p=0, const G4String &nam="BetheHeitler")
static const G4double emax
G4double G4Log(G4double x)
Definition: G4Log.hh:230
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
double G4double
Definition: G4Types.hh:76
#define DBL_MAX
Definition: templates.hh:83