Geant4-11
G4hCoulombScatteringModel.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: G4hCoulombScatteringModel
33//
34// Author: Vladimir Ivanchenko
35//
36// Creation date: 08.06.2012 from G4eCoulombScatteringModel
37//
38// Modifications:
39//
40// Class Description:
41//
42// Implementation of Coulomb Scattering of a charge particle
43// on Atomic Nucleus for interval of scattering anles in Lab system
44// thetaMin - ThetaMax.
45// The model based on analysis of J.M.Fernandez-Varea et al.
46// NIM B73(1993)447 originated from G.Wentzel Z.Phys. 40(1927)590 with
47// screening parameter from H.A.Bethe Phys. Rev. 89 (1953) 1256.
48//
49
50// -------------------------------------------------------------------
51//
52
53#ifndef G4hCoulombScatteringModel_h
54#define G4hCoulombScatteringModel_h 1
55
56#include "G4VEmModel.hh"
57#include "globals.hh"
60
63class G4IonTable;
64class G4NistManager;
65
67{
68
69public:
70
71 explicit G4hCoulombScatteringModel(G4bool combined = true);
72
74
75 void Initialise(const G4ParticleDefinition*, const G4DataVector&) override;
76
78 G4VEmModel* masterModel) override;
79
82 G4double kinEnergy,
83 G4double Z,
84 G4double A,
85 G4double cut,
86 G4double emax) override;
87
88 void SampleSecondaries(std::vector<G4DynamicParticle*>*,
90 const G4DynamicParticle*,
91 G4double tmin,
92 G4double maxEnergy) override;
93
96 G4double) final;
97
98 // defines low energy limit of the model
100
101 // user definition of low-energy threshold of recoil
102 inline void SetRecoilThreshold(G4double eth);
103
104 // defines low energy limit on energy transfer to atomic electron
105 inline void SetFixedCut(G4double);
106
107 // low energy limit on energy transfer to atomic electron
108 inline G4double GetFixedCut() const;
109
110 // hide assignment operator
112 (const G4hCoulombScatteringModel &right) = delete;
114
115protected:
116
117 inline void DefineMaterial(const G4MaterialCutsCouple*);
118
119 inline void SetupParticle(const G4ParticleDefinition*);
120
121private:
122
127
130 const std::vector<G4double>* pCuts;
131
135
141
143
145};
146
147//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
148
149inline
151{
152 if(cup != currentCouple) {
153 currentCouple = cup;
156 }
157}
158
159//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
160
161inline
163{
164 // Initialise mass and charge
165 if(p != particle) {
166 particle = p;
169 }
170}
171
172//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
173
175{
176 recoilThreshold = eth;
177}
178
179//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
180
182{
183 fixedCut = val;
184}
185
186//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
187
189{
190 return fixedCut;
191}
192
193//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
194
195#endif
static const G4double emax
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 SetupParticle(const G4ParticleDefinition *)
G4hCoulombScatteringModel(G4bool combined=true)
const G4MaterialCutsCouple * currentCouple
G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A, G4double cut, G4double emax) override
void SetLowEnergyThreshold(G4double val)
const G4ParticleDefinition * theProton
void DefineMaterial(const G4MaterialCutsCouple *)
void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
const std::vector< G4double > * pCuts
const G4ParticleDefinition * particle
G4ParticleChangeForGamma * fParticleChange
G4hCoulombScatteringModel(const G4hCoulombScatteringModel &)=delete
G4double MinPrimaryEnergy(const G4Material *, const G4ParticleDefinition *, G4double) final
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel) override