Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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 98737 2016-08-09 12:51:38Z 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"
65 #include "G4ThreeVector.hh"
66 #include "G4Pow.hh"
67 
69 
70 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
71 
73 {
74 
75 public:
76 
77  explicit G4WentzelOKandVIxSection(G4bool combined = true);
78 
79  virtual ~G4WentzelOKandVIxSection();
80 
81  void Initialise(const G4ParticleDefinition*, G4double CosThetaLim);
82 
83  void SetupParticle(const G4ParticleDefinition*) ;
84 
85  // return cos(ThetaMax) for msc and cos(thetaMin) for single scattering
86  // cut = DBL_MAX means no scattering off electrons
88 
90 
92  G4double CosThetaMax,
93  G4double elecRatio = 0.0);
94 
96 
97  inline G4double ComputeNuclearCrossSection(G4double CosThetaMin,
98  G4double CosThetaMax);
99 
100  inline G4double ComputeElectronCrossSection(G4double CosThetaMin,
101  G4double CosThetaMax);
102 
103  inline G4double SetupKinematic(G4double kinEnergy, const G4Material* mat);
104 
105  inline void SetTargetMass(G4double value);
106 
107  inline G4double GetMomentumSquare() const;
108 
109  inline G4double GetCosThetaNuc() const;
110 
111  inline G4double GetCosThetaElec() const;
112 
113 private:
114 
115  void ComputeMaxElectronScattering(G4double cut);
116 
117  inline G4double FlatFormfactor(G4double x);
118 
119  // hide assignment operator
120  G4WentzelOKandVIxSection & operator=
121  (const G4WentzelOKandVIxSection &right) = delete;
123 
124  const G4ParticleDefinition* theProton;
125  const G4ParticleDefinition* theElectron;
126  const G4ParticleDefinition* thePositron;
127  const G4Material* currentMaterial;
128 
129  G4NistManager* fNistManager;
130  G4Pow* fG4pow;
131 
132  G4ThreeVector temp;
133 
134  G4double numlimit;
135 
136  // integer parameters
137  G4int nwarnings;
138  G4int nwarnlimit;
139 
140  G4NuclearFormfactorType fNucFormfactor;
141 
142  G4bool isCombined;
143 
144  // single scattering parameters
145  G4double coeff;
146  G4double cosTetMaxElec;
147  G4double cosTetMaxNuc;
148  G4double cosThetaMax;
149  G4double alpha2;
150 
151  // projectile
152  const G4ParticleDefinition* particle;
153 
154  G4double chargeSquare;
155  G4double charge3;
156  G4double spin;
157  G4double mass;
158  G4double tkin;
159  G4double mom2;
160  G4double momCM2;
161  G4double invbeta2;
162  G4double kinFactor;
163  G4double etag;
164  G4double ecut;
165  G4double lowEnergyLimit;
166 
167  // target
168  G4int targetZ;
169  G4double targetMass;
170  G4double screenZ;
171  G4double formfactA;
172  G4double factorA2;
173  G4double factB;
174  G4double factB1;
175  G4double factD;
176  G4double gam0pcmp;
177  G4double pcmp2;
178 
179  static G4double ScreenRSquareElec[100];
180  static G4double ScreenRSquare[100];
181  static G4double FormFactor[100];
182 };
183 
184 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
185 
186 inline G4double
188 {
189  if(ekin != tkin || mat != currentMaterial) {
190  currentMaterial = mat;
191  tkin = ekin;
192  mom2 = tkin*(tkin + 2.0*mass);
193  invbeta2 = 1.0 + mass*mass/mom2;
194  factB = spin/invbeta2;
195  cosTetMaxNuc = cosThetaMax;
196  if(isCombined) {
197  cosTetMaxNuc = std::max(cosTetMaxNuc,
198  1.-factorA2*mat->GetIonisation()->GetInvA23()/mom2);
199  }
200  }
201  return cosTetMaxNuc;
202 }
203 
204 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
205 
207 {
208  targetMass = value;
209  factD = std::sqrt(mom2)/value;
210 }
211 
212 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
213 
215 {
216  return mom2;
217 }
218 
219 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
220 
222 {
223  return cosTetMaxNuc;
224 }
225 
226 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
227 
229 {
230  return cosTetMaxElec;
231 }
232 
233 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
234 
235 inline G4double
237  G4double cosTMax)
238 {
239  G4double xsec = 0.0;
240  if(cosTMax < cosTMin) {
241  xsec = targetZ*kinFactor*(cosTMin - cosTMax)/
242  ((1.0 - cosTMin + screenZ)*(1.0 - cosTMax + screenZ));
243  }
244  return xsec;
245 }
246 
247 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
248 
249 inline G4double
251  G4double cosTMax)
252 {
253  G4double xsec = 0.0;
254  G4double cost1 = std::max(cosTMin,cosTetMaxElec);
255  G4double cost2 = std::max(cosTMax,cosTetMaxElec);
256  if(cost1 > cost2) {
257  xsec = kinFactor*(cost1 - cost2)/
258  ((1.0 - cost1 + screenZ)*(1.0 - cost2 + screenZ));
259  }
260  return xsec;
261 }
262 
263 inline G4double G4WentzelOKandVIxSection::FlatFormfactor(G4double x)
264 {
265  return 3*(std::sin(x) - x*std::cos(x))/(x*x*x);
266 }
267 
268 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
269 
270 #endif
271 
void SetupParticle(const G4ParticleDefinition *)
G4IonisParamMat * GetIonisation() const
Definition: G4Material.hh:226
G4WentzelOKandVIxSection(G4bool combined=true)
Definition: G4Pow.hh:56
G4double ComputeSecondTransportMoment(G4double CosThetaMax)
tuple x
Definition: test.py:50
int G4int
Definition: G4Types.hh:78
G4double SetupTarget(G4int Z, G4double cut=DBL_MAX)
G4double ComputeTransportCrossSectionPerAtom(G4double CosThetaMax)
const XML_Char int const XML_Char * value
Definition: expat.h:331
bool G4bool
Definition: G4Types.hh:79
G4double GetInvA23() const
G4double SetupKinematic(G4double kinEnergy, const G4Material *mat)
G4double ComputeElectronCrossSection(G4double CosThetaMin, G4double CosThetaMax)
T max(const T t1, const T t2)
brief Return the largest of the two arguments
G4NuclearFormfactorType
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)