Geant4-11
G4PAIxSection.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// G4PAIxSection.hh -- header file
27//
28// GEANT 4 class header file --- Copyright CERN 1995
29// CERB Geneva Switzerland
30//
31// for information related to this code, please, contact
32// CERN, CN Division, ASD Group
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// 28.10.11, V. Ivanchenko: Migration of exceptions to the new design
41// 19.10.03, V. Grichine: Integral dEdx was added for G4PAIModel class
42// 13.05.03, V. Grichine: Numerical instability was fixed in SumOverInterval/Border
43// functions
44// 10.02.02, V. Grichine: New functions and arrays/gets for Cerenkov and
45// plasmon collisions dN/dx
46// 27.10.99, V. Grichine: Bug fixed in constructors, 3rd constructor and
47// GetStepEnergyLoss(step) were added, fDelta = 0.005
48// 30.11.97, V. Grichine: 2nd version
49// 11.06.97, V. Grichine: 1st version
50
51#ifndef G4PAIXSECTION_HH
52#define G4PAIXSECTION_HH
53
54#include "G4ios.hh"
55#include "globals.hh"
56#include "Randomize.hh"
57
58#include"G4SandiaTable.hh"
59
61class G4Sandiatable;
62
63
65{
66public:
67 // Constructors
70
71 G4PAIxSection( G4int materialIndex,
72 G4double maxEnergyTransfer );
73
74 G4PAIxSection( G4int materialIndex, // for proton loss table
75 G4double maxEnergyTransfer,
76 G4double betaGammaSq ,
77 G4double** photoAbsCof, G4int intNumber );
78
79 G4PAIxSection( G4int materialIndex, // test constructor
80 G4double maxEnergyTransfer,
81 G4double betaGammaSq );
82
84
85 void Initialize(const G4Material* material, G4double maxEnergyTransfer,
86 G4double betaGammaSq, G4SandiaTable*);
87
88 // General control functions
89
92
93 void InitPAI();
94
95 void NormShift( G4double betaGammaSq );
96
97 void SplainPAI( G4double betaGammaSq );
98
99 // Physical methods
100
101 G4double RutherfordIntegral( G4int intervalNumber,
102 G4double limitLow,
103 G4double limitHigh );
104
105 G4double ImPartDielectricConst( G4int intervalNumber,
107
110
112
113 G4double DifPAIxSection( G4int intervalNumber,
114 G4double betaGammaSq );
115
116 G4double PAIdNdxCerenkov( G4int intervalNumber,
117 G4double betaGammaSq );
118 G4double PAIdNdxMM( G4int intervalNumber,
119 G4double betaGammaSq );
120
121 G4double PAIdNdxPlasmon( G4int intervalNumber,
122 G4double betaGammaSq );
123
124 G4double PAIdNdxResonance( G4int intervalNumber,
125 G4double betaGammaSq );
126
127
128 void IntegralPAIxSection();
129 void IntegralCerenkov();
130 void IntegralMM();
131 void IntegralPlasmon();
132 void IntegralResonance();
133
134 G4double SumOverInterval(G4int intervalNumber);
135 G4double SumOverIntervaldEdx(G4int intervalNumber);
136 G4double SumOverInterCerenkov(G4int intervalNumber);
137 G4double SumOverInterMM(G4int intervalNumber);
138 G4double SumOverInterPlasmon(G4int intervalNumber);
139 G4double SumOverInterResonance(G4int intervalNumber);
140
141 G4double SumOverBorder( G4int intervalNumber,
143 G4double SumOverBorderdEdx( G4int intervalNumber,
145 G4double SumOverBordCerenkov( G4int intervalNumber,
147 G4double SumOverBordMM( G4int intervalNumber,
149 G4double SumOverBordPlasmon( G4int intervalNumber,
151 G4double SumOverBordResonance( G4int intervalNumber,
153
159
166
167 // Inline access functions
168
170
172
174
176
182
185 G4double GetMeanMMLoss() const {return fIntegralMM[0]; }
188
190
192
194
195 inline void SetVerbose(G4int v) { fVerbose=v; };
196
197
198 inline G4double GetPAItable(G4int i,G4int j) const;
199
200 inline G4double GetSplineEnergy(G4int i) const;
201
202 inline G4double GetIntegralPAIxSection(G4int i) const;
203 inline G4double GetIntegralPAIdEdx(G4int i) const;
204 inline G4double GetIntegralCerenkov(G4int i) const;
205 inline G4double GetIntegralMM(G4int i) const;
206 inline G4double GetIntegralPlasmon(G4int i) const;
207 inline G4double GetIntegralResonance(G4int i) const;
208
209 G4PAIxSection & operator=(const G4PAIxSection &right) = delete;
210 G4PAIxSection(const G4PAIxSection&) = delete;
211
212private :
213
214 void CallError(G4int i, const G4String& methodName) const;
215
216 // Local class constants
217
218 static const G4double fDelta; // energy shift from interval border = 0.001
219 static const G4double fError; // error in lin-log approximation = 0.005
220
221 static G4int fNumberOfGammas; // = 111;
222 static const G4double fLorentzFactor[112]; // static gamma array
223
224 static
225 const G4int fRefGammaNumber; //The number of gamma for creation of spline (15)
226
227 G4int fIntervalNumber ; // The number of energy intervals
228 G4double fNormalizationCof; // Normalization cof for PhotoAbsorptionXsection
229
230 // G4double fBetaGammaSq; // (beta*gamma)^2
231
232 G4int fMaterialIndex; // current material index
233 G4double fDensity; // Current density
234 G4double fElectronDensity; // Current electron (number) density
235 G4double fLowEnergyCof; // Correction cof for low energy region
236 G4int fSplineNumber; // Current size of spline
237 G4int fVerbose; // verbose flag
238
239 // Arrays of Sandia coefficients
240
242
244
250
251 static
252 const G4int fMaxSplineSize ; // Max size of output splain arrays = 500
253
254 G4DataVector fSplineEnergy; // energy points of splain
255 G4DataVector fRePartDielectricConst; // Real part of dielectric const
256 G4DataVector fImPartDielectricConst; // Imaginary part of dielectric const
257 G4DataVector fIntegralTerm; // Integral term in PAI cross section
258 G4DataVector fDifPAIxSection; // Differential PAI cross section
259 G4DataVector fdNdxCerenkov; // dNdx of Cerenkov collisions
260 G4DataVector fdNdxPlasmon; // dNdx of Plasmon collisions
261 G4DataVector fdNdxMM; // dNdx of MM-Cerenkov collisions
262 G4DataVector fdNdxResonance; // dNdx of Resonance collisions
263
264 G4DataVector fIntegralPAIxSection; // Integral PAI cross section ?
265 G4DataVector fIntegralPAIdEdx; // Integral PAI dEdx ?
266 G4DataVector fIntegralCerenkov; // Integral Cerenkov N>omega ?
267 G4DataVector fIntegralPlasmon; // Integral Plasmon N>omega ?
268 G4DataVector fIntegralMM; // Integral MM N>omega ?
269 G4DataVector fIntegralResonance; // Integral resonance N>omega ?
270
271 G4double fPAItable[500][112]; // Output array
272
273};
274
276//
277
279{
280 return fPAItable[i][j];
281}
282
284{
285 if(i < 1 || i > fSplineNumber) { CallError(i, "GetSplineEnergy"); }
286 return fSplineEnergy[i];
287}
288
290{
291 if(i < 1 || i > fSplineNumber) { CallError(i, "GetIntegralPAIxSection"); }
292 return fIntegralPAIxSection[i];
293}
294
296{
297 if(i < 1 || i > fSplineNumber) { CallError(i, "GetIntegralPAIdEdx"); }
298 return fIntegralPAIdEdx[i];
299}
300
302{
303 if(i < 1 || i > fSplineNumber) { CallError(i, "GetIntegralCerenkov"); }
304 return fIntegralCerenkov[i];
305}
306
308{
309 if(i < 1 || i > fSplineNumber) { CallError(i, "GetIntegralMM"); }
310 return fIntegralMM[i];
311}
312
314{
315 if(i < 1 || i > fSplineNumber) { CallError(i, "GetIntegralPlasmon"); }
316 return fIntegralPlasmon[i];
317}
318
320{
321 if(i < 1 || i > fSplineNumber) { CallError(i, "GetIntegralResonance"); }
322 return fIntegralResonance[i];
323}
324
325#endif
326
327// ----------------- end of G4PAIxSection header file -------------------
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
G4double GetIntegralCerenkov(G4int i) const
G4int GetSplineSize() const
G4double RutherfordIntegral(G4int intervalNumber, G4double limitLow, G4double limitHigh)
G4double ImPartDielectricConst(G4int intervalNumber, G4double energy)
G4double GetPlasmonEnergyTransfer()
G4double SumOverIntervaldEdx(G4int intervalNumber)
void IntegralPlasmon()
G4DataVector fImPartDielectricConst
G4DataVector fRePartDielectricConst
G4double SumOverInterCerenkov(G4int intervalNumber)
G4DataVector fdNdxPlasmon
G4double GetStepMMLoss(G4double step)
G4double GetIntegralPAIdEdx(G4int i) const
G4double GetDifPAIxSection(G4int i)
G4double GetMeanResonanceLoss() const
G4double SumOverBordCerenkov(G4int intervalNumber, G4double energy)
G4double fPAItable[500][112]
G4double GetStepPlasmonLoss(G4double step)
G4double GetRutherfordEnergyTransfer()
G4int GetNumberOfGammas() const
G4double fElectronDensity
G4double SumOverBordResonance(G4int intervalNumber, G4double energy)
G4DataVector fIntegralPAIdEdx
G4DataVector fA3
G4double GetPAIdNdxPlasmon(G4int i)
G4DataVector fIntegralMM
void ComputeLowEnergyCof()
G4DataVector fIntegralPlasmon
G4double SumOverInterMM(G4int intervalNumber)
G4double GetEnergyInterval(G4int i)
G4double RePartDielectricConst(G4double energy)
static const G4double fDelta
G4DataVector fdNdxCerenkov
G4DataVector fDifPAIxSection
G4double GetElectronRange(G4double energy)
G4double GetIntegralPlasmon(G4int i) const
static const G4double fError
static const G4int fMaxSplineSize
G4double GetPAItable(G4int i, G4int j) const
G4double SumOverInterPlasmon(G4int intervalNumber)
G4double GetMMEnergyTransfer()
G4double PAIdNdxMM(G4int intervalNumber, G4double betaGammaSq)
G4double GetPAIdNdxCerenkov(G4int i)
G4DataVector fdNdxMM
G4double SumOverBorderdEdx(G4int intervalNumber, G4double energy)
static const G4double fLorentzFactor[112]
G4double SumOverBordMM(G4int intervalNumber, G4double energy)
G4double GetLowEnergyCof() const
G4double DifPAIxSection(G4int intervalNumber, G4double betaGammaSq)
void SetVerbose(G4int v)
G4DataVector fEnergyInterval
G4double GetPAIdNdxMM(G4int i)
G4double GetMeanPlasmonLoss() const
G4double GetPAIdNdxResonance(G4int i)
G4double PAIdNdxCerenkov(G4int intervalNumber, G4double betaGammaSq)
G4DataVector fA2
static const G4int fRefGammaNumber
G4double GetStepResonanceLoss(G4double step)
static G4int fNumberOfGammas
G4double GetResonanceEnergyTransfer()
G4double PAIdNdxResonance(G4int intervalNumber, G4double betaGammaSq)
G4PAIxSection & operator=(const G4PAIxSection &right)=delete
void SplainPAI(G4double betaGammaSq)
G4double GetEnergyTransfer()
G4OrderedTable * fMatSandiaMatrix
G4double GetCerenkovEnergyTransfer()
G4double GetMeanCerenkovLoss() const
G4double GetSplineEnergy(G4int i) const
G4double GetIntegralMM(G4int i) const
G4double fDensity
void NormShift(G4double betaGammaSq)
G4PAIxSection(const G4PAIxSection &)=delete
G4int GetIntervalNumber() const
G4double PAIdNdxPlasmon(G4int intervalNumber, G4double betaGammaSq)
G4DataVector fIntegralResonance
G4DataVector fIntegralCerenkov
G4double GetPhotonRange(G4double energy)
void IntegralPAIxSection()
G4double GetNormalizationCof() const
G4double SumOverInterval(G4int intervalNumber)
G4double SumOverBorder(G4int intervalNumber, G4double energy)
G4double fNormalizationCof
G4double GetStepCerenkovLoss(G4double step)
G4double GetIntegralPAIxSection(G4int i) const
G4SandiaTable * fSandia
G4DataVector fdNdxResonance
void CallError(G4int i, const G4String &methodName) const
G4double GetIntegralResonance(G4int i) const
G4double GetMeanEnergyLoss() const
G4DataVector fSplineEnergy
G4double GetStepEnergyLoss(G4double step)
G4double GetMeanMMLoss() const
G4DataVector fA1
void IntegralCerenkov()
G4double fLowEnergyCof
G4double SumOverInterResonance(G4int intervalNumber)
void IntegralResonance()
G4DataVector fIntegralPAIxSection
G4DataVector fA4
G4DataVector fIntegralTerm
G4double SumOverBordPlasmon(G4int intervalNumber, G4double energy)
G4double GetLorentzFactor(G4int i) const
void Initialize(const G4Material *material, G4double maxEnergyTransfer, G4double betaGammaSq, G4SandiaTable *)
G4double energy(const ThreeVector &p, const G4double m)
string material
Definition: eplot.py:19