Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
G4BetheHeitlerModel.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//
28// -------------------------------------------------------------------
29//
30// GEANT4 Class header file
31//
32//
33// File name: G4BetheHeitlerModel
34//
35// Author: Vladimir Ivanchenko on base of Michel Maire code
36//
37// Creation date: 19.04.2005
38//
39// Modifications:
40// 02-02-06 Remove InitialiseCrossSectionPerAtom();
41//
42// Class Description:
43//
44// Implementation of gamma convertion to e+e- in the field of a nucleus
45//
46
47// -------------------------------------------------------------------
48//
49
50#ifndef G4BetheHeitlerModel_h
51#define G4BetheHeitlerModel_h 1
52
53#include "G4VEmModel.hh"
54#include "G4PhysicsTable.hh"
55
57
59{
60
61public:
62
64 const G4String& nam = "BetheHeitler");
65
66 virtual ~G4BetheHeitlerModel();
67
68 virtual void Initialise(const G4ParticleDefinition*, const G4DataVector&);
69
72 G4double kinEnergy,
73 G4double Z,
74 G4double A=0.,
75 G4double cut=0.,
76 G4double emax=DBL_MAX);
77
78 virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
80 const G4DynamicParticle*,
81 G4double tmin,
82 G4double maxEnergy);
83
84private:
85
86 G4double ScreenFunction1(G4double ScreenVariable);
87
88 G4double ScreenFunction2(G4double ScreenVariable);
89
90 // hide assignment operator
91 G4BetheHeitlerModel & operator=(const G4BetheHeitlerModel &right);
93
94 G4ParticleDefinition* theGamma;
95 G4ParticleDefinition* theElectron;
96 G4ParticleDefinition* thePositron;
97 G4ParticleChangeForGamma* fParticleChange;
98};
99
100//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
101
102inline G4double G4BetheHeitlerModel::ScreenFunction1(G4double ScreenVariable)
103
104// compute the value of the screening function 3*PHI1 - PHI2
105
106{
107 G4double screenVal;
108
109 if (ScreenVariable > 1.)
110 screenVal = 42.24 - 8.368*std::log(ScreenVariable+0.952);
111 else
112 screenVal = 42.392 - ScreenVariable*(7.796 - 1.961*ScreenVariable);
113
114 return screenVal;
115}
116
117//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
118
119inline G4double G4BetheHeitlerModel::ScreenFunction2(G4double ScreenVariable)
120
121// compute the value of the screening function 1.5*PHI1 - 0.5*PHI2
122
123{
124 G4double screenVal;
125
126 if (ScreenVariable > 1.)
127 screenVal = 42.24 - 8.368*std::log(ScreenVariable+0.952);
128 else
129 screenVal = 41.405 - ScreenVariable*(5.828 - 0.8945*ScreenVariable);
130
131 return screenVal;
132}
133
134//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
135
136#endif
double G4double
Definition: G4Types.hh:64
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cut=0., G4double emax=DBL_MAX)
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
#define DBL_MAX
Definition: templates.hh:83