source: trunk/source/processes/electromagnetic/standard/include/G4eBremsstrahlungModel.hh@ 830

Last change on this file since 830 was 819, checked in by garnier, 17 years ago

import all except CVS

File size: 7.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// $Id: G4eBremsstrahlungModel.hh,v 1.22 2007/05/23 08:47:34 vnivanch Exp $
27// GEANT4 tag $Name: $
28//
29// -------------------------------------------------------------------
30//
31// GEANT4 Class header file
32//
33//
34// File name: G4eBremsstrahlungModel
35//
36// Author: Vladimir Ivanchenko on base of Laszlo Urban code
37//
38// Creation date: 07.01.2002
39//
40// Modifications:
41//
42// 23-12-02 Change interface in order to move to cut per region (V.Ivanchenko)
43// 24-01-03 Make models region aware (V.Ivanchenko)
44// 13-02-03 Add name (V.Ivanchenko)
45// 08-04-05 Major optimisation of internal interfaces (V.Ivantchenko)
46// 07-02-06 public function ComputeCrossSectionPerAtom() (mma)
47//
48//
49// Class Description:
50//
51// Implementation of energy loss for gamma emission by electrons and
52// positrons
53
54// -------------------------------------------------------------------
55//
56
57#ifndef G4eBremsstrahlungModel_h
58#define G4eBremsstrahlungModel_h 1
59
60#include "G4VEmModel.hh"
61
62class G4Element;
63class G4ParticleChangeForLoss;
64
65class G4eBremsstrahlungModel : public G4VEmModel
66{
67
68public:
69
70 G4eBremsstrahlungModel(const G4ParticleDefinition* p = 0,
71 const G4String& nam = "StandBrem");
72
73 virtual ~G4eBremsstrahlungModel();
74
75 virtual void Initialise(const G4ParticleDefinition*, const G4DataVector&);
76
77 G4double MinEnergyCut(const G4ParticleDefinition*,
78 const G4MaterialCutsCouple*);
79
80 virtual G4double ComputeDEDXPerVolume(const G4Material*,
81 const G4ParticleDefinition*,
82 G4double kineticEnergy,
83 G4double cutEnergy);
84
85 virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition*,
86 G4double tkin,
87 G4double Z, G4double,
88 G4double cut,
89 G4double maxE = DBL_MAX);
90
91 virtual G4double CrossSectionPerVolume(const G4Material*,
92 const G4ParticleDefinition*,
93 G4double kineticEnergy,
94 G4double cutEnergy,
95 G4double maxEnergy);
96
97 virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
98 const G4MaterialCutsCouple*,
99 const G4DynamicParticle*,
100 G4double tmin,
101 G4double maxEnergy);
102
103 inline void SetLPMflag(G4bool val);
104 inline G4bool LPMflag() const;
105
106 inline void SetEnergyThreshold(G4double val);
107 inline G4double EnergyThreshold() const;
108
109protected:
110
111 inline G4double MaxSecondaryEnergy(const G4ParticleDefinition*,
112 G4double kineticEnergy);
113
114 const G4Element* SelectRandomAtom(const G4MaterialCutsCouple* couple);
115
116private:
117
118 void SetParticle(const G4ParticleDefinition* p);
119
120 G4double ComputeBremLoss(G4double Z, G4double tkin, G4double cut);
121
122 G4double PositronCorrFactorLoss(G4double Z, G4double tkin, G4double cut);
123
124 G4double PositronCorrFactorSigma(G4double Z, G4double tkin, G4double cut);
125
126 G4DataVector* ComputePartialSumSigma(const G4Material* material,
127 G4double tkin, G4double cut);
128
129 G4double SupressionFunction(const G4Material* material, G4double tkin,
130 G4double gammaEnergy);
131
132 inline G4double ScreenFunction1(G4double ScreenVariable);
133
134 inline G4double ScreenFunction2(G4double ScreenVariable);
135
136 // hide assignment operator
137 G4eBremsstrahlungModel & operator=(const G4eBremsstrahlungModel &right);
138 G4eBremsstrahlungModel(const G4eBremsstrahlungModel&);
139
140protected:
141
142 const G4ParticleDefinition* particle;
143 G4ParticleDefinition* theGamma;
144 G4ParticleChangeForLoss* fParticleChange;
145
146 G4double minThreshold;
147 G4bool isElectron;
148
149private:
150
151 G4double highKinEnergy;
152 G4double lowKinEnergy;
153 G4double probsup;
154 G4double MigdalConstant;
155 G4double LPMconstant;
156 G4double highEnergyTh;
157 G4bool theLPMflag;
158 G4bool isInitialised;
159
160 std::vector<G4DataVector*> partialSumSigma;
161
162};
163
164//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
165
166inline G4double G4eBremsstrahlungModel::ScreenFunction1(G4double ScreenVariable)
167
168// compute the value of the screening function 3*PHI1 - PHI2
169
170{
171 G4double screenVal;
172
173 if (ScreenVariable > 1.)
174 screenVal = 42.24 - 8.368*std::log(ScreenVariable+0.952);
175 else
176 screenVal = 42.392 - ScreenVariable* (7.796 - 1.961*ScreenVariable);
177
178 return screenVal;
179}
180
181//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
182
183inline
184G4double G4eBremsstrahlungModel::ScreenFunction2(G4double ScreenVariable)
185
186// compute the value of the screening function 1.5*PHI1 - 0.5*PHI2
187
188{
189 G4double screenVal;
190
191 if (ScreenVariable > 1.)
192 screenVal = 42.24 - 8.368*std::log(ScreenVariable+0.952);
193 else
194 screenVal = 41.734 - ScreenVariable* (6.484 - 1.250*ScreenVariable);
195
196 return screenVal;
197}
198
199//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
200
201inline
202G4double G4eBremsstrahlungModel::MaxSecondaryEnergy(
203 const G4ParticleDefinition*,
204 G4double kineticEnergy)
205{
206 return kineticEnergy;
207}
208
209//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
210
211inline
212void G4eBremsstrahlungModel::SetLPMflag(G4bool val)
213{
214 theLPMflag = val;
215}
216
217//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
218
219inline
220G4bool G4eBremsstrahlungModel::LPMflag() const
221{
222 return theLPMflag;
223}
224
225//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
226
227inline
228void G4eBremsstrahlungModel::SetEnergyThreshold(G4double val)
229{
230 highEnergyTh = val;
231}
232
233//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
234
235inline
236G4double G4eBremsstrahlungModel::EnergyThreshold() const
237{
238 return highEnergyTh;
239}
240
241//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
242
243#endif
Note: See TracBrowser for help on using the repository browser.