Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/electromagnetic/standard/include/G4eBremsstrahlungRelModel.hh

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

  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 // -------------------------------------------------------------------
 28 //
 29 // GEANT4 Class header file
 30 //
 31 //
 32 // File name:     G4eBremsstrahlungRelModel
 33 //                extention of standard G4eBremsstrahlungModel
 34 //
 35 // Author:        Andreas Schaelicke
 36 //
 37 // Creation date: 28.03.2008
 38 //
 39 // Modifications:
 40 //
 41 // 15.07.18  introduced data structures to store LPM functions and element depen
 42 //           dent data for faster run-time computation (see more in .cc M.Novak)
 43 //
 44 // Class Description:
 45 //
 46 // Implementation of energy loss for gamma emission by electrons and
 47 // positrons including an improved version of the LPM effect
 48 
 49 // -------------------------------------------------------------------
 50 //
 51 
 52 #ifndef G4eBremsstrahlungRelModel_h
 53 #define G4eBremsstrahlungRelModel_h 1
 54 
 55 #include "G4VEmModel.hh"
 56 #include <memory>
 57 
 58 class G4ParticleChangeForLoss;
 59 
 60 class G4eBremsstrahlungRelModel : public G4VEmModel {
 61 
 62 public:
 63 
 64   explicit G4eBremsstrahlungRelModel(const G4ParticleDefinition* p=nullptr,
 65                                      const G4String& nam="eBremLPM");
 66 
 67   ~G4eBremsstrahlungRelModel() override;
 68 
 69   void Initialise(const G4ParticleDefinition*, const G4DataVector&) override;
 70 
 71   void InitialiseLocal(const G4ParticleDefinition*,
 72                        G4VEmModel* masterModel) override;
 73 
 74   G4double ComputeDEDXPerVolume(const G4Material*,
 75                                 const G4ParticleDefinition*,
 76                                 G4double ekin,
 77                                 G4double cutEnergy) override;
 78 
 79   G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition*,
 80                                       G4double ekin,
 81                                       G4double zet,
 82                                       G4double,
 83                                       G4double cutEnergy,
 84                                       G4double maxEnergy = DBL_MAX) override;
 85 
 86   void SampleSecondaries(std::vector<G4DynamicParticle*>*,
 87                          const G4MaterialCutsCouple*,
 88                          const G4DynamicParticle*,
 89                          G4double cutEnergy,
 90                          G4double maxEnergy) override;
 91 
 92   void SetupForMaterial(const G4ParticleDefinition*,
 93                         const G4Material*, G4double) override;
 94 
 95   G4double MinPrimaryEnergy(const G4Material*,
 96                             const G4ParticleDefinition*,
 97                             G4double cutEnergy) override;
 98 
 99 protected:
100 
101   virtual G4double ComputeDXSectionPerAtom(G4double gammaEnergy);
102 
103   void SetParticle(const G4ParticleDefinition* p);
104 
105 private:
106 
107   G4double ComputeBremLoss(G4double cutEnergy);
108 
109   G4double ComputeXSectionPerAtom(G4double cutEnergy);
110 
111   G4double ComputeRelDXSectionPerAtom(G4double gammaEnergy);
112 
113   // init special data per element i.e. per Z
114   void InitialiseElementData();
115 
116   // methods for initialisation and run-time evaluation of LPM functions:
117   void InitLPMFunctions();
118 
119   void ComputeLPMfunctions(G4double& funcXiS,
120                            G4double& funcGS,
121                            G4double& funcPhiS,
122                            const G4double egamma);
123 
124   void GetLPMFunctions(G4double& lpmGs,
125                        G4double& lpmPhis,
126                        const G4double ss);
127 
128   void ComputeLPMGsPhis(G4double& funcGS,
129                         G4double& funcPhiS,
130                         const G4double varShat);
131 
132   // for evaluating screening related functions
133   void ComputeScreeningFunctions(G4double& phi1,
134                                  G4double& phi1m2,
135                                  G4double& psi1,
136                                  G4double& psi1m2,
137                                  const G4double gam,
138                                  const G4double eps);
139   // hide assignment operator and cctr
140   G4eBremsstrahlungRelModel& operator=
141   (const G4eBremsstrahlungRelModel& right) = delete;
142   G4eBremsstrahlungRelModel(const  G4eBremsstrahlungRelModel&) = delete;
143 
144 private:
145 
146   G4bool                      fIsUseCompleteScreening = false;
147   G4bool                      fIsInitializer = false;
148   G4bool                      fUseLPM = true;
149 
150 protected:
151 
152   G4bool                      fIsElectron = true;
153   G4bool                      fIsScatOffElectron = false;
154   G4bool                      fIsLPMActive = false;
155   //
156   G4int                       fCurrentIZ = 0;
157   const G4ParticleDefinition* fPrimaryParticle = nullptr;
158   G4ParticleDefinition*       fGammaParticle = nullptr;
159   G4ParticleChangeForLoss*    fParticleChange = nullptr;
160   // cash
161   G4double                    fPrimaryParticleMass = 0.;
162   G4double                    fPrimaryKinEnergy = 0.;
163   G4double                    fPrimaryTotalEnergy = 0.;
164   G4double                    fDensityFactor = 0.;
165   G4double                    fDensityCorr = 0.;
166   G4double                    fLowestKinEnergy;
167   // scattering off electrons
168   G4double                    fNucTerm = 0.;
169   G4double                    fSumTerm = 0.;
170 
171 private:
172 
173   // LPM related members
174   G4double                    fLPMEnergyThreshold;
175   G4double                    fLPMEnergy;
176 
177 protected:
178 
179   static const G4double       gBremFactor;
180   static const G4double       gMigdalConstant;
181 
182 private:
183 
184   static const G4int          gMaxZet;
185   //
186   static const G4double       gLPMconstant;
187   //
188   static const G4double       gXGL[8];
189   static const G4double       gWGL[8];
190   static const G4double       gFelLowZet[8];
191   static const G4double       gFinelLowZet[8];
192   //
193   struct ElementData {
194     /** @brief \f$ \ln(Z) \f$  */
195     G4double  fLogZ;
196     /** @brief \f$ \ln(Z)/3 + f_c \f$  */
197     G4double  fFz;
198     /** @brief \f$ ((Fel-fc)+Finel*invZ)\f$  */
199     G4double  fZFactor1;
200     /** @brief \f$ (Fel-fc)\f$  */
201     G4double  fZFactor11;
202     /** @brief \f$ (1.0+invZ)/12  \f$  */
203     G4double  fZFactor2;
204     // LPM variables
205     G4double  fVarS1;
206     G4double  fILVarS1;
207     G4double  fILVarS1Cond;
208     // constant factors to the screening function evaluations
209     G4double  fGammaFactor;
210     G4double  fEpsilonFactor;
211   };
212   //
213   struct LPMFuncs {
214     LPMFuncs() : fIsInitialized(false), fISDelta(100.), fSLimit(2.) {}
215     G4bool                 fIsInitialized;
216     G4double               fISDelta;
217     G4double               fSLimit;
218     std::vector<G4double>  fLPMFuncG;
219     std::vector<G4double>  fLPMFuncPhi;
220   };
221   //
222   static std::shared_ptr<LPMFuncs> gLPMFuncs();
223   static std::shared_ptr<std::vector<ElementData*>> gElementData();
224   std::shared_ptr<LPMFuncs> fLPMFuncs;
225   std::shared_ptr<std::vector<ElementData*>> fElementData;
226 };
227 
228 #endif
229