Geant4  10.00.p02
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: G4PAIySection.hh 70170 2013-05-24 13:32:33Z gcosmo $
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, G4double betaGammaSq, G4SandiaTable*);
62 
63  void ComputeLowEnergyCof(const G4Material* material);
64 
65  void InitPAI();
66 
67  void NormShift( G4double betaGammaSq );
68 
69  void SplainPAI( G4double betaGammaSq );
70 
71  // Physical methods
72  G4double RutherfordIntegral( G4int intervalNumber,
73  G4double limitLow,
74  G4double limitHigh );
75 
76  G4double ImPartDielectricConst( G4int intervalNumber,
77  G4double energy );
78 
80 
81  G4double DifPAIySection( G4int intervalNumber,
82  G4double betaGammaSq );
83 
84  G4double PAIdNdxCerenkov( G4int intervalNumber,
85  G4double betaGammaSq );
86 
87  G4double PAIdNdxPlasmon( G4int intervalNumber,
88  G4double betaGammaSq );
89 
90  void IntegralPAIySection();
91  void IntegralCerenkov();
92  void IntegralPlasmon();
93 
94  G4double SumOverInterval(G4int intervalNumber);
95  G4double SumOverIntervaldEdx(G4int intervalNumber);
96  G4double SumOverInterCerenkov(G4int intervalNumber);
97  G4double SumOverInterPlasmon(G4int intervalNumber);
98 
99  G4double SumOverBorder( G4int intervalNumber,
100  G4double energy );
101  G4double SumOverBorderdEdx( G4int intervalNumber,
102  G4double energy );
103  G4double SumOverBordCerenkov( G4int intervalNumber,
104  G4double energy );
105  G4double SumOverBordPlasmon( G4int intervalNumber,
106  G4double energy );
107 
111 
112  // Inline access functions
113 
114  inline G4int GetNumberOfGammas() const { return fNumberOfGammas; }
115 
116  inline G4int GetSplineSize() const { return fSplineNumber; }
117 
118  inline G4int GetIntervalNumber() const { return fIntervalNumber; }
119 
121 
124  inline G4double GetPAIdNdxPlasmon(G4int i){ return fdNdxPlasmon[i]; }
125 
126  inline G4double GetMeanEnergyLoss() const {return fIntegralPAIySection[0]; }
127  inline G4double GetMeanCerenkovLoss() const {return fIntegralCerenkov[0]; }
128  inline G4double GetMeanPlasmonLoss() const {return fIntegralPlasmon[0]; }
129 
130  inline G4double GetNormalizationCof() const { return fNormalizationCof; }
131 
132  inline G4double GetPAItable(G4int i,G4int j) const;
133 
134  inline G4double GetLorentzFactor(G4int i) const;
135 
136  inline G4double GetSplineEnergy(G4int i) const;
137 
138  inline G4double GetIntegralPAIySection(G4int i) const;
139  inline G4double GetIntegralPAIdEdx(G4int i) const;
140  inline G4double GetIntegralCerenkov(G4int i) const;
141  inline G4double GetIntegralPlasmon(G4int i) const;
142 
143  inline void SetVerbose(G4int v) { fVerbose = v; };
144 
145 private :
146 
147  void CallError(G4int i, const G4String& methodName) const;
148 
149  // Local class constants
150 
151  static const G4double fDelta; // energy shift from interval border = 0.001
152  static const G4double fError; // error in lin-log approximation = 0.005
153 
154  static G4int fNumberOfGammas; // = 111;
155  static const G4double fLorentzFactor[112]; // static gamma array
156 
157  static
158  const G4int fRefGammaNumber ; // The number of gamma for creation of spline (15)
159 
160  G4int fIntervalNumber ; // The number of energy intervals
161  G4double fNormalizationCof; // Normalization cof for PhotoAbsorptionXsection
162 
163  G4double fDensity; // Current density
164  G4double fElectronDensity; // Current electron (number) density
165  G4double fLowEnergyCof; // Correction cof for low energy region
166  G4int fSplineNumber; // Current size of spline
167  G4int fVerbose; // verbose flag
168 
170 
176 
177  static
178  const G4int fMaxSplineSize ; // Max size of output splain arrays = 500
179 
180  G4DataVector fSplineEnergy; // energy points of splain
181  G4DataVector fRePartDielectricConst; // Real part of dielectric const
182  G4DataVector fImPartDielectricConst; // Imaginary part of dielectric const
183  G4DataVector fIntegralTerm; // Integral term in PAI cross section
184  G4DataVector fDifPAIySection; // Differential PAI cross section
185  G4DataVector fdNdxCerenkov; // dNdx of Cerenkov collisions
186  G4DataVector fdNdxPlasmon; // dNdx of Plasmon collisions
187 
188  G4DataVector fIntegralPAIySection; // Integral PAI cross section ?
189  G4DataVector fIntegralPAIdEdx; // Integral PAI dEdx ?
190  G4DataVector fIntegralCerenkov; // Integral Cerenkov N>omega ?
191  G4DataVector fIntegralPlasmon; // Integral Plasmon N>omega ?
192 
193  G4double fPAItable[500][112]; // Output array
194 
195 };
196 
198 //
199 
200 
202 {
203  return fPAItable[i][j];
204 }
205 
207 {
208  return fLorentzFactor[j];
209 }
210 
212 {
213  if(i < 1 || i > fSplineNumber) { CallError(i, "GetSplineEnergy"); }
214  return fSplineEnergy[i];
215 }
216 
218 {
219  if(i < 1 || i > fSplineNumber) { CallError(i, "GetIntegralPAIySection"); }
220  return fIntegralPAIySection[i];
221 }
222 
224 {
225  if(i < 1 || i > fSplineNumber) { CallError(i, "GetIntegralPAIdEdx"); }
226  return fIntegralPAIdEdx[i];
227 }
228 
230 {
231  if(i < 1 || i > fSplineNumber) { CallError(i, "GetIntegralCerenkov"); }
232  return fIntegralCerenkov[i];
233 }
234 
236 {
237  if(i < 1 || i > fSplineNumber) { CallError(i, "GetIntegralPlasmon"); }
238  return fIntegralPlasmon[i];
239 }
240 
241 #endif
242 
243 // ----------------- end of G4PAIySection header file -------------------
G4DataVector fDifPAIySection
G4double SumOverIntervaldEdx(G4int intervalNumber)
G4DataVector fA3
G4double RutherfordIntegral(G4int intervalNumber, G4double limitLow, G4double limitHigh)
static const G4double fDelta
G4double GetDifPAIySection(G4int i)
G4double GetPAIdNdxCrenkov(G4int i)
static const G4int fMaxSplineSize
G4double GetNormalizationCof() const
void NormShift(G4double betaGammaSq)
static const G4double fLorentzFactor[112]
G4double GetMeanEnergyLoss() const
G4double ImPartDielectricConst(G4int intervalNumber, G4double energy)
G4DataVector fImPartDielectricConst
G4double SumOverInterCerenkov(G4int intervalNumber)
G4DataVector fRePartDielectricConst
G4double fDensity
G4double GetMeanPlasmonLoss() const
void CallError(G4int i, const G4String &methodName) const
G4double PAIdNdxCerenkov(G4int intervalNumber, G4double betaGammaSq)
void ComputeLowEnergyCof(const G4Material *material)
void IntegralCerenkov()
G4double fNormalizationCof
G4int GetNumberOfGammas() const
G4DataVector fIntegralPlasmon
G4double SumOverBorderdEdx(G4int intervalNumber, G4double energy)
int G4int
Definition: G4Types.hh:78
void IntegralPlasmon()
G4DataVector fIntegralTerm
G4double fPAItable[500][112]
G4double GetIntegralPAIdEdx(G4int i) const
G4double GetIntegralPlasmon(G4int i) const
G4DataVector fA2
G4DataVector fA1
static G4int fNumberOfGammas
G4double SumOverBordPlasmon(G4int intervalNumber, G4double energy)
G4double GetStepEnergyLoss(G4double step)
void Initialize(const G4Material *material, G4double maxEnergyTransfer, G4double betaGammaSq, G4SandiaTable *)
G4double fLowEnergyCof
G4double GetSplineEnergy(G4int i) const
G4DataVector fdNdxPlasmon
G4double GetIntegralPAIySection(G4int i) const
G4SandiaTable * fSandia
G4double GetIntegralCerenkov(G4int i) const
G4double GetStepCerenkovLoss(G4double step)
G4int GetIntervalNumber() const
static const G4int fRefGammaNumber
void SetVerbose(G4int v)
G4double SumOverBorder(G4int intervalNumber, G4double energy)
G4double energy(const ThreeVector &p, const G4double m)
G4double GetPAItable(G4int i, G4int j) const
void SplainPAI(G4double betaGammaSq)
static const G4double fError
G4DataVector fIntegralPAIySection
G4DataVector fIntegralPAIdEdx
G4DataVector fSplineEnergy
G4DataVector fEnergyInterval
G4double fElectronDensity
G4double SumOverBordCerenkov(G4int intervalNumber, G4double energy)
G4int GetSplineSize() const
double G4double
Definition: G4Types.hh:76
G4double RePartDielectricConst(G4double energy)
G4DataVector fdNdxCerenkov
G4DataVector fA4
G4DataVector fIntegralCerenkov
G4double GetStepPlasmonLoss(G4double step)
G4double GetEnergyInterval(G4int i)
G4double DifPAIySection(G4int intervalNumber, G4double betaGammaSq)
G4double SumOverInterval(G4int intervalNumber)
G4double PAIdNdxPlasmon(G4int intervalNumber, G4double betaGammaSq)
G4double GetMeanCerenkovLoss() const
void IntegralPAIySection()
G4double GetPAIdNdxPlasmon(G4int i)
G4double SumOverInterPlasmon(G4int intervalNumber)
G4double GetLorentzFactor(G4int i) const