Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4PolarizedComptonCrossSection.cc
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// GEANT4 Class file
29//
30//
31// File name: G4PolarizedComptonCrossSection
32//
33// Author: Andreas Schaelicke
34//
35// Creation date: 15.05.2005
36//
37// Modifications:
38//
39// Class Description:
40// determine the polarization of the final state
41// in a Compton scattering process employing the differential
42// cross section by F.W.Lipps & H.A.Tolhoek
43// ( Physica 20 (1954) 395 )
44// recalculated by P.Starovoitov
45//
48#include "Randomize.hh"
49
50//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
52 gammaPol0(false), electronPol1(false), gammaPol2(false), electronPol3(false),
53 epsilon(0.)
54 {
55 SetYmin(0.);
56
57 // G4cout<<"G4PolarizedComptonCrossSection() init\n";
58
59 re2 = classic_electr_radius * classic_electr_radius * sqr(4*pi/hbarc);
60 // G4double unit_conversion = hbarc_squared ;
61 // G4cout<<" (keV)^2* m^2 ="<<unit_conversion<<"\n";
62 phi0 = 0.; polXS = 0.; unpXS = 0.;
63 phi2 = G4ThreeVector(0., 0., 0.);
64 phi3 = G4ThreeVector(0., 0., 0.);
65 polxx = polyy = polzz = polxz = polzx = polyz = polzy = polxy = polyx = 0.;
66 diffXSFactor = 1.;
67 totalXSFactor = 1.;
68}
69
70//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
72{}
73
74//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
76 const G4StokesVector & pol0,
77 const G4StokesVector & pol1,
78 G4int flag)
79{
80 G4double cosT = 1. - (1./eps - 1.)/X;
81 if(cosT > 1.+1.e-8) cosT = 1.;
82 if(cosT < -1.-1.e-8) cosT = -1.;
83 G4double cosT2 = cosT*cosT;
84 G4double cosT3 = cosT2*cosT;
85 G4double sinT2 = 1. - cosT2;
86 if(sinT2 > 1. + 1.e-8) sinT2 = 1.;
87 if(sinT2 < 0.) sinT2 = 0.;
88 G4double sinT = std::sqrt(sinT2);
89 G4double cos2T = 2.*cosT2 - 1.;
90 G4double sin2T = 2.*sinT*cosT;
91 G4double eps2 = sqr(eps);
92 DefineCoefficients(pol0,pol1);
93 diffXSFactor = re2/(4.*X);
94
95 // unpolarized Cross Section
96 unpXS = (eps2 + 1. - eps*sinT2)/(2.*eps);
97 // initial polarization dependence
98 polXS = -sinT2*pol0.x() + (1. - eps)*sinT*polzx + ((eps2 - 1.)/eps)*cosT*polzz;
99 polXS *= 0.5;
100
101 phi0 = unpXS + polXS;
102
103 if (flag == 2 ){
104 // polarization of outgoing photon
105 G4double PHI21 = -sinT2 + 0.5*(cos2T + 3.)*pol0.x() - ((1. - eps)/eps)*sinT*polzx;
106 PHI21 *= 0.5;
107 G4double PHI22 = cosT*pol0.y() + ((1. - eps)/(2.*eps))*sinT*polzy;
108 G4double PHI23 = ((eps2 + 1.)/eps)*cosT*pol0.z() - ((1. - eps)/eps)*(eps*cosT2 + 1.)*pol1.z();
109 PHI23 += 0.5*(1. - eps)*sin2T*pol1.x();
110 PHI23 += (eps - 1.)*(-sinT2*polxz + sinT*polyy - 0.5*sin2T*polxx);
111 PHI23 *= 0.5;
112 phi2 = G4ThreeVector(PHI21, PHI22, PHI23);
113
114 // polarization of outgoing electron
115 G4double PHI32 = -sinT2*polxy + ((1. - eps)/eps)*sinT*polyz + 0.5*(cos2T + 3.)*pol1.y();
116 PHI32 *= 0.5;
117
118 G4double PHI31 = 0., PHI31add = 0., PHI33 = 0., PHI33add = 0.;
119
120 if ((1. - eps) > 1.e-12){
121 G4double helpVar = std::sqrt(eps2 - 2.*cosT*eps + 1.);
122
123 PHI31 = (1. - eps)*(1. + cosT)*sinT*pol0.z();
124 PHI31 += (-eps*cosT3 + eps*cosT2 + (eps - 2.)*cosT + eps)*pol1.x();
125 PHI31 += -(eps*cosT2 - eps*cosT + cosT + 1.)*sinT*pol1.z();
126 PHI31 /= 2.*helpVar;
127
128 PHI31add = -eps*sqr(1. - cosT)*(1. + cosT)*polxx;
129 PHI31add += (1. - eps)*sinT2*polyy;
130 PHI31add += -(-eps2 + cosT*(cosT*eps - eps + 1.)*eps + eps - 1.)*sinT*polxz/eps;
131 PHI31add /= 2.*helpVar;
132
133 PHI33 = ((1. - eps)/eps)*(-eps*cosT2 + eps*(eps + 1.)*cosT - 1.)*pol0.z();
134 PHI33 += -(eps*cosT2 + (1. - eps)*eps*cosT + 1.)*sinT*pol1.x();
135 PHI33 += -(-eps2*cosT3 + eps*(eps2 - eps + 1.)*cosT2 - cosT + eps2)*pol1.z()/eps;
136 PHI33 /= -2.*helpVar;
137
138 PHI33add = (eps*(eps - cosT - 1.)*cosT + 1.)*sinT*polxx;
139 PHI33add += -(-eps2 + cosT*eps + eps - 1.)*sinT2*polxz;
140 PHI33add += (eps - 1.)*(cosT - eps)*sinT*polyy;
141 PHI33add /= -2.*helpVar;
142 }else{
143 PHI31 = -pol1.z() - (X - 1.)*std::sqrt(1. - eps)*pol1.x()/std::sqrt(2.*X);
144 PHI31add = -(-X*X*pol1.z() - 2.*X*(2.*pol0.z() - pol1.z()) - (4.*pol0.x() + 5.)*pol1.z())*(1. - eps)/(4.*X);
145
146 PHI33 = pol1.x() - (X - 1.)*std::sqrt(1. - eps)*pol1.z()/std::sqrt(2.*X);
147 PHI33add = -(X*X - 2.*X + 4.*pol0.x() + 5.)*(1. - eps)*pol1.x()/(4.*X);
148 }
149 phi3 = G4ThreeVector(PHI31 + PHI31add, PHI32, PHI33 + PHI33add);
150
151 }
152 unpXS *= diffXSFactor;
153 polXS *= diffXSFactor;
154 phi0 *= diffXSFactor;
155 phi2 *= diffXSFactor;
156 phi3 *= diffXSFactor;
157
158}
159
161{
162 gammaPol2 = !(pol2==G4StokesVector::ZERO);
163 electronPol3 = !(pol3==G4StokesVector::ZERO);
164
165 G4double phi = 0.;
166 // polarization independent part
167 phi += phi0;
168
169
170 if (gammaPol2) {
171 // part depending on the polarization of the final photon
172 phi += phi2*pol2;
173 }
174
175 if (electronPol3) {
176 // part depending on the polarization of the final electron
177 phi += phi3*pol3;
178 }
179
180 // return cross section.
181 return phi;
182}
183
184//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
185
186
188 const G4StokesVector & pol0,
189 const G4StokesVector & pol1)
190{
191
192 // G4double k0 = gammaEnergy / electron_mass_c2 ;
193 G4double k1 = 1 + 2*k0 ;
194
195// // pi*re^2
196// G4double re=2.81794e-15; //m
197// G4double barn=1.e-28; //m^2
198 G4double Z=theZ;
199
200 G4double unit = Z*pi*classic_electr_radius
201 * classic_electr_radius ; // *1./barn;
202
203 G4double pre = unit/(sqr(k0)*sqr(1.+2.*k0));
204
205 G4double xs_0 = ((k0 - 2.)*k0 -2.)*sqr(k1)*std::log(k1) + 2.*k0*(k0*(k0 + 1.)*(k0 + 8.) + 2.);
206 G4double xs_pol = (k0 + 1.)*sqr(k1)*std::log(k1) - 2.*k0*(5.*sqr(k0) + 4.*k0 + 1.);
207
208 return pre*(xs_0/k0 + pol0.p3()*pol1.z()*xs_pol);
209}
210
211
212//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
213
214
216{
217 // Note, mean polarization can not contain correlation
218 // effects.
219 return 1./phi0 * phi2;
220}
221
222//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
223
225{
226 // Note, mean polarization can not contain correlation
227 // effects.
228 return 1./phi0 * phi3;
229}
230
231
232
233//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
234void G4PolarizedComptonCrossSection::DefineCoefficients(const G4StokesVector & pol0,
235 const G4StokesVector & pol1)
236{
237 polxx=pol0.x()*pol1.x();
238 polyy=pol0.y()*pol1.y();
239 polzz=pol0.z()*pol1.z();
240
241 polxz=pol0.x()*pol1.z();
242 polzx=pol0.z()*pol1.x();
243
244 polyz=pol0.y()*pol1.z();
245 polzy=pol0.z()*pol1.y();
246
247 polxy=pol0.x()*pol1.y();
248 polyx=pol0.y()*pol1.x();
249}
250
251
252//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
double z() const
double x() const
double y() const
virtual void Initialize(G4double eps, G4double X, G4double phi, const G4StokesVector &p0, const G4StokesVector &p1, G4int flag=0)
G4double XSection(const G4StokesVector &pol2, const G4StokesVector &pol3)
G4double TotalXSection(G4double xmin, G4double xmax, G4double y, const G4StokesVector &pol0, const G4StokesVector &pol1)
G4double p3() const
static const G4StokesVector ZERO
T sqr(const T &x)
Definition: templates.hh:145