ECCE @ EIC Software
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4PAIxSection.hh
Go to the documentation of this file. Or view the newest version in sPHENIX GitHub for file G4PAIxSection.hh
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 //
28 //
29 // G4PAIxSection.hh -- header file
30 //
31 // GEANT 4 class header file --- Copyright CERN 1995
32 // CERB Geneva Switzerland
33 //
34 // for information related to this code, please, contact
35 // CERN, CN Division, ASD Group
36 //
37 // Preparation of ionizing collision cross section according to Photo Absorption
38 // Ionization (PAI) model for simulation of ionization energy losses in very thin
39 // absorbers. Author: Vladimir.Grichine@cern.ch
40 //
41 // History:
42 //
43 // 28.10.11, V. Ivanchenko: Migration of exceptions to the new design
44 // 19.10.03, V. Grichine: Integral dEdx was added for G4PAIModel class
45 // 13.05.03, V. Grichine: Numerical instability was fixed in SumOverInterval/Border
46 // functions
47 // 10.02.02, V. Grichine: New functions and arrays/gets for Cerenkov and
48 // plasmon collisions dN/dx
49 // 27.10.99, V. Grichine: Bug fixed in constructors, 3rd constructor and
50 // GetStepEnergyLoss(step) were added, fDelta = 0.005
51 // 30.11.97, V. Grichine: 2nd version
52 // 11.06.97, V. Grichine: 1st version
53 
54 #ifndef G4PAIXSECTION_HH
55 #define G4PAIXSECTION_HH
56 
57 #include "G4ios.hh"
58 #include "globals.hh"
59 #include "Randomize.hh"
60 
61 #include"G4SandiaTable.hh"
62 
64 class G4Sandiatable;
65 
66 
68 {
69 public:
70  // Constructors
71  G4PAIxSection();
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 
87 
88  void Initialize(const G4Material* material, G4double maxEnergyTransfer,
89  G4double betaGammaSq, G4SandiaTable*);
90 
91  // General control functions
92 
93  void ComputeLowEnergyCof(const G4Material* material);
94  void ComputeLowEnergyCof();
95 
96  void InitPAI();
97 
98  void NormShift( G4double betaGammaSq );
99 
100  void SplainPAI( G4double betaGammaSq );
101 
102  // Physical methods
103 
104  G4double RutherfordIntegral( G4int intervalNumber,
105  G4double limitLow,
106  G4double limitHigh );
107 
108  G4double ImPartDielectricConst( G4int intervalNumber,
109  G4double energy );
110 
111  G4double GetPhotonRange( G4double energy );
113 
115 
116  G4double DifPAIxSection( G4int intervalNumber,
117  G4double betaGammaSq );
118 
119  G4double PAIdNdxCerenkov( G4int intervalNumber,
120  G4double betaGammaSq );
121  G4double PAIdNdxMM( G4int intervalNumber,
122  G4double betaGammaSq );
123 
124  G4double PAIdNdxPlasmon( G4int intervalNumber,
125  G4double betaGammaSq );
126 
127  G4double PAIdNdxResonance( G4int intervalNumber,
128  G4double betaGammaSq );
129 
130 
131  void IntegralPAIxSection();
132  void IntegralCerenkov();
133  void IntegralMM();
134  void IntegralPlasmon();
135  void IntegralResonance();
136 
137  G4double SumOverInterval(G4int intervalNumber);
138  G4double SumOverIntervaldEdx(G4int intervalNumber);
139  G4double SumOverInterCerenkov(G4int intervalNumber);
140  G4double SumOverInterMM(G4int intervalNumber);
141  G4double SumOverInterPlasmon(G4int intervalNumber);
142  G4double SumOverInterResonance(G4int intervalNumber);
143 
144  G4double SumOverBorder( G4int intervalNumber,
145  G4double energy );
146  G4double SumOverBorderdEdx( G4int intervalNumber,
147  G4double energy );
148  G4double SumOverBordCerenkov( G4int intervalNumber,
149  G4double energy );
150  G4double SumOverBordMM( G4int intervalNumber,
151  G4double energy );
152  G4double SumOverBordPlasmon( G4int intervalNumber,
153  G4double energy );
154  G4double SumOverBordResonance( G4int intervalNumber,
155  G4double energy );
156 
162 
169 
170  // Inline access functions
171 
173 
174  G4int GetSplineSize() const { return fSplineNumber; }
175 
177 
179 
182  G4double GetPAIdNdxMM(G4int i){ return fdNdxMM[i]; }
185 
188  G4double GetMeanMMLoss() const {return fIntegralMM[0]; }
191 
193 
195 
197 
198  inline void SetVerbose(G4int v) { fVerbose=v; };
199 
200 
201  inline G4double GetPAItable(G4int i,G4int j) 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 
212 private :
213 
214  void CallError(G4int i, const G4String& methodName) const;
215 
216  G4PAIxSection & operator=(const G4PAIxSection &right) = delete;
217  G4PAIxSection(const G4PAIxSection&) = delete;
218 
219  // Local class constants
220 
221  static const G4double fDelta; // energy shift from interval border = 0.001
222  static const G4double fError; // error in lin-log approximation = 0.005
223 
224  static G4int fNumberOfGammas; // = 111;
225  static const G4double fLorentzFactor[112]; // static gamma array
226 
227  static
228  const G4int fRefGammaNumber; //The number of gamma for creation of spline (15)
229 
230  G4int fIntervalNumber ; // The number of energy intervals
231  G4double 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  G4double fLowEnergyCof; // Correction cof for low energy region
239  G4int fSplineNumber; // Current size of spline
240  G4int fVerbose; // verbose flag
241 
242  // Arrays of Sandia coefficients
243 
245 
247 
253 
254  static
255  const G4int fMaxSplineSize ; // Max size of output splain arrays = 500
256 
257  G4DataVector fSplineEnergy; // energy points of splain
258  G4DataVector fRePartDielectricConst; // Real part of dielectric const
259  G4DataVector fImPartDielectricConst; // Imaginary part of dielectric const
260  G4DataVector fIntegralTerm; // Integral term in PAI cross section
261  G4DataVector fDifPAIxSection; // Differential PAI cross section
262  G4DataVector fdNdxCerenkov; // dNdx of Cerenkov collisions
263  G4DataVector fdNdxPlasmon; // dNdx of Plasmon collisions
264  G4DataVector fdNdxMM; // dNdx of MM-Cerenkov collisions
265  G4DataVector fdNdxResonance; // dNdx of Resonance collisions
266 
267  G4DataVector fIntegralPAIxSection; // Integral PAI cross section ?
268  G4DataVector fIntegralPAIdEdx; // Integral PAI dEdx ?
269  G4DataVector fIntegralCerenkov; // Integral Cerenkov N>omega ?
270  G4DataVector fIntegralPlasmon; // Integral Plasmon N>omega ?
271  G4DataVector fIntegralMM; // Integral MM N>omega ?
272  G4DataVector fIntegralResonance; // Integral resonance N>omega ?
273 
274  G4double fPAItable[500][112]; // Output array
275 
276 };
277 
279 //
280 
282 {
283  return fPAItable[i][j];
284 }
285 
287 {
288  if(i < 1 || i > fSplineNumber) { CallError(i, "GetSplineEnergy"); }
289  return fSplineEnergy[i];
290 }
291 
293 {
294  if(i < 1 || i > fSplineNumber) { CallError(i, "GetIntegralPAIxSection"); }
295  return fIntegralPAIxSection[i];
296 }
297 
299 {
300  if(i < 1 || i > fSplineNumber) { CallError(i, "GetIntegralPAIdEdx"); }
301  return fIntegralPAIdEdx[i];
302 }
303 
305 {
306  if(i < 1 || i > fSplineNumber) { CallError(i, "GetIntegralCerenkov"); }
307  return fIntegralCerenkov[i];
308 }
309 
311 {
312  if(i < 1 || i > fSplineNumber) { CallError(i, "GetIntegralMM"); }
313  return fIntegralMM[i];
314 }
315 
317 {
318  if(i < 1 || i > fSplineNumber) { CallError(i, "GetIntegralPlasmon"); }
319  return fIntegralPlasmon[i];
320 }
321 
323 {
324  if(i < 1 || i > fSplineNumber) { CallError(i, "GetIntegralResonance"); }
325  return fIntegralResonance[i];
326 }
327 
328 #endif
329 
330 // ----------------- end of G4PAIxSection header file -------------------