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