Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/electromagnetic/standard/include/G4PAIxSection.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 // 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 
 60 class G4MaterialCutsCouple;
 61 class G4Sandiatable;
 62 
 63 
 64 class G4PAIxSection
 65 {
 66 public:
 67     // Constructors
 68   G4PAIxSection();
 69   G4PAIxSection( G4MaterialCutsCouple* matCC);
 70     
 71   G4PAIxSection( G4int materialIndex, G4double maxEnergyTransfer );
 72     
 73   G4PAIxSection( G4int materialIndex,           // for proton loss table
 74              G4double maxEnergyTransfer,
 75              G4double betaGammaSq ,
 76                          G4double** photoAbsCof, G4int intNumber         );
 77 
 78   G4PAIxSection( G4int materialIndex,           // test constructor
 79              G4double maxEnergyTransfer,
 80              G4double betaGammaSq          );
 81     
 82   ~G4PAIxSection();
 83 
 84   void Initialize(const G4Material* material, G4double maxEnergyTransfer, 
 85       G4double betaGammaSq, G4SandiaTable*);
 86     
 87   // General control functions
 88 
 89   void     ComputeLowEnergyCof(const G4Material* material);
 90   void     ComputeLowEnergyCof();
 91           
 92   void InitPAI();
 93 
 94   void NormShift( G4double betaGammaSq );
 95 
 96   void SplainPAI( G4double betaGammaSq );
 97         
 98   // Physical methods
 99           
100   G4double RutherfordIntegral( G4int intervalNumber,
101              G4double limitLow,
102              G4double limitHigh );
103 
104   G4double ImPartDielectricConst( G4int intervalNumber,
105           G4double energy );
106 
107   G4double GetPhotonRange( G4double energy );
108   G4double GetElectronRange( G4double energy );
109 
110   G4double RePartDielectricConst(G4double energy);
111 
112   G4double DifPAIxSection( G4int intervalNumber,
113          G4double betaGammaSq );
114 
115   G4double PAIdNdxCerenkov( G4int intervalNumber,
116           G4double betaGammaSq );
117   G4double PAIdNdxMM( G4int intervalNumber,
118           G4double betaGammaSq );
119 
120   G4double PAIdNdxPlasmon( G4int intervalNumber,
121          G4double betaGammaSq );
122 
123   G4double PAIdNdxResonance( G4int intervalNumber,
124            G4double betaGammaSq );
125 
126 
127   void     IntegralPAIxSection();
128   void     IntegralCerenkov();
129   void     IntegralMM();
130   void     IntegralPlasmon();
131   void     IntegralResonance();
132 
133   G4double SumOverInterval(G4int intervalNumber);
134   G4double SumOverIntervaldEdx(G4int intervalNumber);
135   G4double SumOverInterCerenkov(G4int intervalNumber);
136   G4double SumOverInterMM(G4int intervalNumber);
137   G4double SumOverInterPlasmon(G4int intervalNumber);
138   G4double SumOverInterResonance(G4int intervalNumber);
139 
140   G4double SumOverBorder( G4int intervalNumber,
141         G4double energy          );
142   G4double SumOverBorderdEdx( G4int intervalNumber,
143             G4double energy          );
144   G4double SumOverBordCerenkov( G4int intervalNumber,
145         G4double energy          );
146   G4double SumOverBordMM( G4int intervalNumber,
147         G4double energy          );
148   G4double SumOverBordPlasmon( G4int intervalNumber,
149              G4double energy          );
150   G4double SumOverBordResonance( G4int intervalNumber,
151          G4double energy          );
152 
153   G4double GetStepEnergyLoss( G4double step );
154   G4double GetStepCerenkovLoss( G4double step );
155   G4double GetStepMMLoss( G4double step );
156   G4double GetStepPlasmonLoss( G4double step );
157   G4double GetStepResonanceLoss( G4double step );
158    
159   G4double GetEnergyTransfer();
160   G4double GetCerenkovEnergyTransfer();
161   G4double GetMMEnergyTransfer();
162   G4double GetPlasmonEnergyTransfer();
163   G4double GetResonanceEnergyTransfer();
164   G4double GetRutherfordEnergyTransfer();
165    
166   // Inline access functions
167 
168   G4int GetNumberOfGammas() const { return fNumberOfGammas; }
169     
170   G4int GetSplineSize() const { return fSplineNumber; }
171     
172   G4int GetIntervalNumber() const { return fIntervalNumber; }
173 
174   G4double GetEnergyInterval(G4int i){ return fEnergyInterval[i]; } 
175 
176   G4double GetDifPAIxSection(G4int i){ return fDifPAIxSection[i]; } 
177   G4double GetPAIdNdxCerenkov(G4int i){ return fdNdxCerenkov[i]; } 
178   G4double GetPAIdNdxMM(G4int i){ return fdNdxMM[i]; } 
179   G4double GetPAIdNdxPlasmon(G4int i){ return fdNdxPlasmon[i]; } 
180   G4double GetPAIdNdxResonance(G4int i){ return fdNdxResonance[i]; } 
181     
182   G4double GetMeanEnergyLoss() const {return fIntegralPAIxSection[0]; }
183   G4double GetMeanCerenkovLoss() const {return fIntegralCerenkov[0]; }
184   G4double GetMeanMMLoss() const {return fIntegralMM[0]; }
185   G4double GetMeanPlasmonLoss() const {return fIntegralPlasmon[0]; }
186   G4double GetMeanResonanceLoss() const {return fIntegralResonance[0]; }
187 
188   G4double GetNormalizationCof() const { return fNormalizationCof; }
189 
190   G4double GetLowEnergyCof() const { return fLowEnergyCof; }
191 
192   G4double GetLorentzFactor(G4int i) const;
193 
194   inline void SetVerbose(G4int v) { fVerbose=v; };
195   
196         
197   inline G4double GetPAItable(G4int i,G4int j) const;
198           
199   inline G4double GetSplineEnergy(G4int i) const;
200     
201   inline G4double GetIntegralPAIxSection(G4int i) const;
202   inline G4double GetIntegralPAIdEdx(G4int i) const;
203   inline G4double GetIntegralCerenkov(G4int i) const;
204   inline G4double GetIntegralMM(G4int i) const;
205   inline G4double GetIntegralPlasmon(G4int i) const;
206   inline G4double GetIntegralResonance(G4int i) const;
207 
208   G4PAIxSection & operator=(const G4PAIxSection &right) = delete;
209   G4PAIxSection(const G4PAIxSection&) = delete;
210 
211 private :
212 
213   void CallError(G4int i, const G4String& methodName) const;
214 
215   // Local class constants
216  
217   static const G4double fDelta; // energy shift from interval border = 0.001
218   static const G4double fError; // error in lin-log approximation = 0.005
219 
220   static G4int fNumberOfGammas; // = 111;
221   static const G4double fLorentzFactor[112]; // static gamma array
222 
223   static 
224   const G4int fRefGammaNumber; // The number of gamma for creation of spline (15)
225 
226   G4int    fIntervalNumber ;   // The number of energy intervals
227   G4double fNormalizationCof;  // Normalization cof for PhotoAbsorptionXsection
228 
229   G4int fMaterialIndex;      // current material index
230   G4double fDensity;         // Current density
231   G4double fElectronDensity; // Current electron (number) density
232   G4double fLowEnergyCof;    // Correction cof for low energy region
233   G4int    fSplineNumber;    // Current size of spline
234   G4int    fVerbose;         // verbose flag
235 
236   // Arrays of Sandia coefficients
237 
238   G4OrderedTable* fMatSandiaMatrix;
239 
240   G4SandiaTable*  fSandia;
241 
242   G4DataVector fEnergyInterval;
243   G4DataVector fA1; 
244   G4DataVector fA2;
245   G4DataVector fA3; 
246   G4DataVector fA4;
247 
248   static
249   const G4int  fMaxSplineSize ; // Max size of output splain arrays = 500
250 
251   G4DataVector          fSplineEnergy;     // energy points of splain
252   G4DataVector   fRePartDielectricConst;   // Real part of dielectric const
253   G4DataVector   fImPartDielectricConst;   // Imaginary part of dielectric const
254   G4DataVector            fIntegralTerm;   // Integral term in PAI cross section
255   G4DataVector          fDifPAIxSection;   // Differential PAI cross section
256   G4DataVector            fdNdxCerenkov;   // dNdx of Cerenkov collisions
257   G4DataVector            fdNdxPlasmon;    // dNdx of Plasmon collisions
258   G4DataVector            fdNdxMM;         // dNdx of MM-Cerenkov collisions
259   G4DataVector            fdNdxResonance;  // dNdx of Resonance collisions
260 
261   G4DataVector     fIntegralPAIxSection;   // Integral PAI cross section  ?
262   G4DataVector     fIntegralPAIdEdx;       // Integral PAI dEdx  ?
263   G4DataVector     fIntegralCerenkov;      // Integral Cerenkov N>omega  ?
264   G4DataVector     fIntegralPlasmon;       // Integral Plasmon N>omega  ?
265   G4DataVector     fIntegralMM;            // Integral MM N>omega  ?
266   G4DataVector     fIntegralResonance;     // Integral resonance N>omega  ?
267 
268   G4double fPAItable[500][112]; // Output array
269 
270 };    
271 
272 ////////////////  Inline methods //////////////////////////////////
273 //
274 
275 inline G4double G4PAIxSection::GetPAItable(G4int i, G4int j) const
276 {
277    return fPAItable[i][j];
278 }
279 
280 inline G4double G4PAIxSection::GetSplineEnergy(G4int i) const 
281 {
282   if(i < 1 || i > fSplineNumber) { CallError(i, "GetSplineEnergy"); }
283   return fSplineEnergy[i];
284 }
285     
286 inline G4double G4PAIxSection::GetIntegralPAIxSection(G4int i) const 
287 {
288   if(i < 1 || i > fSplineNumber) { CallError(i, "GetIntegralPAIxSection"); }
289   return fIntegralPAIxSection[i];
290 }
291 
292 inline G4double G4PAIxSection::GetIntegralPAIdEdx(G4int i) const 
293 {
294   if(i < 1 || i > fSplineNumber) { CallError(i, "GetIntegralPAIdEdx"); }
295   return fIntegralPAIdEdx[i];
296 }
297 
298 inline G4double G4PAIxSection::GetIntegralCerenkov(G4int i) const 
299 {
300   if(i < 1 || i > fSplineNumber) { CallError(i, "GetIntegralCerenkov"); }
301   return fIntegralCerenkov[i];
302 }
303 
304 inline G4double G4PAIxSection::GetIntegralMM(G4int i) const 
305 {
306   if(i < 1 || i > fSplineNumber) { CallError(i, "GetIntegralMM"); }
307   return fIntegralMM[i];
308 }
309 
310 inline G4double G4PAIxSection::GetIntegralPlasmon(G4int i) const 
311 {
312   if(i < 1 || i > fSplineNumber) { CallError(i, "GetIntegralPlasmon"); }
313   return fIntegralPlasmon[i];
314 }
315 
316 inline G4double G4PAIxSection::GetIntegralResonance(G4int i) const 
317 {
318   if(i < 1 || i > fSplineNumber) { CallError(i, "GetIntegralResonance"); }
319   return fIntegralResonance[i];
320 }
321 
322 #endif   
323 
324 // -----------------   end of G4PAIxSection header file    -------------------
325