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