Back to home page

EIC code displayed by LXR

 
 

    


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

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 //
0027 // -------------------------------------------------------------------
0028 //
0029 // GEANT4 Class header file
0030 //
0031 //
0032 // File name:     G4eDPWAElasticDCS
0033 //
0034 // Author:        Mihaly Novak
0035 //
0036 // Creation date: 02.07.2020
0037 //
0038 // Modifications:
0039 //
0040 // Class Description:
0041 //
0042 // Contains numerical Differential Cross Sections (DCS) for e-/e+ Coulomb
0043 // scattering computed by Dirac Partial Wave Analysis (DPWA) [1]:
0044 // - electrostatic interaction, with a local exchange correction in the case of
0045 //   electrons (using Dirac-Fock e- densities; finite nuclear size with Fermi
0046 //   charge distribution; exchange potential with Furness and McCarthy for e-)[2]
0047 // - correlation-polarization (projectiles cause the polarization of the charge
0048 //   cloud of the target atom and the induced dipole moment acts back on the
0049 //   projectile) was accounted by using Local-Density Approximation (LDA) [2]
0050 // - absorption: not included since it's an inelastic channel [2] (the cor-
0051 //   responding excitations needs to be modelled by a separate, independent,
0052 //   inelastic model).
0053 // Using the above mentioned DPWA computation with a free atom approximation
0054 // might lead to questionable results below few hundred [eV] where possible
0055 // solid state or bounding effects might start to affect the potential.
0056 // Nevertheless, the lower energy was set to 10 eV in order to provide(at least)
0057 // some model even at low energies (with this caution). The highest projectile
0058 // kinetic energy is 100 [MeV].
0059 //
0060 // The class provides interface methods for elastic, first-, second-transport
0061 // cross section computations as well as for sampling cosine of polar angular
0062 // deflections. These interface methods are also available for resricted cross
0063 // section computations and angular deflection sampling.
0064 //
0065 // References:
0066 //
0067 // [1] Salvat, F., Jablonski, A. and Powell, C.J., 2005. ELSEPA—Dirac partial-
0068 //     wave calculation of elastic scattering of electrons and positrons by
0069 //     atoms, positive ions and molecules. Computer physics communications,
0070 //     165(2), pp.157-190.
0071 // [2] Salvat, F., 2003. Optical-model potential for electron and positron
0072 //     elastic scattering by atoms. Physical Review A, 68(1), p.012708.
0073 // [3] Benedito, E., Fernández-Varea, J.M. and Salvat, F.,2001. Mixed simulation
0074 //     of the multiple elastic scattering of electrons and positrons using
0075 //     partial-wave differential cross-sections. Nuclear Instruments and Methods
0076 //     in Physics Research Section B: Beam Interactions with Materials and Atoms,
0077 //     174(1-2), pp.91-110.
0078 //
0079 // -------------------------------------------------------------------
0080 
0081 #ifndef G4eDPWAElasticDCS_h
0082 #define G4eDPWAElasticDCS_h 1
0083 
0084 
0085 #include <vector>
0086 #include <fstream>
0087 #include <iomanip>
0088 #include <sstream>
0089 
0090 #include "globals.hh"
0091 #include "G4String.hh"
0092 #include "G4Physics2DVector.hh"
0093 
0094 
0095 #include "G4MaterialTable.hh"
0096 #include "G4Material.hh"
0097 #include "G4Element.hh"
0098 #include "G4MaterialCutsCouple.hh"
0099 #include "G4ProductionCutsTable.hh"
0100 
0101 
0102 #include "G4Log.hh"
0103 #include "G4Exp.hh"
0104 
0105 
0106 class G4eDPWAElasticDCS {
0107 
0108 public:
0109 
0110   // CTR:
0111   // - iselectron   : data for e- (for e+ otherwise)
0112   // - isrestricted : sampling of angular deflection on restricted interavl is
0113   //                  required (i.e. in case of mixed-simulation models)
0114   G4eDPWAElasticDCS(G4bool iselectron=true, G4bool isrestricted=false);
0115 
0116    // DTR
0117  ~G4eDPWAElasticDCS();
0118 
0119   // initialise for a given 'iz' atomic number:
0120   //  - nothing happens if it has already been initialised for that Z.
0121   void InitialiseForZ(std::size_t iz);
0122 
0123   // Computes the elastic, first and second cross sections for the given kinetic
0124   // energy and target atom.
0125   // Cross sections are zero ff ekin is below/above the kinetic energy grid
0126   void  ComputeCSPerAtom(G4int iz, G4double ekin, G4double& elcs, G4double& tr1cs,
0127                         G4double& tr2cs, G4double mumin=0.0, G4double mumax=1.0);
0128 
0129 
0130   // samples cos(theta) i.e. cosine of the polar angle of scattering in elastic
0131   // interaction (Coulomb scattering) of the projectile (e- or e+ depending on
0132   // fIsElectron) with kinetic energy of exp('lekin'), target atom with atomic
0133   // muber of 'iz'. See the 'SampleCosineThetaRestricted' for obtain samples on
0134   // a restricted inteval.
0135   G4double SampleCosineTheta(std::size_t iz, G4double lekin, G4double r1,
0136                              G4double r2, G4double r3);
0137 
0138   // samples cos(theta) i.e. cosine of the polar angle of scattering in elastic
0139   // interaction (Coulomb scattering) of the projectile (e- or e+ depending on
0140   // fIsElectron) with kinetic energy of exp('lekin'), target atom with atomic
0141   // muber of 'iz'.
0142   // The cosine theta will be in the [costMin, costMax] interval where costMin
0143   // corresponds to a maximum allowed polar scattering angle thetaMax while
0144   // costMin corresponds to minimum allowed polar scatterin angle thetaMin.
0145   // See the 'SampleCosineTheta' for obtain samples on the entire [-1,1] range.
0146   G4double SampleCosineThetaRestricted(std::size_t iz, G4double lekin,
0147                                        G4double r1, G4double r2,
0148                                        G4double costMax, G4double costMin);
0149 
0150   // interpolate scattering power correction form table buit at init.
0151   G4double ComputeScatteringPowerCorrection(const G4MaterialCutsCouple *matcut,
0152                                             G4double ekin);
0153 
0154   // build scattering power correction table at init.
0155   void     InitSCPCorrection(G4double lowEnergyLimit, G4double highEnergyLimit);
0156 
0157 
0158 private:
0159 
0160   // data structure to store one sampling table: combined Alias + RatIn
0161   // NOTE: when Alias is used, sampling on a resctricted interval is not possible
0162   //       However, Alias makes possible faster sampling. Alias is used in case
0163   //       of single scattering model while it's not used in case of mixed-model
0164   //       when restricted interval sampling is needed. This is controlled by
0165   //       the fIsRestrictedSamplingRequired flag (false by default).
0166   struct OneSamplingTable {
0167     OneSamplingTable () = default;
0168     void SetSize(std::size_t nx, G4bool useAlias)  {
0169       fN = nx;
0170       // Alias
0171       if (useAlias) {
0172         fW.resize(nx);
0173         fI.resize(nx);
0174       }
0175       // Ratin
0176       fCum.resize(nx);
0177       fA.resize(nx);
0178       fB.resize(nx);
0179     }
0180 
0181     // members
0182     std::size_t           fN;            // # data points
0183     G4double              fScreenParA;   // the screening parameter
0184     std::vector<G4double> fW;
0185     std::vector<G4double> fCum;
0186     std::vector<G4double> fA;
0187     std::vector<G4double> fB;
0188     std::vector<G4int>    fI;
0189   };
0190 
0191 
0192   // loads the kinetic energy and theta grids for the DCS data (first init step)
0193   // should be called only by the master
0194   void LoadGrid();
0195 
0196   // load DCS data for a given Z
0197   void LoadDCSForZ(G4int iz);
0198 
0199   // loads sampling table for the given Z over the enrgy grid
0200   void BuildSmplingTableForZ(G4int iz);
0201 
0202   G4double SampleMu(std::size_t izet, std::size_t ie, G4double r1, G4double r2);
0203 
0204   G4double FindCumValue(G4double u, const OneSamplingTable& stable,
0205                         const std::vector<G4double>& uvect);
0206 
0207   // muMin and muMax : no checks on these
0208   G4double SampleMu(std::size_t izet, std::size_t ie, G4double r1, G4double muMin,
0209                     G4double muMax);
0210 
0211   // set the DCS data directory path
0212   const G4String& FindDirectoryPath();
0213 
0214   // uncompress one data file into the input string stream
0215   void ReadCompressedFile(G4String fname, std::istringstream &iss);
0216 
0217   // compute Molier material dependent parameters
0218   void ComputeMParams(const G4Material* mat, G4double& theBc, G4double& theXc2);
0219 
0220 
0221 // members
0222 private:
0223 
0224   // indicates if the object is for mixed-simulation (single scatterin otherwise)
0225   G4bool                           fIsRestrictedSamplingRequired;
0226   // indicates if the object is for e- (for e+ otherwise)
0227   G4bool                           fIsElectron;
0228   // indicates if the ekin, mu grids has already been loaded (only once)
0229   static G4bool                    gIsGridLoaded;
0230   // data directory
0231   static G4String                  gDataDirectory;
0232   // max atomic number (Z) for which DCS has been computed (103)
0233   static constexpr std::size_t     gMaxZ = 103;
0234   // energy and theta grid(s) relaed variables: loaded from gridinfo by LoadGrid
0235   static std::size_t               gNumEnergies;
0236   static std::size_t               gIndxEnergyLim;// the energy index just above 2 [keV]
0237   static std::size_t               gNumThetas1;   // used for e- below 2 [keV]
0238   static std::size_t               gNumThetas2;   // used for e+ and for e- bove 2 [keV]
0239   static std::vector<G4double>     gTheEnergies;  // log-kinetic energy grid
0240   static std::vector<G4double>     gTheMus1;      // mu(theta) = 0.5[1-cos(theta)]
0241   static std::vector<G4double>     gTheMus2;
0242   static std::vector<G4double>     gTheU1;        // u(mu; A'=0.01) = (A'+1)mu/(mu+A')
0243   static std::vector<G4double>     gTheU2;
0244   static G4double                  gLogMinEkin;   // log(gTheEnergies[0])
0245   static G4double                  gInvDelLogEkin;// 1./log(gTheEnergies[i+1]/gTheEnergies[i])
0246   // abscissas and weights of an 8 point Gauss-Legendre quadrature
0247   // for numerical integration on [0,1]
0248   static const G4double            gXGL[8];
0249   static const G4double            gWGL[8];
0250   //
0251   std::vector<G4Physics2DVector*>  fDCS;    // log(DCS) data per Z
0252   std::vector<G4Physics2DVector*>  fDCSLow; // only for e- E < 2keV
0253   // sampling tables: only one of the followings will be utilized
0254   std::vector< std::vector<OneSamplingTable>* > fSamplingTables;
0255   //
0256   // scattering power correction: to account sub-threshold inelastic deflections
0257   const G4int                      fNumSPCEbinPerDec = 3;
0258   struct SCPCorrection {
0259     G4bool   fIsUse;               //
0260     G4double fPrCut;               // sec. e- production cut energy
0261     G4double fLEmin;               // log min energy
0262     G4double fILDel;               // inverse log delta kinetic energy
0263     std::vector<G4double> fVSCPC;  // scattering power correction vector
0264   };
0265   std::vector<SCPCorrection*>      fSCPCPerMatCuts;
0266 
0267 
0268 };
0269 
0270 #endif