Geant4  10.01.p02
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: G4eCoulombScatteringModel.hh 80656 2014-05-06 08:31:39Z gcosmo $
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"
70 #include "G4MaterialCutsCouple.hh"
72 
75 class G4IonTable;
76 class G4NistManager;
77 
79 {
80 
81 public:
82 
83  G4eCoulombScatteringModel(G4bool combined = true);
84 
86 
87  virtual void Initialise(const G4ParticleDefinition*, const G4DataVector&);
88 
89  virtual void InitialiseLocal(const G4ParticleDefinition*,
90  G4VEmModel* masterModel);
91 
93  const G4ParticleDefinition*,
94  G4double kinEnergy,
95  G4double Z,
96  G4double A,
97  G4double cut,
98  G4double emax);
99 
100  virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
101  const G4MaterialCutsCouple*,
102  const G4DynamicParticle*,
103  G4double tmin,
104  G4double maxEnergy);
105 
106  virtual G4double MinPrimaryEnergy(const G4Material*,
107  const G4ParticleDefinition*,
108  G4double);
109 
110  // defines low energy limit of the model
111  inline void SetLowEnergyThreshold(G4double val);
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 
122 protected:
123 
124  inline void DefineMaterial(const G4MaterialCutsCouple*);
125 
126  inline void SetupParticle(const G4ParticleDefinition*);
127 
128 private:
129 
130  // hide assignment operator
133 
134  //protected:
135 
140 
141  const std::vector<G4double>* pCuts;
142 
146 
152 
154 
155  // projectile
158 
160 
162 
163  //private:
164 
166 };
167 
168 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
169 
170 inline
172 {
173  if(cup != currentCouple) {
174  currentCouple = cup;
175  currentMaterial = cup->GetMaterial();
177  }
178 }
179 
180 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
181 
182 inline
184 {
185  // Initialise mass and charge
186  if(p != particle) {
187  particle = p;
188  mass = particle->GetPDGMass();
189  wokvi->SetupParticle(p);
190  }
191 }
192 
193 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
194 
196 {
197  lowEnergyThreshold = val;
198 }
199 
200 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
201 
203 {
204  recoilThreshold = eth;
205 }
206 
207 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
208 
210 {
211  fixedCut = val;
212 }
213 
214 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
215 
217 {
218  return fixedCut;
219 }
220 
221 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
222 
223 #endif
void SetupParticle(const G4ParticleDefinition *)
const G4ParticleDefinition * particle
const std::vector< G4double > * pCuts
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
const G4MaterialCutsCouple * currentCouple
int G4int
Definition: G4Types.hh:78
bool G4bool
Definition: G4Types.hh:79
void DefineMaterial(const G4MaterialCutsCouple *)
G4eCoulombScatteringModel(G4bool combined=true)
static const G4double A[nN]
virtual G4double MinPrimaryEnergy(const G4Material *, const G4ParticleDefinition *, G4double)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
G4double GetPDGMass() const
const G4ParticleDefinition * theProton
virtual void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel)
G4eCoulombScatteringModel & operator=(const G4eCoulombScatteringModel &right)
double G4double
Definition: G4Types.hh:76
G4WentzelOKandVIxSection * wokvi
void SetupParticle(const G4ParticleDefinition *)
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A, G4double cut, G4double emax)
G4ParticleChangeForGamma * fParticleChange
const G4Material * GetMaterial() const