Geant4  10.00.p02
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 // $Id: G4hCoulombScatteringModel.hh 70370 2013-05-29 15:15:40Z gcosmo $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class header file
31 //
32 //
33 // File name: G4hCoulombScatteringModel
34 //
35 // Author: Vladimir Ivanchenko
36 //
37 // Creation date: 08.06.2012 from G4eCoulombScatteringModel
38 //
39 // Modifications:
40 //
41 // Class Description:
42 //
43 // Implementation of Coulomb Scattering of pointlike charge particle
44 // on Atomic Nucleus for interval of scattering anles in Lab system
45 // thetaMin - ThetaMax, nucleus recoil is neglected.
46 // The model based on analysis of J.M.Fernandez-Varea et al.
47 // NIM B73(1993)447 originated from G.Wentzel Z.Phys. 40(1927)590 with
48 // screening parameter from H.A.Bethe Phys. Rev. 89 (1953) 1256.
49 //
50 
51 // -------------------------------------------------------------------
52 //
53 
54 #ifndef G4hCoulombScatteringModel_h
55 #define G4hCoulombScatteringModel_h 1
56 
57 #include "G4VEmModel.hh"
58 #include "globals.hh"
59 #include "G4MaterialCutsCouple.hh"
61 
64 class G4IonTable;
65 class G4NistManager;
66 
68 {
69 
70 public:
71 
72  G4hCoulombScatteringModel(const G4String& nam = "eCoulombScattering");
73 
75 
76  virtual void Initialise(const G4ParticleDefinition*, const G4DataVector&);
77 
79  const G4ParticleDefinition*,
80  G4double kinEnergy,
81  G4double Z,
82  G4double A,
83  G4double cut,
84  G4double emax);
85 
86  virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
87  const G4MaterialCutsCouple*,
88  const G4DynamicParticle*,
89  G4double tmin,
90  G4double maxEnergy);
91 
92  // defines low energy limit of the model
93  inline void SetLowEnergyThreshold(G4double val);
94 
95  // user definition of low-energy threshold of recoil
96  inline void SetRecoilThreshold(G4double eth);
97 
98 protected:
99 
100  inline void DefineMaterial(const G4MaterialCutsCouple*);
101 
102  inline void SetupParticle(const G4ParticleDefinition*);
103 
104 private:
105 
106  // hide assignment operator
109 
110 protected:
111 
116 
117  const std::vector<G4double>* pCuts;
118 
122 
130 
131  // projectile
134 
136 
137 private:
138 
140 };
141 
142 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
143 
144 inline
146 {
147  if(cup != currentCouple) {
148  currentCouple = cup;
149  currentMaterial = cup->GetMaterial();
151  }
152 }
153 
154 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
155 
156 inline
158 {
159  // Initialise mass and charge
160  if(p != particle) {
161  particle = p;
162  mass = particle->GetPDGMass();
163  wokvi->SetupParticle(p);
164  }
165 }
166 
167 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
168 
170 {
171  lowEnergyThreshold = val;
172 }
173 
174 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
175 
177 {
178  recoilThreshold = eth;
179 }
180 
181 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
182 
183 #endif
const G4ParticleDefinition * theProton
G4hCoulombScatteringModel(const G4String &nam="eCoulombScattering")
void DefineMaterial(const G4MaterialCutsCouple *)
const std::vector< G4double > * pCuts
void SetupParticle(const G4ParticleDefinition *)
int G4int
Definition: G4Types.hh:78
const G4MaterialCutsCouple * currentCouple
void SetupParticle(const G4ParticleDefinition *)
const G4ParticleDefinition * particle
G4hCoulombScatteringModel & operator=(const G4hCoulombScatteringModel &right)
bool G4bool
Definition: G4Types.hh:79
G4ParticleChangeForGamma * fParticleChange
static const G4double A[nN]
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A, G4double cut, G4double emax)
G4double GetPDGMass() const
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
double G4double
Definition: G4Types.hh:76
const G4Material * GetMaterial() const