Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
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//
27// -------------------------------------------------------------------
28//
29// GEANT4 Class header file
30//
31//
32// File name: G4eBremParametrizedModel
33// extention of standard G4eBremsstrahlungModel
34//
35// Author: Andreas Schaelicke
36//
37// Creation date: 28.03.2008
38//
39// Modifications:
40//
41//
42// Class Description:
43//
44// Implementation of energy loss for gamma emission by electrons and
45// positrons including an improved version of the LPM effect
46
47// -------------------------------------------------------------------
48//
49
50#ifndef G4eBremParametrizedModel_h
51#define G4eBremParametrizedModel_h 1
52
53#include "G4VEmModel.hh"
54#include "G4NistManager.hh"
55
57class G4PhysicsVector;
58
60{
61
62public:
63
64 explicit G4eBremParametrizedModel(const G4ParticleDefinition* p = nullptr,
65 const G4String& nam = "eBremParam");
66
68
69 void Initialise(const G4ParticleDefinition*, const G4DataVector&) override;
70
72 G4VEmModel* masterModel) override;
73
75 const G4MaterialCutsCouple*) override;
76
79 G4double kineticEnergy,
80 G4double cutEnergy) override;
81
83 G4double tkin,
85 G4double cutEnergy,
86 G4double maxEnergy = DBL_MAX) override;
87
88 void SampleSecondaries(std::vector<G4DynamicParticle*>*,
90 const G4DynamicParticle*,
91 G4double cutEnergy,
92 G4double maxEnergy) override;
93
95 const G4Material*,G4double) override;
96
97 // hide assignment operator
100
101private:
102
103 void InitialiseConstants();
104
105 G4double ComputeBremLoss(G4double cutEnergy);
106
107 G4double ComputeXSectionPerAtom(G4double cutEnergy);
108
109 G4double ComputeDXSectionPerAtom(G4double gammaEnergy);
110
111 G4double ComputeParametrizedDXSectionPerAtom(G4double kineticEnergy,
112 G4double gammaEnergy,
113 G4double Z);
114
115 void SetParticle(const G4ParticleDefinition* p);
116
117 G4double ScreenFunction1(G4double ScreenVariable);
118
119 G4double ScreenFunction2(G4double ScreenVariable);
120
121 inline void SetCurrentElement(const G4double);
122
123protected:
124
129
130 static const G4double xgi[8], wgi[8];
131
133
134 // cash
145
146private:
147
148 G4double lowKinEnergy;
149 G4double fMigdalConstant;
150 G4double bremFactor;
151
152 G4bool isInitialised;
153
154protected:
155
157
158};
159
160//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
161
162inline void G4eBremParametrizedModel::SetCurrentElement(const G4double Z)
163{
164 if(Z != currentZ) {
165 currentZ = Z;
166
167 G4int iz = G4lrint(Z);
168 z13 = nist->GetZ13(iz);
169 z23 = z13*z13;
170 lnZ = nist->GetLOGZ(iz);
171
172 Fel = facFel - lnZ/3. ;
173 Finel = facFinel - 2.*lnZ/3. ;
174
176 fMax = Fel-fCoulomb + Finel/currentZ + (1.+1./currentZ)/12.;
177 }
178}
179
180//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
181
182
183#endif
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
G4double GetfCoulomb() const
Definition: G4Element.hh:190
G4double GetZ13(G4double Z) const
G4double GetLOGZ(G4int Z) const
const G4Element * GetCurrentElement(const G4Material *mat=nullptr) const
Definition: G4VEmModel.cc:242
void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel) override
static const G4double wgi[8]
~G4eBremParametrizedModel() override
G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy) override
G4double MinEnergyCut(const G4ParticleDefinition *, const G4MaterialCutsCouple *) override
G4ParticleChangeForLoss * fParticleChange
G4eBremParametrizedModel & operator=(const G4eBremParametrizedModel &right)=delete
const G4ParticleDefinition * particle
void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double cutEnergy, G4double maxEnergy) override
G4eBremParametrizedModel(const G4eBremParametrizedModel &)=delete
G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double tkin, G4double Z, G4double, G4double cutEnergy, G4double maxEnergy=DBL_MAX) override
static const G4double xgi[8]
void SetupForMaterial(const G4ParticleDefinition *, const G4Material *, G4double) override
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
int G4lrint(double ad)
Definition: templates.hh:134
#define DBL_MAX
Definition: templates.hh:62