Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-01-18 09:58:49

0001 //
0002 // ********************************************************************
0003 // * License and Disclaimer                                           *
0004 // *                                                                  *
0005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
0006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
0007 // * conditions of the Geant4 Software License,  included in the file *
0008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
0009 // * include a list of copyright holders.                             *
0010 // *                                                                  *
0011 // * Neither the authors of this software system, nor their employing *
0012 // * institutes,nor the agencies providing financial support for this *
0013 // * work  make  any representation or  warranty, express or implied, *
0014 // * regarding  this  software system or assume any liability for its *
0015 // * use.  Please see the license in the file  LICENSE  and URL above *
0016 // * for the full disclaimer and the limitation of liability.         *
0017 // *                                                                  *
0018 // * This  code  implementation is the result of  the  scientific and *
0019 // * technical work of the GEANT4 collaboration.                      *
0020 // * By using,  copying,  modifying or  distributing the software (or *
0021 // * any work based  on the software)  you  agree  to acknowledge its *
0022 // * use  in  resulting  scientific  publications,  and indicate your *
0023 // * acceptance of all terms of the Geant4 Software license.          *
0024 // ********************************************************************
0025 // 
0026 // G4PAIxSection.hh -- header file
0027 //
0028 // GEANT 4 class header file --- Copyright CERN 1995
0029 // CERB Geneva Switzerland
0030 //
0031 // for information related to this code, please, contact
0032 // CERN, CN Division, ASD Group
0033 //
0034 // Preparation of ionizing collision cross section according to Photo Absorption 
0035 // Ionization (PAI) model for simulation of ionization energy losses in very thin
0036 // absorbers. Author: Vladimir.Grichine@cern.ch
0037 //
0038 // History:
0039 //
0040 // 28.10.11, V. Ivanchenko: Migration of exceptions to the new design 
0041 // 19.10.03, V. Grichine: Integral dEdx was added for G4PAIModel class 
0042 // 13.05.03, V. Grichine: Numerical instability was fixed in SumOverInterval/Border 
0043 //                        functions
0044 // 10.02.02, V. Grichine: New functions and arrays/gets for Cerenkov and 
0045 //                        plasmon collisions dN/dx
0046 // 27.10.99, V. Grichine: Bug fixed in constructors, 3rd constructor and 
0047 //                        GetStepEnergyLoss(step) were added, fDelta = 0.005
0048 // 30.11.97, V. Grichine: 2nd version 
0049 // 11.06.97, V. Grichine: 1st version
0050 
0051 #ifndef G4PAIXSECTION_HH
0052 #define G4PAIXSECTION_HH
0053 
0054 #include "G4ios.hh"
0055 #include "globals.hh"
0056 #include "Randomize.hh"
0057 
0058 #include "G4SandiaTable.hh"
0059 
0060 class G4MaterialCutsCouple;
0061 class G4Sandiatable;
0062 
0063 
0064 class G4PAIxSection
0065 {
0066 public:
0067       // Constructors
0068   G4PAIxSection();
0069   G4PAIxSection( G4MaterialCutsCouple* matCC);
0070       
0071   G4PAIxSection( G4int materialIndex, G4double maxEnergyTransfer );
0072       
0073   G4PAIxSection( G4int materialIndex,           // for proton loss table
0074                  G4double maxEnergyTransfer,
0075                  G4double betaGammaSq ,
0076                          G4double** photoAbsCof, G4int intNumber         );
0077 
0078   G4PAIxSection( G4int materialIndex,           // test constructor
0079                  G4double maxEnergyTransfer,
0080                  G4double betaGammaSq          );
0081       
0082   ~G4PAIxSection();
0083 
0084   void Initialize(const G4Material* material, G4double maxEnergyTransfer, 
0085           G4double betaGammaSq, G4SandiaTable*);
0086       
0087   // General control functions
0088 
0089   void     ComputeLowEnergyCof(const G4Material* material);
0090   void     ComputeLowEnergyCof();
0091           
0092   void InitPAI();
0093 
0094   void NormShift( G4double betaGammaSq );
0095 
0096   void SplainPAI( G4double betaGammaSq );
0097           
0098   // Physical methods
0099           
0100   G4double RutherfordIntegral( G4int intervalNumber,
0101                    G4double limitLow,
0102                    G4double limitHigh );
0103 
0104   G4double ImPartDielectricConst( G4int intervalNumber,
0105                   G4double energy );
0106 
0107   G4double GetPhotonRange( G4double energy );
0108   G4double GetElectronRange( G4double energy );
0109 
0110   G4double RePartDielectricConst(G4double energy);
0111 
0112   G4double DifPAIxSection( G4int intervalNumber,
0113                G4double betaGammaSq );
0114 
0115   G4double PAIdNdxCerenkov( G4int intervalNumber,
0116                 G4double betaGammaSq );
0117   G4double PAIdNdxMM( G4int intervalNumber,
0118               G4double betaGammaSq );
0119 
0120   G4double PAIdNdxPlasmon( G4int intervalNumber,
0121                G4double betaGammaSq );
0122 
0123   G4double PAIdNdxResonance( G4int intervalNumber,
0124                  G4double betaGammaSq );
0125 
0126 
0127   void     IntegralPAIxSection();
0128   void     IntegralCerenkov();
0129   void     IntegralMM();
0130   void     IntegralPlasmon();
0131   void     IntegralResonance();
0132 
0133   G4double SumOverInterval(G4int intervalNumber);
0134   G4double SumOverIntervaldEdx(G4int intervalNumber);
0135   G4double SumOverInterCerenkov(G4int intervalNumber);
0136   G4double SumOverInterMM(G4int intervalNumber);
0137   G4double SumOverInterPlasmon(G4int intervalNumber);
0138   G4double SumOverInterResonance(G4int intervalNumber);
0139 
0140   G4double SumOverBorder( G4int intervalNumber,
0141               G4double energy          );
0142   G4double SumOverBorderdEdx( G4int intervalNumber,
0143                   G4double energy          );
0144   G4double SumOverBordCerenkov( G4int intervalNumber,
0145                 G4double energy          );
0146   G4double SumOverBordMM( G4int intervalNumber,
0147               G4double energy          );
0148   G4double SumOverBordPlasmon( G4int intervalNumber,
0149                    G4double energy          );
0150   G4double SumOverBordResonance( G4int intervalNumber,
0151                  G4double energy          );
0152 
0153   G4double GetStepEnergyLoss( G4double step );
0154   G4double GetStepCerenkovLoss( G4double step );
0155   G4double GetStepMMLoss( G4double step );
0156   G4double GetStepPlasmonLoss( G4double step );
0157   G4double GetStepResonanceLoss( G4double step );
0158      
0159   G4double GetEnergyTransfer();
0160   G4double GetCerenkovEnergyTransfer();
0161   G4double GetMMEnergyTransfer();
0162   G4double GetPlasmonEnergyTransfer();
0163   G4double GetResonanceEnergyTransfer();
0164   G4double GetRutherfordEnergyTransfer();
0165      
0166   // Inline access functions
0167 
0168   G4int GetNumberOfGammas() const { return fNumberOfGammas; }
0169       
0170   G4int GetSplineSize() const { return fSplineNumber; }
0171       
0172   G4int GetIntervalNumber() const { return fIntervalNumber; }
0173 
0174   G4double GetEnergyInterval(G4int i){ return fEnergyInterval[i]; } 
0175 
0176   G4double GetDifPAIxSection(G4int i){ return fDifPAIxSection[i]; } 
0177   G4double GetPAIdNdxCerenkov(G4int i){ return fdNdxCerenkov[i]; } 
0178   G4double GetPAIdNdxMM(G4int i){ return fdNdxMM[i]; } 
0179   G4double GetPAIdNdxPlasmon(G4int i){ return fdNdxPlasmon[i]; } 
0180   G4double GetPAIdNdxResonance(G4int i){ return fdNdxResonance[i]; } 
0181       
0182   G4double GetMeanEnergyLoss() const {return fIntegralPAIxSection[0]; }
0183   G4double GetMeanCerenkovLoss() const {return fIntegralCerenkov[0]; }
0184   G4double GetMeanMMLoss() const {return fIntegralMM[0]; }
0185   G4double GetMeanPlasmonLoss() const {return fIntegralPlasmon[0]; }
0186   G4double GetMeanResonanceLoss() const {return fIntegralResonance[0]; }
0187 
0188   G4double GetNormalizationCof() const { return fNormalizationCof; }
0189 
0190   G4double GetLowEnergyCof() const { return fLowEnergyCof; }
0191 
0192   G4double GetLorentzFactor(G4int i) const;
0193 
0194   inline void SetVerbose(G4int v) { fVerbose=v; };
0195   
0196         
0197   inline G4double GetPAItable(G4int i,G4int j) const;
0198           
0199   inline G4double GetSplineEnergy(G4int i) const;
0200       
0201   inline G4double GetIntegralPAIxSection(G4int i) const;
0202   inline G4double GetIntegralPAIdEdx(G4int i) const;
0203   inline G4double GetIntegralCerenkov(G4int i) const;
0204   inline G4double GetIntegralMM(G4int i) const;
0205   inline G4double GetIntegralPlasmon(G4int i) const;
0206   inline G4double GetIntegralResonance(G4int i) const;
0207 
0208   G4PAIxSection & operator=(const G4PAIxSection &right) = delete;
0209   G4PAIxSection(const G4PAIxSection&) = delete;
0210 
0211 private :
0212 
0213   void CallError(G4int i, const G4String& methodName) const;
0214 
0215   // Local class constants
0216  
0217   static const G4double fDelta; // energy shift from interval border = 0.001
0218   static const G4double fError; // error in lin-log approximation = 0.005
0219 
0220   static G4int fNumberOfGammas; // = 111;
0221   static const G4double fLorentzFactor[112]; // static gamma array
0222 
0223   static 
0224   const G4int fRefGammaNumber; // The number of gamma for creation of spline (15)
0225 
0226   G4int    fIntervalNumber ;   // The number of energy intervals
0227   G4double fNormalizationCof;  // Normalization cof for PhotoAbsorptionXsection
0228 
0229   G4int fMaterialIndex;      // current material index
0230   G4double fDensity;         // Current density
0231   G4double fElectronDensity; // Current electron (number) density
0232   G4double fLowEnergyCof;    // Correction cof for low energy region
0233   G4int    fSplineNumber;    // Current size of spline
0234   G4int    fVerbose;         // verbose flag
0235 
0236   // Arrays of Sandia coefficients
0237 
0238   G4OrderedTable* fMatSandiaMatrix;
0239 
0240   G4SandiaTable*  fSandia;
0241 
0242   G4DataVector fEnergyInterval;
0243   G4DataVector fA1; 
0244   G4DataVector fA2;
0245   G4DataVector fA3; 
0246   G4DataVector fA4;
0247 
0248   static
0249   const G4int  fMaxSplineSize ; // Max size of output splain arrays = 500
0250 
0251   G4DataVector          fSplineEnergy;     // energy points of splain
0252   G4DataVector   fRePartDielectricConst;   // Real part of dielectric const
0253   G4DataVector   fImPartDielectricConst;   // Imaginary part of dielectric const
0254   G4DataVector            fIntegralTerm;   // Integral term in PAI cross section
0255   G4DataVector          fDifPAIxSection;   // Differential PAI cross section
0256   G4DataVector            fdNdxCerenkov;   // dNdx of Cerenkov collisions
0257   G4DataVector            fdNdxPlasmon;    // dNdx of Plasmon collisions
0258   G4DataVector            fdNdxMM;         // dNdx of MM-Cerenkov collisions
0259   G4DataVector            fdNdxResonance;  // dNdx of Resonance collisions
0260 
0261   G4DataVector     fIntegralPAIxSection;   // Integral PAI cross section  ?
0262   G4DataVector     fIntegralPAIdEdx;       // Integral PAI dEdx  ?
0263   G4DataVector     fIntegralCerenkov;      // Integral Cerenkov N>omega  ?
0264   G4DataVector     fIntegralPlasmon;       // Integral Plasmon N>omega  ?
0265   G4DataVector     fIntegralMM;            // Integral MM N>omega  ?
0266   G4DataVector     fIntegralResonance;     // Integral resonance N>omega  ?
0267 
0268   G4double fPAItable[500][112]; // Output array
0269 
0270 };    
0271 
0272 ////////////////  Inline methods //////////////////////////////////
0273 //
0274 
0275 inline G4double G4PAIxSection::GetPAItable(G4int i, G4int j) const
0276 {
0277    return fPAItable[i][j];
0278 }
0279 
0280 inline G4double G4PAIxSection::GetSplineEnergy(G4int i) const 
0281 {
0282   if(i < 1 || i > fSplineNumber) { CallError(i, "GetSplineEnergy"); }
0283   return fSplineEnergy[i];
0284 }
0285       
0286 inline G4double G4PAIxSection::GetIntegralPAIxSection(G4int i) const 
0287 {
0288   if(i < 1 || i > fSplineNumber) { CallError(i, "GetIntegralPAIxSection"); }
0289   return fIntegralPAIxSection[i];
0290 }
0291 
0292 inline G4double G4PAIxSection::GetIntegralPAIdEdx(G4int i) const 
0293 {
0294   if(i < 1 || i > fSplineNumber) { CallError(i, "GetIntegralPAIdEdx"); }
0295   return fIntegralPAIdEdx[i];
0296 }
0297 
0298 inline G4double G4PAIxSection::GetIntegralCerenkov(G4int i) const 
0299 {
0300   if(i < 1 || i > fSplineNumber) { CallError(i, "GetIntegralCerenkov"); }
0301   return fIntegralCerenkov[i];
0302 }
0303 
0304 inline G4double G4PAIxSection::GetIntegralMM(G4int i) const 
0305 {
0306   if(i < 1 || i > fSplineNumber) { CallError(i, "GetIntegralMM"); }
0307   return fIntegralMM[i];
0308 }
0309 
0310 inline G4double G4PAIxSection::GetIntegralPlasmon(G4int i) const 
0311 {
0312   if(i < 1 || i > fSplineNumber) { CallError(i, "GetIntegralPlasmon"); }
0313   return fIntegralPlasmon[i];
0314 }
0315 
0316 inline G4double G4PAIxSection::GetIntegralResonance(G4int i) const 
0317 {
0318   if(i < 1 || i > fSplineNumber) { CallError(i, "GetIntegralResonance"); }
0319   return fIntegralResonance[i];
0320 }
0321 
0322 #endif   
0323 
0324 // -----------------   end of G4PAIxSection header file    -------------------