ECCE @ EIC Software
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4PAIySection.hh
Go to the documentation of this file. Or view the newest version in sPHENIX GitHub for file G4PAIySection.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 // G4PAIySection.hh -- header file
30 //
31 //
32 // Preparation of ionizing collision cross section according to Photo Absorption
33 // Ionization (PAI) model for simulation of ionization energy losses in very thin
34 // absorbers. Author: Vladimir.Grichine@cern.ch
35 //
36 // History:
37 //
38 // 01.10.07, V.Ivanchenko create using V.Grichine G4PAIxSection class
39 // 21.11.10, V.Grichine fVerbose and SetVerbose added
40 // 28.10.11, V.Ivanchenko Migration of exceptions to the new design
41 
42 #ifndef G4PAIYSECTION_HH
43 #define G4PAIYSECTION_HH
44 
45 #include "G4ios.hh"
46 #include "globals.hh"
47 #include "Randomize.hh"
48 
49 #include "G4SandiaTable.hh"
50 
52 {
53 public:
54 
55  explicit G4PAIySection();
56 
57  ~G4PAIySection() = default;
58 
59  void Initialize(const G4Material* material, G4double maxEnergyTransfer,
60  G4double betaGammaSq, G4SandiaTable*);
61 
62  void ComputeLowEnergyCof(const G4Material* material);
63 
64  void InitPAI();
65 
66  void NormShift( G4double betaGammaSq );
67 
68  void SplainPAI( G4double betaGammaSq );
69 
70  // Physical methods
71  G4double RutherfordIntegral( G4int intervalNumber,
72  G4double limitLow,
73  G4double limitHigh );
74 
75  G4double ImPartDielectricConst( G4int intervalNumber,
76  G4double energy );
77 
79 
80  G4double DifPAIySection( G4int intervalNumber,
81  G4double betaGammaSq );
82 
83  G4double PAIdNdxCerenkov( G4int intervalNumber,
84  G4double betaGammaSq );
85 
86  G4double PAIdNdxPlasmon( G4int intervalNumber,
87  G4double betaGammaSq );
88 
89  void IntegralPAIySection();
90  void IntegralCerenkov();
91  void IntegralPlasmon();
92 
93  G4double SumOverInterval(G4int intervalNumber);
94  G4double SumOverIntervaldEdx(G4int intervalNumber);
95  G4double SumOverInterCerenkov(G4int intervalNumber);
96  G4double SumOverInterPlasmon(G4int intervalNumber);
97 
98  G4double SumOverBorder( G4int intervalNumber,
99  G4double energy );
100  G4double SumOverBorderdEdx( G4int intervalNumber,
101  G4double energy );
102  G4double SumOverBordCerenkov( G4int intervalNumber,
103  G4double energy );
104  G4double SumOverBordPlasmon( G4int intervalNumber,
105  G4double energy );
106 
110 
112 
113  // Inline access functions
114 
115  inline G4int GetNumberOfGammas() const { return fNumberOfGammas; }
116 
117  inline G4int GetSplineSize() const { return fSplineNumber; }
118 
119  inline G4int GetIntervalNumber() const { return fIntervalNumber; }
120 
122 
125  inline G4double GetPAIdNdxPlasmon(G4int i){ return fdNdxPlasmon[i]; }
126 
127  inline G4double GetMeanEnergyLoss() const {return fIntegralPAIySection[0]; }
128  inline G4double GetMeanCerenkovLoss() const {return fIntegralCerenkov[0]; }
129  inline G4double GetMeanPlasmonLoss() const {return fIntegralPlasmon[0]; }
130 
131  inline G4double GetNormalizationCof() const { return fNormalizationCof; }
132 
133  inline G4double GetPAItable(G4int i,G4int j) const;
134 
135  inline G4double GetSplineEnergy(G4int i) const;
136 
137  inline G4double GetIntegralPAIySection(G4int i) const;
138  inline G4double GetIntegralPAIdEdx(G4int i) const;
139  inline G4double GetIntegralCerenkov(G4int i) const;
140  inline G4double GetIntegralPlasmon(G4int i) const;
141 
142  inline void SetVerbose(G4int v) { fVerbose = v; };
143 
144 private :
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 
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 
172 
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 -------------------