Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4eBremsstrahlungModel.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$
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class header file
31 //
32 //
33 // File name: G4eBremsstrahlungModel
34 //
35 // Author: Vladimir Ivanchenko on base of Laszlo Urban code
36 //
37 // Creation date: 07.01.2002
38 //
39 // Modifications:
40 //
41 // 23-12-02 Change interface in order to move to cut per region (V.Ivanchenko)
42 // 24-01-03 Make models region aware (V.Ivanchenko)
43 // 13-02-03 Add name (V.Ivanchenko)
44 // 08-04-05 Major optimisation of internal interfaces (V.Ivantchenko)
45 // 07-02-06 public function ComputeCrossSectionPerAtom() (mma)
46 //
47 //
48 // Class Description:
49 //
50 // Implementation of energy loss for gamma emission by electrons and
51 // positrons
52 
53 // -------------------------------------------------------------------
54 //
55 
56 #ifndef G4eBremsstrahlungModel_h
57 #define G4eBremsstrahlungModel_h 1
58 
59 #include "G4VEmModel.hh"
60 
61 class G4Element;
63 
65 {
66 
67 public:
68 
70  const G4String& nam = "eBrem");
71 
72  virtual ~G4eBremsstrahlungModel();
73 
74  virtual void Initialise(const G4ParticleDefinition*, const G4DataVector&);
75 
77  const G4ParticleDefinition*,
78  G4double kineticEnergy,
79  G4double cutEnergy);
80 
82  G4double tkin,
84  G4double cut,
85  G4double maxE = DBL_MAX);
86 
88  const G4ParticleDefinition*,
89  G4double kineticEnergy,
90  G4double cutEnergy,
91  G4double maxEnergy);
92 
93  virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
94  const G4MaterialCutsCouple*,
95  const G4DynamicParticle*,
96  G4double tmin,
97  G4double maxEnergy);
98 
99 protected:
100 
101  const G4Element* SelectRandomAtom(const G4MaterialCutsCouple* couple);
102 
103 private:
104 
105  void SetParticle(const G4ParticleDefinition* p);
106 
107  G4double ComputeBremLoss(G4double Z, G4double tkin, G4double cut);
108 
109  G4double PositronCorrFactorLoss(G4double Z, G4double tkin, G4double cut);
110 
111  G4double PositronCorrFactorSigma(G4double Z, G4double tkin, G4double cut);
112 
113  G4DataVector* ComputePartialSumSigma(const G4Material* material,
114  G4double tkin, G4double cut);
115 
116  G4double SupressionFunction(const G4Material* material, G4double tkin,
117  G4double gammaEnergy);
118 
119  inline G4double ScreenFunction1(G4double ScreenVariable);
120 
121  inline G4double ScreenFunction2(G4double ScreenVariable);
122 
123  // hide assignment operator
126 
127 protected:
128 
132 
134 
135 private:
136 
137  G4double highKinEnergy;
138  G4double lowKinEnergy;
139  G4double probsup;
140  G4double MigdalConstant;
141  G4double LPMconstant;
142  G4bool isInitialised;
143 
144  std::vector<G4DataVector*> partialSumSigma;
145 
146 };
147 
148 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
149 
150 inline G4double G4eBremsstrahlungModel::ScreenFunction1(G4double ScreenVariable)
151 
152 // compute the value of the screening function 3*PHI1 - PHI2
153 
154 {
155  G4double screenVal;
156 
157  if (ScreenVariable > 1.)
158  screenVal = 42.24 - 8.368*std::log(ScreenVariable+0.952);
159  else
160  screenVal = 42.392 - ScreenVariable* (7.796 - 1.961*ScreenVariable);
161 
162  return screenVal;
163 }
164 
165 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
166 
167 inline
168 G4double G4eBremsstrahlungModel::ScreenFunction2(G4double ScreenVariable)
169 
170 // compute the value of the screening function 1.5*PHI1 - 0.5*PHI2
171 
172 {
173  G4double screenVal;
174 
175  if (ScreenVariable > 1.)
176  screenVal = 42.24 - 8.368*std::log(ScreenVariable+0.952);
177  else
178  screenVal = 41.734 - ScreenVariable* (6.484 - 1.250*ScreenVariable);
179 
180  return screenVal;
181 }
182 
183 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
184 
185 #endif