Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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"
72
75class G4ParticleTable;
76class G4NistManager;
77
79{
80
81public:
82
83 G4eCoulombScatteringModel(const G4String& nam = "eCoulombScattering");
84
86
87 virtual void Initialise(const G4ParticleDefinition*, const G4DataVector&);
88
91 G4double kinEnergy,
92 G4double Z,
93 G4double A,
94 G4double cut,
95 G4double emax);
96
97 virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
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
109protected:
110
111 inline void DefineMaterial(const G4MaterialCutsCouple*);
112
113 inline void SetupParticle(const G4ParticleDefinition*);
114
115private:
116
117 // hide assignment operator
120
121protected:
122
127
128 const std::vector<G4double>* pCuts;
129
133
141
142 // projectile
145
147
148private:
149
150 G4bool isInitialised;
151};
152
153//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
154
155inline
157{
158 if(cup != currentCouple) {
159 currentCouple = cup;
162 }
163}
164
165//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
166
167inline
169{
170 // Initialise mass and charge
171 if(p != particle) {
172 particle = p;
175 }
176}
177
178//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
179
181{
182 lowEnergyThreshold = val;
183}
184
185//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
186
188{
189 recoilThreshold = eth;
190}
191
192//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
193
194#endif
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
const G4Material * GetMaterial() const
void SetupParticle(const G4ParticleDefinition *)
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
G4WentzelOKandVIxSection * wokvi
const G4ParticleDefinition * particle
const G4MaterialCutsCouple * currentCouple
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
G4ParticleChangeForGamma * fParticleChange
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A, G4double cut, G4double emax)
void DefineMaterial(const G4MaterialCutsCouple *)
const std::vector< G4double > * pCuts
const G4ParticleDefinition * theProton
void SetupParticle(const G4ParticleDefinition *)