source: trunk/source/processes/electromagnetic/standard/include/G4PAIySection.hh@ 1036

Last change on this file since 1036 was 1007, checked in by garnier, 17 years ago

update to geant4.9.2

File size: 8.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//
27// $Id: G4PAIySection.hh,v 1.1 2007/10/01 17:45:14 vnivanch Exp $
28// GEANT4 tag $Name: geant4-09-02 $
29//
30//
31// G4PAIySection.hh -- header file
32//
33//
34// Preparation of ionizing collision cross section according to Photo Absorption
35// Ionization (PAI) model for simulation of ionization energy losses in very thin
36// absorbers. Author: Vladimir.Grichine@cern.ch
37//
38// History:
39//
40// 01.10.07, V.Ivanchenko create using V.Grichine G4PAIxSection class
41
42#ifndef G4PAIYSECTION_HH
43#define G4PAIYSECTION_HH
44
45#include "G4ios.hh"
46#include "globals.hh"
47#include "Randomize.hh"
48
49#include "G4SandiaTable.hh"
50
51class G4PAIySection
52{
53public:
54
55 G4PAIySection();
56
57 ~G4PAIySection();
58
59 void Initialize(const G4Material* material,
60 G4double maxEnergyTransfer,
61 G4double betaGammaSq);
62
63 void InitPAI() ;
64
65 void NormShift( G4double betaGammaSq ) ;
66
67 void SplainPAI( G4double betaGammaSq ) ;
68
69 // Physical methods
70 G4double RutherfordIntegral( G4int intervalNumber,
71 G4double limitLow,
72 G4double limitHigh ) ;
73
74 G4double ImPartDielectricConst( G4int intervalNumber,
75 G4double energy ) ;
76
77 G4double RePartDielectricConst(G4double energy) ;
78
79 G4double DifPAIySection( G4int intervalNumber,
80 G4double betaGammaSq ) ;
81
82 G4double PAIdNdxCerenkov( G4int intervalNumber,
83 G4double betaGammaSq ) ;
84
85 G4double PAIdNdxPlasmon( G4int intervalNumber,
86 G4double betaGammaSq ) ;
87
88 void IntegralPAIySection() ;
89 void IntegralCerenkov() ;
90 void IntegralPlasmon() ;
91
92 G4double SumOverInterval(G4int intervalNumber) ;
93 G4double SumOverIntervaldEdx(G4int intervalNumber) ;
94 G4double SumOverInterCerenkov(G4int intervalNumber) ;
95 G4double SumOverInterPlasmon(G4int intervalNumber) ;
96
97 G4double SumOverBorder( G4int intervalNumber,
98 G4double energy ) ;
99 G4double SumOverBorderdEdx( G4int intervalNumber,
100 G4double energy ) ;
101 G4double SumOverBordCerenkov( G4int intervalNumber,
102 G4double energy ) ;
103 G4double SumOverBordPlasmon( G4int intervalNumber,
104 G4double energy ) ;
105
106 G4double GetStepEnergyLoss( G4double step ) ;
107 G4double GetStepCerenkovLoss( G4double step ) ;
108 G4double GetStepPlasmonLoss( G4double step ) ;
109
110 // Inline access functions
111
112 G4int GetNumberOfGammas() const { return fNumberOfGammas ; }
113
114 G4int GetSplineSize() const { return fSplineNumber ; }
115
116 G4int GetIntervalNumber() const { return fIntervalNumber ; }
117
118 G4double GetEnergyInterval(G4int i){ return fEnergyInterval[i] ; }
119
120 G4double GetDifPAIySection(G4int i){ return fDifPAIySection[i] ; }
121 G4double GetPAIdNdxCrenkov(G4int i){ return fdNdxCerenkov[i] ; }
122 G4double GetPAIdNdxPlasmon(G4int i){ return fdNdxPlasmon[i] ; }
123
124 G4double GetMeanEnergyLoss() const {return fIntegralPAIySection[0] ; }
125 G4double GetMeanCerenkovLoss() const {return fIntegralCerenkov[0] ; }
126 G4double GetMeanPlasmonLoss() const {return fIntegralPlasmon[0] ; }
127
128 G4double GetNormalizationCof() const { return fNormalizationCof ; }
129
130 inline G4double GetPAItable(G4int i,G4int j) const ;
131
132 inline G4double GetLorentzFactor(G4int i) const ;
133
134 inline G4double GetSplineEnergy(G4int i) const ;
135
136 inline G4double GetIntegralPAIySection(G4int i) const ;
137 inline G4double GetIntegralPAIdEdx(G4int i) const ;
138 inline G4double GetIntegralCerenkov(G4int i) const ;
139 inline G4double GetIntegralPlasmon(G4int i) const ;
140
141private :
142
143// Local class constants
144
145 static const G4double fDelta ; // energy shift from interval border = 0.001
146 static const G4double fError ; // error in lin-log approximation = 0.005
147
148 static G4int fNumberOfGammas ; // = 111 ;
149 static const G4double fLorentzFactor[112] ; // static gamma array
150
151 static
152 const G4int fRefGammaNumber ; // The number of gamma for creation of spline (15)
153
154 G4int fIntervalNumber ; // The number of energy intervals
155 G4double fNormalizationCof ; // Normalization cof for PhotoAbsorptionXsection
156
157 G4double fDensity ; // Current density
158 G4double fElectronDensity ; // Current electron (number) density
159 G4int fSplineNumber ; // Current size of spline
160
161 G4SandiaTable* fSandia;
162
163 G4double fEnergyInterval[500] ;
164 G4double fA1[500] ;
165 G4double fA2[500] ;
166 G4double fA3[500] ;
167 G4double fA4[500] ;
168
169 static
170 const G4int fMaxSplineSize ; // Max size of output splain arrays = 500
171
172 G4double fSplineEnergy[500] ; // energy points of splain
173 G4double fRePartDielectricConst[500] ; // Real part of dielectric const
174 G4double fImPartDielectricConst[500] ; // Imaginary part of dielectric const
175 G4double fIntegralTerm[500] ; // Integral term in PAI cross section
176 G4double fDifPAIySection[500] ; // Differential PAI cross section
177 G4double fdNdxCerenkov[500] ; // dNdx of Cerenkov collisions
178 G4double fdNdxPlasmon[500] ; // dNdx of Plasmon collisions
179
180 G4double fIntegralPAIySection[500] ; // Integral PAI cross section ?
181 G4double fIntegralPAIdEdx[500] ; // Integral PAI dEdx ?
182 G4double fIntegralCerenkov[500] ; // Integral Cerenkov N>omega ?
183 G4double fIntegralPlasmon[500] ; // Integral Plasmon N>omega ?
184
185 G4double fPAItable[500][112] ; // Output array
186
187} ;
188
189//////////////// Inline methods //////////////////////////////////
190//
191
192
193inline G4double G4PAIySection::GetPAItable(G4int i, G4int j) const
194{
195 return fPAItable[i][j] ;
196}
197
198inline G4double G4PAIySection::GetLorentzFactor(G4int j) const
199{
200 return fLorentzFactor[j] ;
201}
202
203inline G4double G4PAIySection::GetSplineEnergy(G4int i) const
204{
205 if(i < 1 || i > fSplineNumber)
206 {
207 G4Exception("Invalid argument in G4PAIySection::GetSplineEnergy");
208 }
209 return fSplineEnergy[i] ;
210}
211
212inline G4double G4PAIySection::GetIntegralPAIySection(G4int i) const
213{
214 if(i < 1 || i > fSplineNumber)
215 {
216 G4Exception("Invalid argument in G4PAIySection::GetIntegralPAIySection");
217 }
218 return fIntegralPAIySection[i] ;
219}
220
221inline G4double G4PAIySection::GetIntegralPAIdEdx(G4int i) const
222{
223 if(i < 1 || i > fSplineNumber)
224 {
225 G4Exception("Invalid argument in G4PAIySection::GetIntegralPAIySection");
226 }
227 return fIntegralPAIdEdx[i] ;
228}
229
230inline G4double G4PAIySection::GetIntegralCerenkov(G4int i) const
231{
232 if(i < 1 || i > fSplineNumber)
233 {
234 G4Exception("Invalid argument in G4PAIySection::GetIntegralCerenkov");
235 }
236 return fIntegralCerenkov[i] ;
237}
238
239inline G4double G4PAIySection::GetIntegralPlasmon(G4int i) const
240{
241 if(i < 1 || i > fSplineNumber)
242 {
243 G4Exception("Invalid argument in G4PAIySection::GetIntegralPlasmon");
244 }
245 return fIntegralPlasmon[i] ;
246}
247
248#endif
249
250// ----------------- end of G4PAIySection header file -------------------
Note: See TracBrowser for help on using the repository browser.