Geant4
9.6.p02
Main Page
Related Pages
Modules
Namespaces
Classes
Files
File List
File Members
All
Classes
Namespaces
Files
Functions
Variables
Typedefs
Enumerations
Enumerator
Friends
Macros
Groups
Pages
geant4_9_6_p02
source
processes
electromagnetic
standard
include
G4eCoulombScatteringModel.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: G4eCoulombScatteringModel
34
//
35
// Author: Vladimir Ivanchenko
36
//
37
// Creation date: 19.02.2006
38
//
39
// Modifications:
40
// 01.08.06 V.Ivanchenko extend upper limit of table to TeV and review the
41
// logic of building - only elements from G4ElementTable
42
// 08.08.06 V.Ivanchenko build internal table in ekin scale, introduce faclim
43
// 19.08.06 V.Ivanchenko add inline function ScreeningParameter and
44
// make some members protected
45
// 09.10.07 V.Ivanchenko reorganized methods, add cut dependence in scattering off e-
46
// 09.06.08 V.Ivanchenko add SelectIsotope and sampling of the recoil ion
47
// 17.06.09 C.Consoalndi modified SetupTarget method - remove kinFactor
48
// 27.05.10 V.Ivanchenko added G4WentzelOKandVIxSection class to
49
// compute cross sections and sample scattering angle
50
//
51
//
52
// Class Description:
53
//
54
// Implementation of eCoulombScattering of pointlike charge particle
55
// on Atomic Nucleus for interval of scattering anles in Lab system
56
// thetaMin - ThetaMax, nucleus recoil is neglected.
57
// The model based on analysis of J.M.Fernandez-Varea et al.
58
// NIM B73(1993)447 originated from G.Wentzel Z.Phys. 40(1927)590 with
59
// screening parameter from H.A.Bethe Phys. Rev. 89 (1953) 1256.
60
//
61
62
// -------------------------------------------------------------------
63
//
64
65
#ifndef G4eCoulombScatteringModel_h
66
#define G4eCoulombScatteringModel_h 1
67
68
#include "
G4VEmModel.hh
"
69
#include "
globals.hh
"
70
#include "
G4MaterialCutsCouple.hh
"
71
#include "
G4WentzelOKandVIxSection.hh
"
72
73
class
G4ParticleChangeForGamma
;
74
class
G4ParticleDefinition
;
75
class
G4ParticleTable
;
76
class
G4NistManager
;
77
78
class
G4eCoulombScatteringModel
:
public
G4VEmModel
79
{
80
81
public
:
82
83
G4eCoulombScatteringModel
(
const
G4String
& nam =
"eCoulombScattering"
);
84
85
virtual
~G4eCoulombScatteringModel
();
86
87
virtual
void
Initialise
(
const
G4ParticleDefinition
*,
const
G4DataVector
&);
88
89
virtual
G4double
ComputeCrossSectionPerAtom
(
90
const
G4ParticleDefinition
*,
91
G4double
kinEnergy,
92
G4double
Z
,
93
G4double
A,
94
G4double
cut,
95
G4double
emax);
96
97
virtual
void
SampleSecondaries
(std::vector<G4DynamicParticle*>*,
98
const
G4MaterialCutsCouple
*,
99
const
G4DynamicParticle
*,
100
G4double
tmin,
101
G4double
maxEnergy);
102
103
// defines low energy limit of the model
104
inline
void
SetLowEnergyThreshold
(
G4double
val);
105
106
// user definition of low-energy threshold of recoil
107
inline
void
SetRecoilThreshold
(
G4double
eth);
108
109
protected
:
110
111
inline
void
DefineMaterial
(
const
G4MaterialCutsCouple
*);
112
113
inline
void
SetupParticle
(
const
G4ParticleDefinition
*);
114
115
private
:
116
117
// hide assignment operator
118
G4eCoulombScatteringModel
& operator=(
const
G4eCoulombScatteringModel
&
right
);
119
G4eCoulombScatteringModel
(
const
G4eCoulombScatteringModel
&);
120
121
protected
:
122
123
G4ParticleTable
*
theParticleTable
;
124
G4ParticleChangeForGamma
*
fParticleChange
;
125
G4WentzelOKandVIxSection
*
wokvi
;
126
G4NistManager
*
fNistManager
;
127
128
const
std::vector<G4double>*
pCuts
;
129
130
const
G4MaterialCutsCouple
*
currentCouple
;
131
const
G4Material
*
currentMaterial
;
132
G4int
currentMaterialIndex
;
133
134
G4double
cosThetaMin
;
135
G4double
cosThetaMax
;
136
G4double
cosTetMinNuc
;
137
G4double
cosTetMaxNuc
;
138
G4double
recoilThreshold
;
139
G4double
elecRatio
;
140
G4double
mass
;
141
142
// projectile
143
const
G4ParticleDefinition
*
particle
;
144
const
G4ParticleDefinition
*
theProton
;
145
146
G4double
lowEnergyThreshold
;
147
148
private
:
149
150
G4bool
isInitialised;
151
};
152
153
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
154
155
inline
156
void
G4eCoulombScatteringModel::DefineMaterial
(
const
G4MaterialCutsCouple
* cup)
157
{
158
if
(cup !=
currentCouple
) {
159
currentCouple
= cup;
160
currentMaterial
= cup->
GetMaterial
();
161
currentMaterialIndex
=
currentCouple
->
GetIndex
();
162
}
163
}
164
165
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
166
167
inline
168
void
G4eCoulombScatteringModel::SetupParticle
(
const
G4ParticleDefinition
*
p
)
169
{
170
// Initialise mass and charge
171
if
(p !=
particle
) {
172
particle
=
p
;
173
mass
=
particle
->
GetPDGMass
();
174
wokvi
->
SetupParticle
(p);
175
}
176
}
177
178
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
179
180
inline
void
G4eCoulombScatteringModel::SetLowEnergyThreshold
(
G4double
val)
181
{
182
lowEnergyThreshold
= val;
183
}
184
185
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
186
187
inline
void
G4eCoulombScatteringModel::SetRecoilThreshold
(
G4double
eth)
188
{
189
recoilThreshold
= eth;
190
}
191
192
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
193
194
#endif
Generated on Sat May 25 2013 14:33:36 for Geant4 by
1.8.4