Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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$
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
59class G4PhysicsVector;
60
62{
63
64public:
65
67 const G4String& nam = "eBremParam");
68
70
71 virtual void Initialise(const G4ParticleDefinition*, const G4DataVector&);
72
75
78 G4double kineticEnergy,
79 G4double cutEnergy);
80
82 G4double tkin,
84 G4double cutEnergy,
85 G4double maxEnergy = DBL_MAX);
86
87 virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
89 const G4DynamicParticle*,
90 G4double cutEnergy,
91 G4double maxEnergy);
92
93 virtual void SetupForMaterial(const G4ParticleDefinition*,
94 const G4Material*,G4double);
95
96
97private:
98
99 void InitialiseConstants();
100
101 G4double ComputeBremLoss(G4double cutEnergy);
102
103 G4double ComputeXSectionPerAtom(G4double cutEnergy);
104
105 G4double ComputeDXSectionPerAtom(G4double gammaEnergy);
106
107 void SetParticle(const G4ParticleDefinition* p);
108
109 // * fast inline functions *
110 inline void SetCurrentElement(const G4double);
111
112 // hide assignment operator
113 G4eBremParametrizedModel & operator=(const G4eBremParametrizedModel &right);
115
116protected:
117
122
123 static const G4double xgi[8], wgi[8];
124
126
127 // cash
138
140
141private:
142
143 // consts
144 G4double lowKinEnergy;
145 G4double fMigdalConstant;
146 G4double bremFactor;
147
148 G4bool isInitialised;
149};
150
151//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
152
153inline void G4eBremParametrizedModel::SetCurrentElement(const G4double Z)
154{
155 std::cout<<"SetCurrentElement Z="<<Z<<std::endl;
156 if(Z != currentZ) {
157 currentZ = Z;
158
159 G4int iz = G4int(Z);
160 z13 = nist->GetZ13(iz);
161 z23 = z13*z13;
162 lnZ = nist->GetLOGZ(iz);
163
164 Fel = facFel - lnZ/3. ;
165 Finel = facFinel - 2.*lnZ/3. ;
166
168 fMax = Fel-fCoulomb + Finel/currentZ + (1.+1./currentZ)/12.;
169
170 }
171}
172
173//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
174
175
176#endif
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
G4double GetfCoulomb() const
Definition: G4Element.hh:201
G4double GetZ13(G4double Z)
G4double GetLOGZ(G4int Z)
const G4Element * GetCurrentElement() const
Definition: G4VEmModel.hh:391
virtual void SetupForMaterial(const G4ParticleDefinition *, const G4Material *, G4double)
static const G4double wgi[8]
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy)
G4ParticleChangeForLoss * fParticleChange
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double tkin, G4double Z, G4double, G4double cutEnergy, G4double maxEnergy=DBL_MAX)
const G4ParticleDefinition * particle
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double cutEnergy, G4double maxEnergy)
virtual G4double MinEnergyCut(const G4ParticleDefinition *, const G4MaterialCutsCouple *)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
static const G4double xgi[8]
#define DBL_MAX
Definition: templates.hh:83