Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
G4WentzelVIRelXSection.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$
27//
28// -------------------------------------------------------------------
29//
30//
31// GEANT4 Class header file
32//
33//
34// File name: G4WentzelVIRelXSection
35//
36// Authors: V.Ivanchenko
37//
38// Creation date: 08.06.2012 from G4WentzelOKandVIxSection
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 G4WentzelVIRelXSection_h
55#define G4WentzelVIRelXSection_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
74public:
75
77
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 G4double CosThetaMax);
96
98 G4double CosThetaMax);
99
100 inline G4double SetupKinematic(G4double kinEnergy, const G4Material* mat);
101
102 inline void SetTargetMass(G4double value);
103
104 //obsolete method
105 inline void SetRelativisticMass(G4double value);
106
107 inline G4double GetMomentumSquare() const;
108
109 inline G4double GetCosThetaNuc() const;
110
111 inline G4double GetCosThetaElec() const;
112
113private:
114
115 void ComputeMaxElectronScattering(G4double cut);
116
117 // hide assignment operator
118 G4WentzelVIRelXSection & operator=(const G4WentzelVIRelXSection &right);
120
121 const G4ParticleDefinition* theProton;
122 const G4ParticleDefinition* theElectron;
123 const G4ParticleDefinition* thePositron;
124 const G4Material* currentMaterial;
125
126 G4NistManager* fNistManager;
127 G4Pow* fG4pow;
128
129 G4double numlimit;
130
131 G4double elecXSRatio;
132
133 // integer parameters
134 G4int nwarnings;
135 G4int nwarnlimit;
136
137 // single scattering parameters
138 G4double coeff;
139 G4double cosTetMaxElec;
140 G4double cosTetMaxNuc;
141 G4double cosThetaMax;
142 G4double alpha2;
143
144 // projectile
145 const G4ParticleDefinition* particle;
146
147 G4double chargeSquare;
148 G4double charge3;
149 G4double spin;
150 G4double mass;
151 G4double tkin;
152 G4double mom2;
153 G4double momCM2;
154 G4double invbeta2;
155 G4double kinFactor;
156 G4double etag;
157 G4double ecut;
158 G4double lowEnergyLimit;
159
160 // target
161 G4int targetZ;
162 G4double targetMass;
163 G4double screenZ;
164 G4double formfactA;
165 G4double factorA2;
166 G4double factB;
167 G4double factB1;
168 G4double factD;
169 G4double gam0pcmp;
170 G4double pcmp2;
171
172 static G4double ScreenRSquare[100];
173 static G4double FormFactor[100];
174};
175
176//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
177
178inline G4double
180{
181 if(ekin != tkin || mat != currentMaterial) {
182 currentMaterial = mat;
183 tkin = ekin;
184 mom2 = tkin*(tkin + 2.0*mass);
185 invbeta2 = 1.0 + mass*mass/mom2;
186 factB = spin/invbeta2;
187 cosTetMaxNuc =
188 std::max(cosThetaMax,1.-factorA2*mat->GetIonisation()->GetInvA23()/mom2);
189 }
190 return cosTetMaxNuc;
191}
192
193//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
194
196{
197 targetMass = value;
198 factD = std::sqrt(mom2)/value;
199}
200
201//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
202
204{
205 mass = value;
206}
207
208//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
209
211{
212 return mom2;
213}
214
215//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
216
218{
219 return cosTetMaxNuc;
220}
221
222//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
223
225{
226 return cosTetMaxElec;
227}
228
229//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
230
231inline G4double
233 G4double cosTMax)
234{
235 G4double xsec = 0.0;
236 if(cosTMax < cosTMin) {
237 xsec = targetZ*kinFactor*(cosTMin - cosTMax)/
238 ((1.0 - cosTMin + screenZ)*(1.0 - cosTMax + screenZ));
239 }
240 return xsec;
241}
242
243//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
244
245inline G4double
247 G4double cosTMax)
248{
249 G4double xsec = 0.0;
250 G4double cost1 = std::max(cosTMin,cosTetMaxElec);
251 G4double cost2 = std::max(cosTMax,cosTetMaxElec);
252 if(cost1 > cost2) {
253 xsec = kinFactor*(cost1 - cost2)/((1.0 - cost1 + screenZ)*(1.0 - cost2 + screenZ));
254 }
255 return xsec;
256}
257//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
258
259#endif
260
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
G4double GetInvA23() const
G4IonisParamMat * GetIonisation() const
Definition: G4Material.hh:225
Definition: G4Pow.hh:54
void SetTargetMass(G4double value)
G4double ComputeTransportCrossSectionPerAtom(G4double CosThetaMax)
G4double SetupTarget(G4int Z, G4double cut=DBL_MAX)
void SetRelativisticMass(G4double value)
void SetupParticle(const G4ParticleDefinition *)
G4ThreeVector SampleSingleScattering(G4double CosThetaMin, G4double CosThetaMax, G4double elecRatio=0.0)
void Initialise(const G4ParticleDefinition *, G4double CosThetaLim)
G4double SetupKinematic(G4double kinEnergy, const G4Material *mat)
G4double ComputeNuclearCrossSection(G4double CosThetaMin, G4double CosThetaMax)
G4double ComputeElectronCrossSection(G4double CosThetaMin, G4double CosThetaMax)
#define DBL_MAX
Definition: templates.hh:83