Geant4  10.02.p01
G4WentzelOKandVIxSection.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: G4WentzelOKandVIxSection.hh 80656 2014-05-06 08:31:39Z gcosmo $
27 //
28 // -------------------------------------------------------------------
29 //
30 //
31 // GEANT4 Class header file
32 //
33 //
34 // File name: G4WentzelOKandVIxSection
35 //
36 // Authors: V.Ivanchenko and O.Kadri
37 //
38 // Creation date: 21.05.2010
39 //
40 // Modifications:
41 //
42 //
43 // Class Description:
44 //
45 // Implementation of the computation of total and transport cross sections,
46 // sample scattering angle for the single scattering case.
47 // to be used by single and multiple scattering models. References:
48 // 1) G.Wentzel, Z. Phys. 40 (1927) 590.
49 // 2) J.M. Fernandez-Varea et al., NIM B73 (1993) 447.
50 //
51 // -------------------------------------------------------------------
52 //
53 
54 #ifndef G4WentzelOKandVIxSection_h
55 #define G4WentzelOKandVIxSection_h 1
56 
57 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
58 
59 #include "globals.hh"
60 #include "G4Material.hh"
61 #include "G4Element.hh"
62 #include "G4ElementVector.hh"
63 #include "G4NistManager.hh"
64 #include "G4ThreeVector.hh"
65 #include "G4Pow.hh"
66 
68 
69 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
70 
72 {
73 
74 public:
75 
76  G4WentzelOKandVIxSection(G4bool combined = true);
77 
78  virtual ~G4WentzelOKandVIxSection();
79 
80  void Initialise(const G4ParticleDefinition*, G4double CosThetaLim);
81 
83 
84  // return cos(ThetaMax) for msc and cos(thetaMin) for single scattering
85  // cut = DBL_MAX means no scattering off electrons
87 
89 
91  G4double CosThetaMax,
92  G4double elecRatio = 0.0);
93 
95 
96  inline G4double ComputeNuclearCrossSection(G4double CosThetaMin,
97  G4double CosThetaMax);
98 
100  G4double CosThetaMax);
101 
102  inline G4double SetupKinematic(G4double kinEnergy, const G4Material* mat);
103 
104  inline void SetTargetMass(G4double value);
105 
106  inline G4double GetMomentumSquare() const;
107 
108  inline G4double GetCosThetaNuc() const;
109 
110  inline G4double GetCosThetaElec() const;
111 
112 private:
113 
115 
116  // hide assignment operator
119 
124 
127 
129 
131 
132  // integer parameters
135 
137 
138  // single scattering parameters
144 
145  // projectile
147 
160 
161  // target
172 
174  static G4double ScreenRSquare[100];
175  static G4double FormFactor[100];
176 };
177 
178 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
179 
180 inline G4double
182 {
183  if(ekin != tkin || mat != currentMaterial) {
184  currentMaterial = mat;
185  tkin = ekin;
186  mom2 = tkin*(tkin + 2.0*mass);
187  invbeta2 = 1.0 + mass*mass/mom2;
188  factB = spin/invbeta2;
190  if(isCombined) {
192  1.-factorA2*mat->GetIonisation()->GetInvA23()/mom2);
193  }
194  }
195  return cosTetMaxNuc;
196 }
197 
198 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
199 
201 {
202  targetMass = value;
203  factD = std::sqrt(mom2)/value;
204 }
205 
206 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
207 
209 {
210  return mom2;
211 }
212 
213 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
214 
216 {
217  return cosTetMaxNuc;
218 }
219 
220 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
221 
223 {
224  return cosTetMaxElec;
225 }
226 
227 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
228 
229 inline G4double
231  G4double cosTMax)
232 {
233  G4double xsec = 0.0;
234  if(cosTMax < cosTMin) {
235  xsec = targetZ*kinFactor*(cosTMin - cosTMax)/
236  ((1.0 - cosTMin + screenZ)*(1.0 - cosTMax + screenZ));
237  }
238  return xsec;
239 }
240 
241 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
242 
243 inline G4double
245  G4double cosTMax)
246 {
247  G4double xsec = 0.0;
248  G4double cost1 = std::max(cosTMin,cosTetMaxElec);
249  G4double cost2 = std::max(cosTMax,cosTetMaxElec);
250  if(cost1 > cost2) {
251  xsec = kinFactor*(cost1 - cost2)/
252  ((1.0 - cost1 + screenZ)*(1.0 - cost2 + screenZ));
253  }
254  return xsec;
255 }
256 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
257 
258 #endif
259 
const G4ParticleDefinition * theElectron
const G4ParticleDefinition * theProton
void SetupParticle(const G4ParticleDefinition *)
G4IonisParamMat * GetIonisation() const
Definition: G4Material.hh:226
static G4double ScreenRSquareElec[100]
G4WentzelOKandVIxSection(G4bool combined=true)
CLHEP::Hep3Vector G4ThreeVector
Definition: G4Pow.hh:56
G4double ComputeSecondTransportMoment(G4double CosThetaMax)
const G4ParticleDefinition * particle
void ComputeMaxElectronScattering(G4double cut)
int G4int
Definition: G4Types.hh:78
G4double SetupTarget(G4int Z, G4double cut=DBL_MAX)
G4double ComputeTransportCrossSectionPerAtom(G4double CosThetaMax)
bool G4bool
Definition: G4Types.hh:79
G4double GetInvA23() const
G4double SetupKinematic(G4double kinEnergy, const G4Material *mat)
G4double ComputeElectronCrossSection(G4double CosThetaMin, G4double CosThetaMax)
G4WentzelOKandVIxSection & operator=(const G4WentzelOKandVIxSection &right)
const G4ParticleDefinition * thePositron
T max(const T t1, const T t2)
brief Return the largest of the two arguments
G4double ComputeNuclearCrossSection(G4double CosThetaMin, G4double CosThetaMax)
double G4double
Definition: G4Types.hh:76
#define DBL_MAX
Definition: templates.hh:83
void Initialise(const G4ParticleDefinition *, G4double CosThetaLim)
G4ThreeVector & SampleSingleScattering(G4double CosThetaMin, G4double CosThetaMax, G4double elecRatio=0.0)