Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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 void Initialise(const G4ParticleDefinition*, const G4DataVector&) override;
87
89 G4VEmModel* masterModel) override;
90
93 G4double kinEnergy,
94 G4double Z,
95 G4double A,
96 G4double cut,
97 G4double emax) override;
98
99 void SampleSecondaries(std::vector<G4DynamicParticle*>*,
101 const G4DynamicParticle*,
102 G4double tmin,
103 G4double maxEnergy) override;
104
106 G4double) final;
107
108 // defines low energy limit of the model
110
111 // user definition of low-energy threshold of recoil
112 inline void SetRecoilThreshold(G4double eth);
113
114 // defines low energy limit on energy transfer to atomic electron
115 inline void SetFixedCut(G4double);
116
117 // low energy limit on energy transfer to atomic electron
118 inline G4double GetFixedCut() const;
119
120 // hide assignment operator
122 (const G4eCoulombScatteringModel &right) = delete;
124
125protected:
126
127 inline void DefineMaterial(const G4MaterialCutsCouple*);
128
129 inline void SetupParticle(const G4ParticleDefinition*);
130
131private:
132
133 G4IonTable* theIonTable;
134 G4ParticleChangeForGamma* fParticleChange;
136 G4NistManager* fNistManager;
137
138 const std::vector<G4double>* pCuts;
139
140 const G4ParticleDefinition* particle;
141 const G4ParticleDefinition* theProton;
142
143 const G4MaterialCutsCouple* currentCouple;
144 const G4Material* currentMaterial;
145 G4int currentMaterialIndex;
146
147 G4double cosThetaMin;
148 G4double cosThetaMax;
149 G4double recoilThreshold;
150 G4double elecRatio;
151 G4double mass;
152 G4double fixedCut;
153
154 G4bool isCombined;
155};
156
157//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
158
159inline
161{
162 if(cup != currentCouple) {
163 currentCouple = cup;
164 currentMaterial = cup->GetMaterial();
165 currentMaterialIndex = currentCouple->GetIndex();
166 }
167}
168
169//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
170
171inline
173{
174 // Initialise mass and charge
175 if(p != particle) {
176 particle = p;
177 mass = particle->GetPDGMass();
178 wokvi->SetupParticle(p);
179 }
180}
181
182//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
183
185{
186 recoilThreshold = eth;
187}
188
189//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
190
192{
193 fixedCut = val;
194}
195
196//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
197
199{
200 return fixedCut;
201}
202
203//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
204
205#endif
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
const G4double A[17]
const G4Material * GetMaterial() const
void SetupParticle(const G4ParticleDefinition *)
void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A, G4double cut, G4double emax) override
G4double MinPrimaryEnergy(const G4Material *, const G4ParticleDefinition *, G4double) final
void SetLowEnergyThreshold(G4double val)
void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel) override
G4eCoulombScatteringModel(const G4eCoulombScatteringModel &)=delete
void DefineMaterial(const G4MaterialCutsCouple *)
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
void SetupParticle(const G4ParticleDefinition *)