Geant4 10.7.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//
27// -------------------------------------------------------------------
28//
29// GEANT4 Class header file
30//
31//
32// File name: G4eCoulombScatteringModel
33//
34// Author: Vladimir Ivanchenko
35//
36// Creation date: 19.02.2006
37//
38// Modifications:
39// 01.08.06 V.Ivanchenko extend upper limit of table to TeV and review the
40// logic of building - only elements from G4ElementTable
41// 08.08.06 V.Ivanchenko build internal table in ekin scale, introduce faclim
42// 19.08.06 V.Ivanchenko add inline function ScreeningParameter and
43// make some members protected
44// 09.10.07 V.Ivanchenko reorganized methods, add cut dependence in scattering off e-
45// 09.06.08 V.Ivanchenko add SelectIsotope and sampling of the recoil ion
46// 17.06.09 C.Consoalndi modified SetupTarget method - remove kinFactor
47// 27.05.10 V.Ivanchenko added G4WentzelOKandVIxSection class to
48// compute cross sections and sample scattering angle
49//
50//
51// Class Description:
52//
53// Implementation of eCoulombScattering of a charge particle
54// on Atomic Nucleus for interval of scattering anles in Lab system
55// thetaMin - ThetaMax.
56// The model based on analysis of J.M.Fernandez-Varea et al.
57// NIM B73(1993)447 originated from G.Wentzel Z.Phys. 40(1927)590 with
58// screening parameter from H.A.Bethe Phys. Rev. 89 (1953) 1256.
59//
60
61// -------------------------------------------------------------------
62//
63
64#ifndef G4eCoulombScatteringModel_h
65#define G4eCoulombScatteringModel_h 1
66
67#include "G4VEmModel.hh"
68#include "globals.hh"
71
74class G4IonTable;
75class G4NistManager;
76
78{
79
80public:
81
82 explicit G4eCoulombScatteringModel(G4bool combined = true);
83
85
86 virtual void Initialise(const G4ParticleDefinition*,
87 const G4DataVector&) override;
88
89 virtual void InitialiseLocal(const G4ParticleDefinition*,
90 G4VEmModel* masterModel) override;
91
94 G4double kinEnergy,
95 G4double Z,
96 G4double A,
97 G4double cut,
98 G4double emax) override;
99
100 virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
102 const G4DynamicParticle*,
103 G4double tmin,
104 G4double maxEnergy) override;
105
106 virtual G4double MinPrimaryEnergy(const G4Material*,
108 G4double) final;
109
110 // defines low energy limit of the model
112
113 // user definition of low-energy threshold of recoil
114 inline void SetRecoilThreshold(G4double eth);
115
116 // defines low energy limit on energy transfer to atomic electron
117 inline void SetFixedCut(G4double);
118
119 // low energy limit on energy transfer to atomic electron
120 inline G4double GetFixedCut() const;
121
122protected:
123
124 inline void DefineMaterial(const G4MaterialCutsCouple*);
125
126 inline void SetupParticle(const G4ParticleDefinition*);
127
128private:
129
130 // hide assignment operator
131 G4eCoulombScatteringModel & operator=
132 (const G4eCoulombScatteringModel &right) = delete;
134
135 //protected:
136
137 G4IonTable* theIonTable;
138 G4ParticleChangeForGamma* fParticleChange;
140 G4NistManager* fNistManager;
141
142 const std::vector<G4double>* pCuts;
143
144 const G4MaterialCutsCouple* currentCouple;
145 const G4Material* currentMaterial;
146 G4int currentMaterialIndex;
147
148 G4double cosThetaMin;
149 G4double cosThetaMax;
150 G4double recoilThreshold;
151 G4double elecRatio;
152 G4double mass;
153
154 G4double fixedCut;
155
156 // projectile
157 const G4ParticleDefinition* particle;
158 const G4ParticleDefinition* theProton;
159
160 G4bool isCombined;
161};
162
163//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
164
165inline
167{
168 if(cup != currentCouple) {
169 currentCouple = cup;
170 currentMaterial = cup->GetMaterial();
171 currentMaterialIndex = currentCouple->GetIndex();
172 }
173}
174
175//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
176
177inline
179{
180 // Initialise mass and charge
181 if(p != particle) {
182 particle = p;
183 mass = particle->GetPDGMass();
184 wokvi->SetupParticle(p);
185 }
186}
187
188//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
189
191{
192 recoilThreshold = eth;
193}
194
195//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
196
198{
199 fixedCut = val;
200}
201
202//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
203
205{
206 return fixedCut;
207}
208
209//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
210
211#endif
double A(double temperature)
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
const G4Material * GetMaterial() const
void SetupParticle(const G4ParticleDefinition *)
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A, G4double cut, G4double emax) override
virtual G4double MinPrimaryEnergy(const G4Material *, const G4ParticleDefinition *, G4double) final
void SetLowEnergyThreshold(G4double val)
virtual void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel) override
void DefineMaterial(const G4MaterialCutsCouple *)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
void SetupParticle(const G4ParticleDefinition *)