Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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//
27// $Id$
28//
29//
30// G4PAIxSection.hh -- header file
31//
32// GEANT 4 class header file --- Copyright CERN 1995
33// CERB Geneva Switzerland
34//
35// for information related to this code, please, contact
36// CERN, CN Division, ASD Group
37//
38// Preparation of ionizing collision cross section according to Photo Absorption
39// Ionization (PAI) model for simulation of ionization energy losses in very thin
40// absorbers. Author: [email protected]
41//
42// History:
43//
44// 28.10.11, V. Ivanchenko: Migration of exceptions to the new design
45// 19.10.03, V. Grichine: Integral dEdx was added for G4PAIModel class
46// 13.05.03, V. Grichine: Numerical instability was fixed in SumOverInterval/Border
47// functions
48// 10.02.02, V. Grichine: New functions and arrays/gets for Cerenkov and
49// plasmon collisions dN/dx
50// 27.10.99, V. Grichine: Bug fixed in constructors, 3rd constructor and
51// GetStepEnergyLoss(step) were added, fDelta = 0.005
52// 30.11.97, V. Grichine: 2nd version
53// 11.06.97, V. Grichine: 1st version
54
55#ifndef G4PAIXSECTION_HH
56#define G4PAIXSECTION_HH
57
58#include "G4ios.hh"
59#include "globals.hh"
60#include "Randomize.hh"
61
62#include"G4SandiaTable.hh"
63
65class G4Sandiatable;
66
67
69{
70public:
71 // Constructors
73
74 G4PAIxSection( G4int materialIndex,
75 G4double maxEnergyTransfer );
76
77 G4PAIxSection( G4int materialIndex, // for proton loss table
78 G4double maxEnergyTransfer,
79 G4double betaGammaSq ,
80 G4double** photoAbsCof, G4int intNumber );
81
82 G4PAIxSection( G4int materialIndex, // test constructor
83 G4double maxEnergyTransfer,
84 G4double betaGammaSq );
85
86 // G4PAIxSection(const G4PAIxSection& right);
87
88 // Destructor
89
91
92 // Operators
93 // G4PAIxSection& operator=(const G4PAIxSection& right);
94 // G4int operator==(const G4PAIxSection& right)const;
95 // G4int operator!=(const G4PAIxSection& right)const;
96
97 // Methods
98
99 // General control functions
100
101 void InitPAI();
102
103 void NormShift( G4double betaGammaSq );
104
105 void SplainPAI( G4double betaGammaSq );
106
107 // Physical methods
108
109
110 G4double RutherfordIntegral( G4int intervalNumber,
111 G4double limitLow,
112 G4double limitHigh );
113
114 G4double ImPartDielectricConst( G4int intervalNumber,
115 G4double energy );
116
119
121
122 G4double DifPAIxSection( G4int intervalNumber,
123 G4double betaGammaSq );
124
125 G4double PAIdNdxCerenkov( G4int intervalNumber,
126 G4double betaGammaSq );
127 G4double PAIdNdxMM( G4int intervalNumber,
128 G4double betaGammaSq );
129
130 G4double PAIdNdxPlasmon( G4int intervalNumber,
131 G4double betaGammaSq );
132
133 G4double PAIdNdxResonance( G4int intervalNumber,
134 G4double betaGammaSq );
135
136 void IntegralPAIxSection();
137 void IntegralCerenkov();
138 void IntegralMM();
139 void IntegralPlasmon();
140 void IntegralResonance();
141
142 G4double SumOverInterval(G4int intervalNumber);
143 G4double SumOverIntervaldEdx(G4int intervalNumber);
144 G4double SumOverInterCerenkov(G4int intervalNumber);
145 G4double SumOverInterMM(G4int intervalNumber);
146 G4double SumOverInterPlasmon(G4int intervalNumber);
147 G4double SumOverInterResonance(G4int intervalNumber);
148
149 G4double SumOverBorder( G4int intervalNumber,
150 G4double energy );
151 G4double SumOverBorderdEdx( G4int intervalNumber,
152 G4double energy );
153 G4double SumOverBordCerenkov( G4int intervalNumber,
154 G4double energy );
155 G4double SumOverBordMM( G4int intervalNumber,
156 G4double energy );
157 G4double SumOverBordPlasmon( G4int intervalNumber,
158 G4double energy );
159 G4double SumOverBordResonance( G4int intervalNumber,
160 G4double energy );
161
167
174
175 // Inline access functions
176
177 G4int GetNumberOfGammas() const { return fNumberOfGammas; }
178
179 G4int GetSplineSize() const { return fSplineNumber; }
180
181 G4int GetIntervalNumber() const { return fIntervalNumber; }
182
183 G4double GetEnergyInterval(G4int i){ return fEnergyInterval[i]; }
184
185 G4double GetDifPAIxSection(G4int i){ return fDifPAIxSection[i]; }
186 G4double GetPAIdNdxCerenkov(G4int i){ return fdNdxCerenkov[i]; }
187 G4double GetPAIdNdxMM(G4int i){ return fdNdxMM[i]; }
188 G4double GetPAIdNdxPlasmon(G4int i){ return fdNdxPlasmon[i]; }
189 G4double GetPAIdNdxResonance(G4int i){ return fdNdxResonance[i]; }
190
191 G4double GetMeanEnergyLoss() const {return fIntegralPAIxSection[0]; }
192 G4double GetMeanCerenkovLoss() const {return fIntegralCerenkov[0]; }
193 G4double GetMeanMMLoss() const {return fIntegralMM[0]; }
194 G4double GetMeanPlasmonLoss() const {return fIntegralPlasmon[0]; }
195 G4double GetMeanResonanceLoss() const {return fIntegralResonance[0]; }
196
197 G4double GetNormalizationCof() const { return fNormalizationCof; }
198
199 inline G4double GetPAItable(G4int i,G4int j) const;
200
201 inline G4double GetLorentzFactor(G4int i) const;
202
203 inline G4double GetSplineEnergy(G4int i) const;
204
205 inline G4double GetIntegralPAIxSection(G4int i) const;
206 inline G4double GetIntegralPAIdEdx(G4int i) const;
207 inline G4double GetIntegralCerenkov(G4int i) const;
208 inline G4double GetIntegralMM(G4int i) const;
209 inline G4double GetIntegralPlasmon(G4int i) const;
210 inline G4double GetIntegralResonance(G4int i) const;
211
212private :
213
214 void CallError(G4int i, const G4String& methodName) const;
215
216 G4PAIxSection & operator=(const G4PAIxSection &right);
218
219// Local class constants
220
221static const G4double fDelta; // energy shift from interval border = 0.001
222static const G4double fError; // error in lin-log approximation = 0.005
223
224static G4int fNumberOfGammas; // = 111;
225static const G4double fLorentzFactor[112]; // static gamma array
226
227static
228const G4int fRefGammaNumber ; // The number of gamma for creation of spline (15)
229
230G4int fIntervalNumber ; // The number of energy intervals
231G4double fNormalizationCof; // Normalization cof for PhotoAbsorptionXsection
232
233// G4double fBetaGammaSq; // (beta*gamma)^2
234
235 G4int fMaterialIndex; // current material index
236 G4double fDensity; // Current density
237 G4double fElectronDensity; // Current electron (number) density
238 G4int fSplineNumber; // Current size of spline
239
240// Arrays of Sandia coefficients
241
242 G4OrderedTable* fMatSandiaMatrix;
243 G4SandiaTable* fSandia;
244
245G4double* fEnergyInterval;
246G4double* fA1;
247G4double* fA2;
248G4double* fA3;
249G4double* fA4;
250
251static
252const G4int fMaxSplineSize ; // Max size of output splain arrays = 500
253
254/* ******************
255G4double* fSplineEnergy; // energy points of splain
256G4double* fRePartDielectricConst; // Real part of dielectric const
257G4double* fImPartDielectricConst; // Imaginary part of dielectric const
258G4double* fIntegralTerm; // Integral term in PAI cross section
259G4double* fDifPAIxSection; // Differential PAI cross section
260G4double* fIntegralPAIxSection; // Integral PAI cross section ?
261*/ ///////////////
262
263
264G4double fSplineEnergy[500]; // energy points of splain
265G4double fRePartDielectricConst[500]; // Real part of dielectric const
266G4double fImPartDielectricConst[500]; // Imaginary part of dielectric const
267G4double fIntegralTerm[500]; // Integral term in PAI cross section
268G4double fDifPAIxSection[500]; // Differential PAI cross section
269G4double fdNdxCerenkov[500]; // dNdx of Cerenkov collisions
270G4double fdNdxMM[500]; // dNdx of MM-Cerenkov collisions
271G4double fdNdxPlasmon[500]; // dNdx of Plasmon collisions
272G4double fdNdxResonance[500]; // dNdx of resonance collisions
273
274G4double fIntegralPAIxSection[500]; // Integral PAI cross section ?
275G4double fIntegralPAIdEdx[500]; // Integral PAI dEdx ?
276G4double fIntegralCerenkov[500]; // Integral Cerenkov N>omega ?
277G4double fIntegralMM[500]; // Integral MM-Cerenkov N>omega ?
278G4double fIntegralPlasmon[500]; // Integral Plasmon N>omega ?
279G4double fIntegralResonance[500]; // Integral resonance N>omega ?
280
281G4double fPAItable[500][112]; // Output array
282
283};
284
285//////////////// Inline methods //////////////////////////////////
286//
287
289{
290 return fPAItable[i][j];
291}
292
294{
295 return fLorentzFactor[j];
296}
297
299{
300 if(i < 1 || i > fSplineNumber) { CallError(i, "GetSplineEnergy"); }
301 return fSplineEnergy[i];
302}
303
305{
306 if(i < 1 || i > fSplineNumber) { CallError(i, "GetIntegralPAIxSection"); }
307 return fIntegralPAIxSection[i];
308}
309
311{
312 if(i < 1 || i > fSplineNumber) { CallError(i, "GetIntegralPAIdEdx"); }
313 return fIntegralPAIdEdx[i];
314}
315
317{
318 if(i < 1 || i > fSplineNumber) { CallError(i, "GetIntegralCerenkov"); }
319 return fIntegralCerenkov[i];
320}
321
323{
324 if(i < 1 || i > fSplineNumber) { CallError(i, "GetIntegralMM"); }
325 return fIntegralMM[i];
326}
327
329{
330 if(i < 1 || i > fSplineNumber) { CallError(i, "GetIntegralPlasmon"); }
331 return fIntegralPlasmon[i];
332}
333
335{
336 if(i < 1 || i > fSplineNumber) { CallError(i, "GetIntegralResonance"); }
337 return fIntegralResonance[i];
338}
339
340#endif
341
342// ----------------- end of G4PAIxSection header file -------------------
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
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()
G4double SumOverInterCerenkov(G4int intervalNumber)
G4double GetStepMMLoss(G4double step)
G4double GetIntegralPAIdEdx(G4int i) const
G4double GetDifPAIxSection(G4int i)
G4double GetMeanResonanceLoss() const
G4double SumOverBordCerenkov(G4int intervalNumber, G4double energy)
G4double GetStepPlasmonLoss(G4double step)
G4double GetRutherfordEnergyTransfer()
G4int GetNumberOfGammas() const
G4double SumOverBordResonance(G4int intervalNumber, G4double energy)
G4double GetPAIdNdxPlasmon(G4int i)
G4double SumOverInterMM(G4int intervalNumber)
G4double GetEnergyInterval(G4int i)
G4double RePartDielectricConst(G4double energy)
G4double GetElectronRange(G4double energy)
G4double GetIntegralPlasmon(G4int i) const
G4double GetPAItable(G4int i, G4int j) const
G4double SumOverInterPlasmon(G4int intervalNumber)
G4double GetMMEnergyTransfer()
G4double PAIdNdxMM(G4int intervalNumber, G4double betaGammaSq)
G4double GetPAIdNdxCerenkov(G4int i)
G4double SumOverBorderdEdx(G4int intervalNumber, G4double energy)
G4double SumOverBordMM(G4int intervalNumber, G4double energy)
G4double DifPAIxSection(G4int intervalNumber, G4double betaGammaSq)
G4double GetPAIdNdxMM(G4int i)
G4double GetMeanPlasmonLoss() const
G4double GetPAIdNdxResonance(G4int i)
G4double PAIdNdxCerenkov(G4int intervalNumber, G4double betaGammaSq)
G4double GetStepResonanceLoss(G4double step)
G4double GetResonanceEnergyTransfer()
G4double PAIdNdxResonance(G4int intervalNumber, G4double betaGammaSq)
void SplainPAI(G4double betaGammaSq)
G4double GetEnergyTransfer()
G4double GetCerenkovEnergyTransfer()
G4double GetMeanCerenkovLoss() const
G4double GetSplineEnergy(G4int i) const
G4double GetIntegralMM(G4int i) const
void NormShift(G4double betaGammaSq)
G4int GetIntervalNumber() const
G4double PAIdNdxPlasmon(G4int intervalNumber, G4double betaGammaSq)
G4double GetPhotonRange(G4double energy)
void IntegralPAIxSection()
G4double GetNormalizationCof() const
G4double SumOverInterval(G4int intervalNumber)
G4double SumOverBorder(G4int intervalNumber, G4double energy)
G4double GetStepCerenkovLoss(G4double step)
G4double GetIntegralPAIxSection(G4int i) const
G4double GetIntegralResonance(G4int i) const
G4double GetMeanEnergyLoss() const
G4double GetStepEnergyLoss(G4double step)
G4double GetMeanMMLoss() const
void IntegralCerenkov()
G4double SumOverInterResonance(G4int intervalNumber)
void IntegralResonance()
G4double SumOverBordPlasmon(G4int intervalNumber, G4double energy)
G4double GetLorentzFactor(G4int i) const