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

Last change on this file since 1350 was 1337, checked in by garnier, 14 years ago

tag geant4.9.4 beta 1 + modifs locales

File size: 6.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: G4eBremsstrahlungModel.hh,v 1.26 2009/02/20 12:06:37 vnivanch Exp $
27// GEANT4 tag $Name: geant4-09-04-beta-01 $
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 = "eBrem");
72
73  virtual ~G4eBremsstrahlungModel();
74
75  virtual void Initialise(const G4ParticleDefinition*, const G4DataVector&);
76
77  virtual 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
103protected:
104
105  const G4Element* SelectRandomAtom(const G4MaterialCutsCouple* couple);
106
107private:
108
109  void SetParticle(const G4ParticleDefinition* p);
110
111  G4double ComputeBremLoss(G4double Z, G4double tkin, G4double cut);
112
113  G4double PositronCorrFactorLoss(G4double Z, G4double tkin, G4double cut);
114
115  G4double PositronCorrFactorSigma(G4double Z, G4double tkin, G4double cut);
116
117  G4DataVector* ComputePartialSumSigma(const G4Material* material,
118                                             G4double tkin, G4double cut);
119
120  G4double SupressionFunction(const G4Material* material, G4double tkin,
121                                    G4double gammaEnergy);
122
123  inline G4double ScreenFunction1(G4double ScreenVariable);
124
125  inline G4double ScreenFunction2(G4double ScreenVariable);
126
127  // hide assignment operator
128  G4eBremsstrahlungModel & operator=(const  G4eBremsstrahlungModel &right);
129  G4eBremsstrahlungModel(const  G4eBremsstrahlungModel&);
130
131protected:
132
133  const G4ParticleDefinition* particle;
134  G4ParticleDefinition*       theGamma;
135  G4ParticleChangeForLoss*    fParticleChange;
136
137  G4double minThreshold;
138  G4bool   isElectron;
139
140private:
141
142  G4double highKinEnergy;
143  G4double lowKinEnergy;
144  G4double probsup;
145  G4double MigdalConstant;
146  G4double LPMconstant;
147  G4bool   isInitialised;
148
149  std::vector<G4DataVector*> partialSumSigma;
150
151};
152
153//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
154
155inline G4double G4eBremsstrahlungModel::ScreenFunction1(G4double ScreenVariable)
156
157// compute the value of the screening function 3*PHI1 - PHI2
158
159{
160  G4double screenVal;
161
162  if (ScreenVariable > 1.)
163    screenVal = 42.24 - 8.368*std::log(ScreenVariable+0.952);
164  else
165    screenVal = 42.392 - ScreenVariable* (7.796 - 1.961*ScreenVariable);
166
167  return screenVal;
168} 
169
170//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
171
172inline
173G4double G4eBremsstrahlungModel::ScreenFunction2(G4double ScreenVariable)
174
175// compute the value of the screening function 1.5*PHI1 - 0.5*PHI2
176
177{
178  G4double screenVal;
179
180  if (ScreenVariable > 1.)
181    screenVal = 42.24 - 8.368*std::log(ScreenVariable+0.952);
182  else
183    screenVal = 41.734 - ScreenVariable* (6.484 - 1.250*ScreenVariable);
184
185  return screenVal;
186} 
187
188//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
189
190#endif
Note: See TracBrowser for help on using the repository browser.