Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/electromagnetic/standard/include/G4PAIySection.hh

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

  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: Vladimir.Grichine@cern.ch
 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 
 48 class G4PAIySection
 49 {
 50 public:
 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 
 75   G4double RePartDielectricConst(G4double energy);
 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 
 86   void     IntegralPAIySection();
 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 
104   G4double GetStepEnergyLoss( G4double step );
105   G4double GetStepCerenkovLoss( G4double step );
106   G4double GetStepPlasmonLoss( G4double step );
107 
108   G4double GetLorentzFactor(G4int j) const;
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 
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 
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 
198 inline G4double G4PAIySection::GetPAItable(G4int i, G4int j) const
199 {
200    return fPAItable[i][j];
201 }
202 
203 inline G4double G4PAIySection::GetSplineEnergy(G4int i) const 
204 {
205   if(i < 1 || i > fSplineNumber) { CallError(i, "GetSplineEnergy"); }
206   return fSplineEnergy[i];
207 }
208   
209 inline G4double G4PAIySection::GetIntegralPAIySection(G4int i) const 
210 {
211   if(i < 1 || i > fSplineNumber) { CallError(i, "GetIntegralPAIySection"); }
212   return fIntegralPAIySection[i];
213 }
214 
215 inline G4double G4PAIySection::GetIntegralPAIdEdx(G4int i) const 
216 {
217   if(i < 1 || i > fSplineNumber) { CallError(i, "GetIntegralPAIdEdx"); }
218   return fIntegralPAIdEdx[i];
219 }
220 
221 inline G4double G4PAIySection::GetIntegralCerenkov(G4int i) const 
222 {
223   if(i < 1 || i > fSplineNumber) { CallError(i, "GetIntegralCerenkov"); }
224   return fIntegralCerenkov[i];
225 }
226 
227 inline G4double G4PAIySection::GetIntegralPlasmon(G4int i) const 
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    -------------------
236