Back to home page

EIC code displayed by LXR

 
 

    


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

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:     G4BetheHeitlerModel
0033 //
0034 // Author:        Vladimir Ivanchenko on base of Michel Maire code
0035 //
0036 // Creation date: 19.04.2005
0037 //
0038 // Modifications by Vladimir Ivanchenko, Michel Maire, Mihaly Novak
0039 //
0040 // Class Description:
0041 //
0042 // Implementation of gamma conversion to e+e- in the field of a nucleus 
0043 // For details see Physics Reference Manual
0044 
0045 // -------------------------------------------------------------------
0046 //
0047 
0048 #ifndef G4BetheHeitlerModel_h
0049 #define G4BetheHeitlerModel_h 1
0050 
0051 #include "G4VEmModel.hh"
0052 #include "G4PhysicsTable.hh"
0053 #include "G4Log.hh"
0054 
0055 #include <vector>
0056 
0057 class G4ParticleChangeForGamma;
0058 class G4Pow;
0059 
0060 class G4BetheHeitlerModel : public G4VEmModel
0061 {
0062 
0063 public:
0064 
0065   explicit G4BetheHeitlerModel(const G4ParticleDefinition* p = nullptr,
0066                                const G4String& nam = "BetheHeitler");
0067  
0068   ~G4BetheHeitlerModel() override;
0069 
0070   void Initialise(const G4ParticleDefinition*, const G4DataVector&) override;
0071 
0072   void InitialiseLocal(const G4ParticleDefinition*, 
0073                G4VEmModel* masterModel) override;
0074 
0075   G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition*,
0076                       G4double kinEnergy, 
0077                       G4double Z, 
0078                       G4double A=0., 
0079                       G4double cut=0.,
0080                       G4double emax=DBL_MAX) override;
0081 
0082   void SampleSecondaries(std::vector<G4DynamicParticle*>*,
0083              const G4MaterialCutsCouple*,
0084              const G4DynamicParticle*,
0085              G4double tmin,
0086              G4double maxEnergy) override;
0087 
0088   // hide assignment operator
0089   G4BetheHeitlerModel & operator=(const G4BetheHeitlerModel &right) = delete;
0090   G4BetheHeitlerModel(const  G4BetheHeitlerModel&) = delete;
0091 
0092 protected:
0093 
0094   inline G4double ScreenFunction1(const G4double delta);
0095 
0096   inline G4double ScreenFunction2(const G4double delta);
0097 
0098   inline void ScreenFunction12(const G4double delta, G4double &f1, G4double &f2);
0099 
0100   void     InitialiseElementData();
0101  
0102   struct ElementData {
0103     G4double fDeltaMaxLow;
0104     G4double fDeltaMaxHigh;
0105   };
0106   
0107   static const G4int                gMaxZet; 
0108   
0109   G4Pow*                            fG4Calc;
0110   const G4ParticleDefinition*       fTheGamma;
0111   const G4ParticleDefinition*       fTheElectron;
0112   const G4ParticleDefinition*       fThePositron;
0113   G4ParticleChangeForGamma*         fParticleChange;
0114 
0115   G4bool isFirstInstance{false};
0116 
0117   static std::vector<ElementData*>  gElementData;
0118 };
0119 
0120 //
0121 // Bethe screening functions for the elastic (coherent) scattering:
0122 // Bethe's phi1, phi2 coherent screening functions were computed numerically 
0123 // by using (the universal) atomic form factors computed based on the Thomas-
0124 // Fermi model of the atom (using numerical solution of the Thomas-Fermi 
0125 // screening function instead of Moliere's analytical approximation). The 
0126 // numerical results can be well approximated (better than Butcher & Messel 
0127 // especially near the delta=1 limit) by:
0128 // ## if delta <= 1.4 
0129 //  phi1(delta) = 20.806 - delta*(3.190 - 0.5710*delta)   
0130 //  phi2(delta) = 20.234 - delta*(2.126 - 0.0903*delta)
0131 // ## if delta  > 1.4
0132 //  phi1(delta) = phi2(delta) = 21.0190 - 4.145*ln(delta + 0.958)
0133 // with delta = 136mc^2kZ^{-1/3}/[E(Eg-E)] = 136Z^{-1/3}eps0/[eps(1-eps)] where 
0134 // Eg is the initial photon energy, E is the total energy transferred to one of 
0135 // the e-/e+ pair, eps0 = mc^2/Eg and eps = E/Eg.
0136 
0137 // Compute the value of the screening function 3*PHI1(delta) - PHI2(delta):
0138 inline G4double G4BetheHeitlerModel::ScreenFunction1(const G4double delta)
0139 {
0140   return (delta > 1.4) ? 42.038 - 8.29*G4Log(delta + 0.958) 
0141                        : 42.184 - delta*(7.444 - 1.623*delta);
0142 }
0143 
0144 // Compute the value of the screening function 1.5*PHI1(delta) +0.5*PHI2(delta):
0145 inline G4double G4BetheHeitlerModel::ScreenFunction2(const G4double delta)
0146 {
0147   return (delta > 1.4) ? 42.038 - 8.29*G4Log(delta + 0.958)
0148                        : 41.326 - delta*(5.848 - 0.902*delta);
0149 }
0150 
0151 // Same as ScreenFunction1 and ScreenFunction2 but computes them at once
0152 inline void G4BetheHeitlerModel::ScreenFunction12(const G4double delta, 
0153                                                   G4double &f1, G4double &f2)
0154 {
0155   if (delta > 1.4) {
0156     f1 = 42.038 - 8.29*G4Log(delta + 0.958);
0157     f2 = f1;
0158   } else {
0159     f1 = 42.184 - delta*(7.444 - 1.623*delta);
0160     f2 = 41.326 - delta*(5.848 - 0.902*delta); 
0161   }
0162 }
0163 
0164 #endif