Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros 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//
27// -------------------------------------------------------------------
28//
29//
30// GEANT4 Class header file
31//
32//
33// File name: G4WentzelOKandVIxSection
34//
35// Authors: V.Ivanchenko and O.Kadri
36//
37// Creation date: 21.05.2010
38//
39// Modifications:
40//
41//
42// Class Description:
43//
44// Implementation of the computation of total and transport cross sections,
45// sample scattering angle for the single scattering case.
46// to be used by single and multiple scattering models. References:
47// 1) G.Wentzel, Z. Phys. 40 (1927) 590.
48// 2) J.M. Fernandez-Varea et al., NIM B73 (1993) 447.
49//
50// -------------------------------------------------------------------
51//
52
53#ifndef G4WentzelOKandVIxSection_h
54#define G4WentzelOKandVIxSection_h 1
55
56//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
57
58#include "globals.hh"
59#include "G4Material.hh"
60#include "G4Element.hh"
61#include "G4ElementVector.hh"
62#include "G4NistManager.hh"
64#include "G4ThreeVector.hh"
65#include "G4Pow.hh"
66
69
70//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
71
73{
74
75public:
76
77 explicit G4WentzelOKandVIxSection(G4bool comb=true);
78
80
81 void Initialise(const G4ParticleDefinition*, G4double CosThetaLim);
82
84
85 // return cos(ThetaMax) for msc and cos(thetaMin) for single scattering
86 // cut = DBL_MAX means no scattering off electrons
87 virtual G4double SetupKinematic(G4double kinEnergy, const G4Material* mat);
89
91
93 G4double CosThetaMax,
94 G4double elecRatio);
95
97
99 G4double CosThetaMax);
100
102 G4double CosThetaMax);
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 // hide assignment operator
114 (const G4WentzelOKandVIxSection &right) = delete;
116
117protected:
118
120
121 void InitialiseA();
122
124
129 const G4Material* currentMaterial = nullptr;
130
133
135
137
138 // single scattering parameters
143
155
156 // target
166
167 // integer parameters
170
172
174
177 static G4double FormFactor[100];
178};
179
180//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
181
183{
184 targetMass = value;
185 factD = std::sqrt(mom2)/value;
186}
187
188//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
189
191{
192 return mom2;
193}
194
195//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
196
198{
199 return cosTetMaxNuc;
200}
201
202//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
203
205{
206 return cosTetMaxElec;
207}
208
209//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
210
211inline G4double
213 G4double cosTMax)
214{
215 return targetZ*kinFactor*fMottFactor*(cosTMin - cosTMax)/
216 ((1.0 - cosTMin + screenZ)*(1.0 - cosTMax + screenZ));
217}
218
219//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
220
221inline G4double
223 G4double cosTMax)
224{
225 G4double cost1 = std::max(cosTMin,cosTetMaxElec);
226 G4double cost2 = std::max(cosTMax,cosTetMaxElec);
227 return (cost1 <= cost2) ? 0.0 : kinFactor*fMottFactor*(cost1 - cost2)/
228 ((1.0 - cost1 + screenZ)*(1.0 - cost2 + screenZ));
229}
230
232{
233 return 3.0*(std::sin(x) - x*std::cos(x))/(x*x*x);
234}
235
236//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
237
238#endif
239
G4NuclearFormfactorType
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
Definition: G4Pow.hh:49
G4double SetupTarget(G4int Z, G4double cut)
static G4double ScreenRSquareElec[100]
G4WentzelOKandVIxSection(const G4WentzelOKandVIxSection &)=delete
G4double ComputeElectronCrossSection(G4double CosThetaMin, G4double CosThetaMax)
G4double ComputeSecondTransportMoment(G4double CosThetaMax)
void ComputeMaxElectronScattering(G4double cut)
const G4ParticleDefinition * theProton
void Initialise(const G4ParticleDefinition *, G4double CosThetaLim)
const G4ParticleDefinition * thePositron
const G4ParticleDefinition * particle
G4double ComputeTransportCrossSectionPerAtom(G4double CosThetaMax)
void SetupParticle(const G4ParticleDefinition *)
virtual G4double SetupKinematic(G4double kinEnergy, const G4Material *mat)
G4ScreeningMottCrossSection * fMottXSection
G4double ComputeNuclearCrossSection(G4double CosThetaMin, G4double CosThetaMax)
G4NuclearFormfactorType fNucFormfactor
G4ThreeVector & SampleSingleScattering(G4double CosThetaMin, G4double CosThetaMax, G4double elecRatio)
const G4ParticleDefinition * theElectron
#define DBL_MAX
Definition: templates.hh:62