source: trunk/source/processes/electromagnetic/standard/include/G4PAIxSection.hh @ 1315

Last change on this file since 1315 was 1228, checked in by garnier, 14 years ago

update geant4.9.3 tag

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