Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4eBremParametrizedModel.hh
Go to the documentation of this file.
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 // $Id: G4eBremParametrizedModel.hh 98737 2016-08-09 12:51:38Z gcosmo $
27 // GEANT4 tag $Name: geant4-09-04 $
28 //
29 // -------------------------------------------------------------------
30 //
31 // GEANT4 Class header file
32 //
33 //
34 // File name: G4eBremParametrizedModel
35 // extention of standard G4eBremsstrahlungModel
36 //
37 // Author: Andreas Schaelicke
38 //
39 // Creation date: 28.03.2008
40 //
41 // Modifications:
42 //
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 G4eBremParametrizedModel_h
53 #define G4eBremParametrizedModel_h 1
54 
55 #include "G4VEmModel.hh"
56 #include "G4NistManager.hh"
57 
59 class G4PhysicsVector;
60 
62 {
63 
64 public:
65 
66  explicit G4eBremParametrizedModel(const G4ParticleDefinition* p = nullptr,
67  const G4String& nam = "eBremParam");
68 
69  virtual ~G4eBremParametrizedModel();
70 
71  virtual void Initialise(const G4ParticleDefinition*, const G4DataVector&) override;
72 
73  virtual void InitialiseLocal(const G4ParticleDefinition*,
74  G4VEmModel* masterModel) override;
75 
77  const G4MaterialCutsCouple*) override;
78 
80  const G4ParticleDefinition*,
81  G4double kineticEnergy,
82  G4double cutEnergy) override;
83 
85  G4double tkin,
87  G4double cutEnergy,
88  G4double maxEnergy = DBL_MAX) override;
89 
90  virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
91  const G4MaterialCutsCouple*,
92  const G4DynamicParticle*,
93  G4double cutEnergy,
94  G4double maxEnergy) override;
95 
96  virtual void SetupForMaterial(const G4ParticleDefinition*,
97  const G4Material*,G4double) override;
98 
99 
100 private:
101 
102  void InitialiseConstants();
103 
104  G4double ComputeBremLoss(G4double cutEnergy);
105 
106  G4double ComputeXSectionPerAtom(G4double cutEnergy);
107 
108  G4double ComputeDXSectionPerAtom(G4double gammaEnergy);
109 
110  G4double ComputeParametrizedDXSectionPerAtom(G4double kineticEnergy,
111  G4double gammaEnergy,
112  G4double Z);
113 
114  void SetParticle(const G4ParticleDefinition* p);
115 
116  G4double ScreenFunction1(G4double ScreenVariable);
117 
118  G4double ScreenFunction2(G4double ScreenVariable);
119 
120  // * fast inline functions *
121  inline void SetCurrentElement(const G4double);
122 
123  // hide assignment operator
124  G4eBremParametrizedModel & operator=(const G4eBremParametrizedModel &right) = delete;
126 
127 protected:
128 
133 
134  static const G4double xgi[8], wgi[8];
135 
137 
138  // cash
149 
151 
152 private:
153 
154  // consts
155  G4double lowKinEnergy;
156  G4double fMigdalConstant;
157  G4double bremFactor;
158 
159  G4bool isInitialised;
160 };
161 
162 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
163 
164 inline void G4eBremParametrizedModel::SetCurrentElement(const G4double Z)
165 {
166  //std::cout<<"SetCurrentElement Z="<<Z<<std::endl;
167  if(Z != currentZ) {
168  currentZ = Z;
169 
170  G4int iz = G4int(Z);
171  z13 = nist->GetZ13(iz);
172  z23 = z13*z13;
173  lnZ = nist->GetLOGZ(iz);
174 
175  Fel = facFel - lnZ/3. ;
176  Finel = facFinel - 2.*lnZ/3. ;
177 
179  fMax = Fel-fCoulomb + Finel/currentZ + (1.+1./currentZ)/12.;
180 
181  }
182 }
183 
184 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
185 
186 
187 #endif
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double tkin, G4double Z, G4double, G4double cutEnergy, G4double maxEnergy=DBL_MAX) override
G4ParticleChangeForLoss * fParticleChange
const char * p
Definition: xmltok.h:285
G4double GetfCoulomb() const
Definition: G4Element.hh:191
static const G4double xgi[8]
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
virtual void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel) override
int G4int
Definition: G4Types.hh:78
G4double GetZ13(G4double Z) const
bool G4bool
Definition: G4Types.hh:79
G4eBremParametrizedModel(const G4ParticleDefinition *p=nullptr, const G4String &nam="eBremParam")
G4double GetLOGZ(G4int Z) const
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy) override
virtual G4double MinEnergyCut(const G4ParticleDefinition *, const G4MaterialCutsCouple *) override
virtual void SetupForMaterial(const G4ParticleDefinition *, const G4Material *, G4double) override
const G4ParticleDefinition * particle
double G4double
Definition: G4Types.hh:76
static const G4double wgi[8]
#define DBL_MAX
Definition: templates.hh:83
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double cutEnergy, G4double maxEnergy) override
const G4Element * GetCurrentElement() const
Definition: G4VEmModel.hh:466