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 // G4PAIySection.hh -- header file
0027 //
0028 //
0029 // Preparation of ionizing collision cross section according to Photo Absorption 
0030 // Ionization (PAI) model for simulation of ionization energy losses in very thin
0031 // absorbers. Author: Vladimir.Grichine@cern.ch
0032 //
0033 // History:
0034 //
0035 // 01.10.07, V.Ivanchenko create using V.Grichine G4PAIxSection class
0036 // 21.11.10, V.Grichine   fVerbose and SetVerbose added
0037 // 28.10.11, V.Ivanchenko Migration of exceptions to the new design 
0038 
0039 #ifndef G4PAIYSECTION_HH
0040 #define G4PAIYSECTION_HH
0041 
0042 #include "G4ios.hh"
0043 #include "globals.hh"
0044 #include "Randomize.hh"
0045 
0046 #include "G4SandiaTable.hh"
0047 
0048 class G4PAIySection
0049 {
0050 public:
0051 
0052   explicit G4PAIySection();
0053           
0054   ~G4PAIySection() = default;
0055 
0056   void Initialize(const G4Material* material, G4double maxEnergyTransfer, 
0057                   G4double betaGammaSq, G4SandiaTable*);
0058 
0059   void ComputeLowEnergyCof(const G4Material* material);
0060 
0061   void InitPAI();
0062 
0063   void NormShift( G4double betaGammaSq );
0064   
0065   void SplainPAI( G4double betaGammaSq );
0066                     
0067   // Physical methods
0068   G4double RutherfordIntegral( G4int intervalNumber,
0069                                G4double limitLow,
0070                                G4double limitHigh     );
0071 
0072   G4double ImPartDielectricConst( G4int intervalNumber,
0073                                   G4double energy        );
0074 
0075   G4double RePartDielectricConst(G4double energy);
0076 
0077   G4double DifPAIySection( G4int intervalNumber,
0078                            G4double betaGammaSq    );
0079 
0080   G4double PAIdNdxCerenkov( G4int intervalNumber,
0081                             G4double betaGammaSq    );
0082 
0083   G4double PAIdNdxPlasmon( G4int intervalNumber,
0084                            G4double betaGammaSq    );
0085 
0086   void     IntegralPAIySection();
0087   void     IntegralCerenkov();
0088   void     IntegralPlasmon();
0089 
0090   G4double SumOverInterval(G4int intervalNumber);
0091   G4double SumOverIntervaldEdx(G4int intervalNumber);
0092   G4double SumOverInterCerenkov(G4int intervalNumber);
0093   G4double SumOverInterPlasmon(G4int intervalNumber);
0094 
0095   G4double SumOverBorder( G4int intervalNumber,
0096                           G4double energy          );
0097   G4double SumOverBorderdEdx( G4int intervalNumber,
0098                               G4double energy          );
0099   G4double SumOverBordCerenkov( G4int intervalNumber,
0100                                 G4double energy          );
0101   G4double SumOverBordPlasmon( G4int intervalNumber,
0102                                G4double energy          );
0103 
0104   G4double GetStepEnergyLoss( G4double step );
0105   G4double GetStepCerenkovLoss( G4double step );
0106   G4double GetStepPlasmonLoss( G4double step );
0107 
0108   G4double GetLorentzFactor(G4int j) const;
0109          
0110   // Inline access functions
0111 
0112   inline G4int GetNumberOfGammas() const { return fNumberOfGammas; }
0113           
0114   inline G4int GetSplineSize() const { return fSplineNumber; }
0115           
0116   inline G4int GetIntervalNumber() const { return fIntervalNumber; }
0117 
0118   inline G4double GetEnergyInterval(G4int i){ return fEnergyInterval[i]; } 
0119 
0120   inline G4double GetDifPAIySection(G4int i){ return fDifPAIySection[i]; } 
0121   inline G4double GetPAIdNdxCrenkov(G4int i){ return fdNdxCerenkov[i]; } 
0122   inline G4double GetPAIdNdxPlasmon(G4int i){ return fdNdxPlasmon[i]; } 
0123           
0124   inline G4double GetMeanEnergyLoss() const {return fIntegralPAIySection[0]; }
0125   inline G4double GetMeanCerenkovLoss() const {return fIntegralCerenkov[0]; }
0126   inline G4double GetMeanPlasmonLoss() const {return fIntegralPlasmon[0]; }
0127 
0128   inline G4double GetNormalizationCof() const { return fNormalizationCof; }
0129           
0130   inline G4double GetPAItable(G4int i,G4int j) const;
0131                     
0132   inline G4double GetSplineEnergy(G4int i) const;
0133           
0134   inline G4double GetIntegralPAIySection(G4int i) const;
0135   inline G4double GetIntegralPAIdEdx(G4int i) const;
0136   inline G4double GetIntegralCerenkov(G4int i) const;
0137   inline G4double GetIntegralPlasmon(G4int i) const;
0138 
0139   inline void SetVerbose(G4int v) { fVerbose = v; };
0140 
0141   G4PAIySection & operator=(const G4PAIySection &right) = delete;
0142   G4PAIySection(const G4PAIySection&) = delete;
0143 
0144 private :
0145 
0146   void CallError(G4int i, const G4String& methodName) const;
0147 
0148   // Local class constants
0149  
0150   static const G4double fDelta; // energy shift from interval border = 0.001
0151   static const G4double fError; // error in lin-log approximation = 0.005
0152 
0153   static G4int fNumberOfGammas; // = 111;
0154   static const G4double fLorentzFactor[112];  //  static gamma array
0155 
0156   static
0157   const G4int fRefGammaNumber; // The number of gamma for creation of spline (15)
0158 
0159   G4int    fIntervalNumber ;   //  The number of energy intervals
0160   G4double fNormalizationCof;  // Normalization cof for PhotoAbsorptionXsection
0161 
0162   G4double betaBohr;
0163   G4double betaBohr4;
0164 
0165   G4double fDensity;            // Current density
0166   G4double fElectronDensity;    // Current electron (number) density
0167   G4double fLowEnergyCof;       // Correction cof for low energy region
0168   G4int    fSplineNumber;       // Current size of spline
0169   G4int    fVerbose;            // verbose flag
0170 
0171   G4SandiaTable*  fSandia;
0172 
0173   G4DataVector fEnergyInterval;
0174   G4DataVector fA1; 
0175   G4DataVector fA2;
0176   G4DataVector fA3; 
0177   G4DataVector fA4;
0178 
0179   static
0180   const G4int  fMaxSplineSize; // Max size of output splain arrays = 500
0181 
0182   G4DataVector fSplineEnergy;          // energy points of splain
0183   G4DataVector fRePartDielectricConst; // Real part of dielectric const
0184   G4DataVector fImPartDielectricConst; // Imaginary part of dielectric const
0185   G4DataVector fIntegralTerm;          // Integral term in PAI cross section
0186   G4DataVector fDifPAIySection;        // Differential PAI cross section
0187   G4DataVector fdNdxCerenkov;          // dNdx of Cerenkov collisions
0188   G4DataVector fdNdxPlasmon;           // dNdx of Plasmon collisions
0189 
0190   G4DataVector fIntegralPAIySection;   // Integral PAI cross section  ?
0191   G4DataVector fIntegralPAIdEdx;       // Integral PAI dEdx  ?
0192   G4DataVector fIntegralCerenkov;      // Integral Cerenkov N>omega  ?
0193   G4DataVector fIntegralPlasmon;       // Integral Plasmon N>omega  ?
0194 
0195   G4double     fPAItable[500][112];    // Output array
0196 };
0197 
0198 inline G4double G4PAIySection::GetPAItable(G4int i, G4int j) const
0199 {
0200    return fPAItable[i][j];
0201 }
0202 
0203 inline G4double G4PAIySection::GetSplineEnergy(G4int i) const 
0204 {
0205   if(i < 1 || i > fSplineNumber) { CallError(i, "GetSplineEnergy"); }
0206   return fSplineEnergy[i];
0207 }
0208   
0209 inline G4double G4PAIySection::GetIntegralPAIySection(G4int i) const 
0210 {
0211   if(i < 1 || i > fSplineNumber) { CallError(i, "GetIntegralPAIySection"); }
0212   return fIntegralPAIySection[i];
0213 }
0214 
0215 inline G4double G4PAIySection::GetIntegralPAIdEdx(G4int i) const 
0216 {
0217   if(i < 1 || i > fSplineNumber) { CallError(i, "GetIntegralPAIdEdx"); }
0218   return fIntegralPAIdEdx[i];
0219 }
0220 
0221 inline G4double G4PAIySection::GetIntegralCerenkov(G4int i) const 
0222 {
0223   if(i < 1 || i > fSplineNumber) { CallError(i, "GetIntegralCerenkov"); }
0224   return fIntegralCerenkov[i];
0225 }
0226 
0227 inline G4double G4PAIySection::GetIntegralPlasmon(G4int i) const 
0228 {
0229   if(i < 1 || i > fSplineNumber) { CallError(i, "GetIntegralPlasmon"); }
0230   return fIntegralPlasmon[i];
0231 }
0232 
0233 #endif   
0234 
0235 // -----------------   end of G4PAIySection header file    -------------------