Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ScreeningMottCrossSection.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// G4ScreeningMottCrossSection.hh
27//-------------------------------------------------------------------
28//
29// GEANT4 Class header file
30//
31// File name: G4ScreeningMottCrossSection
32//
33// Author: Cristina Consolandi
34//
35// Creation date: 20.10.2011
36//
37// Modifications:
38//
39// Class Description:
40// Computation of electron Coulomb Scattering Cross Section.
41// Suitable for high energy electrons and light target materials.
42//
43// Reference:
44// M.J. Boschini et al.
45// "Non Ionizing Energy Loss induced by Electrons in the Space Environment"
46// Proc. of the 13th Int. Conf. on Particle Physics and Advanced Technology
47// (13th ICPPAT, Como 3-7/10/2011), World Scientific (Singapore).
48// Available at: http://arxiv.org/abs/1111.4042v4
49//
50// 1) Mott Differential Cross Section Approximation:
51// For Target material up to Z=92 (U):
52// As described in http://arxiv.org/abs/1111.4042v4
53// par. 2.1 , eq. (16)-(17)
54// Else (Z>92):
55// W. A. McKinley and H. Fashbach, Phys. Rev. 74, (1948) 1759.
56// 2) Screening coefficient:
57// vomn G. Moliere, Z. Naturforsh A2 (1947), 133-145; A3 (1948), 78.
58// 3) Nuclear Form Factor:
59// A.V. Butkevich et al. Nucl. Instr. Meth. A488 (2002), 282-294.
60//
61// -----------------------------------------------------------------------------
62
63//
64#ifndef G4ScreeningMottCrossSection_h
65#define G4ScreeningMottCrossSection_h 1
66
67//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
68
69#include "globals.hh"
71#include <vector>
72
73//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
74
75static const G4int DIMMOTT = 750;
76
77class G4NistManager;
78class G4Pow;
79
81{
82
83public:
84
86
88
89 void Initialise(const G4ParticleDefinition*, G4double cosThetaLim);
90
91 void SetupKinematic(G4double kinEnergy, G4int Z);
92
95
98
100 inline void SetupParticle(const G4ParticleDefinition*);
101
103 (const G4ScreeningMottCrossSection &right) = delete;
105
106private:
107
108 G4double ComputeAngle(G4int idx, G4double& rand);
109
110 G4double FormFactor2ExpHof(G4double sin2t2);
111 G4double FormFactor2Gauss(G4double sin2t2);
112 G4double FormFactor2UniformHelm(G4double sin2t2);
113 G4double DifferentialXSection(G4int idx, G4int form);
114
115 G4double GetTransitionRandom();
116
117 G4NistManager* fNistManager;
118 G4Pow* fG4pow;
119
120 const G4ParticleDefinition* particle;
121
122 G4double fTotalCross;
123 //cost - min - max
124 G4double cosThetaMin;// def 1.0
125 G4double cosThetaMax;// def -1.0
126
127 G4double cosTetMinNuc;
128 G4double cosTetMaxNuc;
129
130 //energy cut
131 G4double ecut;
132 G4double etag;
133
134 G4double spin;
135 G4double mass;
136
137 //lab of incedent particle
138 G4double tkinLab;
139 G4double momLab2;
140 G4double invbetaLab2;
141
142 //relative system with nucleus
143 G4double mu_rel;
144 G4double tkin;
145 G4double mom2;
146 G4double invbeta2;
147 G4double beta;
148 G4double gamma;
149
150 //constants
151 G4double alpha;
152 G4double htc2;
153 G4double e2;
154
155 // target nucleus
156 G4double targetMass;
157 G4double As;
158 G4int targetZ;
159 G4int targetA;
160
161 // working array
162 std::vector<G4double> cross;
163};
164
165
166//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
167
168inline
170{
171 particle = p;
172 mass = particle->GetPDGMass();
173 spin = particle->GetPDGSpin();
174 if(0.0 != spin) { spin = 0.5; }
175 tkin = 0.0;
176}
177
178//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
179
180#endif
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
Definition: G4Pow.hh:49
G4double NuclearCrossSection(G4int form, G4int fast)
void SetupKinematic(G4double kinEnergy, G4int Z)
void Initialise(const G4ParticleDefinition *, G4double cosThetaLim)
G4ScreeningMottCrossSection(const G4ScreeningMottCrossSection &)=delete
G4double GetScatteringAngle(G4int form, G4int fast)
G4double McFcorrection(G4double sin2t2)
void SetupParticle(const G4ParticleDefinition *)
G4double RatioMottRutherfordCosT(G4double sin2t2)