Geant4-11
Public Member Functions | Protected Member Functions | Protected Attributes | Private Member Functions | Private Attributes | Static Private Attributes
G4PolarizedComptonXS Class Reference

#include <G4PolarizedComptonXS.hh>

Inheritance diagram for G4PolarizedComptonXS:
G4VPolarizedXS

Public Member Functions

 G4PolarizedComptonXS ()
 
 G4PolarizedComptonXS (const G4PolarizedComptonXS &)=delete
 
G4StokesVector GetPol2 () override
 
G4StokesVector GetPol3 () override
 
virtual G4double GetXmax (G4double y)
 
virtual G4double GetXmin (G4double y)
 
G4double GetYmin ()
 
void Initialize (G4double eps, G4double X, G4double phi, const G4StokesVector &p0, const G4StokesVector &p1, G4int flag=0) override
 
G4PolarizedComptonXSoperator= (const G4PolarizedComptonXS &right)=delete
 
void SetMaterial (G4double A, G4double Z, G4double coul)
 
G4double TotalXSection (G4double xmin, G4double xmax, G4double y, const G4StokesVector &pol0, const G4StokesVector &pol1) override
 
G4double XSection (const G4StokesVector &pol2, const G4StokesVector &pol3) override
 
 ~G4PolarizedComptonXS () override
 

Protected Member Functions

void SetXmax (G4double xmax)
 
void SetXmin (G4double xmin)
 
void SetYmin (G4double ymin)
 

Protected Attributes

G4double fA
 
G4double fCoul
 
G4double fXmax
 
G4double fXmin
 
G4double fYmin
 
G4double fZ
 

Private Member Functions

void DefineCoefficients (const G4StokesVector &pol0, const G4StokesVector &pol1)
 

Private Attributes

G4double fPhi0
 
G4ThreeVector fPhi2
 
G4ThreeVector fPhi3
 
G4double fPolXS
 
G4double fUnpXS
 
G4double polxx
 
G4double polxy
 
G4double polxz
 
G4double polyx
 
G4double polyy
 
G4double polyz
 
G4double polzx
 
G4double polzy
 
G4double polzz
 

Static Private Attributes

static constexpr G4double re2
 

Detailed Description

Definition at line 45 of file G4PolarizedComptonXS.hh.

Constructor & Destructor Documentation

◆ G4PolarizedComptonXS() [1/2]

G4PolarizedComptonXS::G4PolarizedComptonXS ( )

Definition at line 41 of file G4PolarizedComptonXS.cc.

42{
43 SetYmin(0.);
44
45 fPhi0 = 0.;
46 fPolXS = 0.;
47 fUnpXS = 0.;
48 fPhi2 = G4ThreeVector(0., 0., 0.);
49 fPhi3 = G4ThreeVector(0., 0., 0.);
50 polxx = polyy = polzz = polxz = polzx = polyz = polzy = polxy = polyx = 0.;
51}
CLHEP::Hep3Vector G4ThreeVector
void SetYmin(G4double ymin)

References fPhi0, fPhi2, fPhi3, fPolXS, fUnpXS, polxx, polxy, polxz, polyx, polyy, polyz, polzx, polzy, polzz, and G4VPolarizedXS::SetYmin().

◆ ~G4PolarizedComptonXS()

G4PolarizedComptonXS::~G4PolarizedComptonXS ( )
override

Definition at line 54 of file G4PolarizedComptonXS.cc.

54{}

◆ G4PolarizedComptonXS() [2/2]

G4PolarizedComptonXS::G4PolarizedComptonXS ( const G4PolarizedComptonXS )
delete

Member Function Documentation

◆ DefineCoefficients()

void G4PolarizedComptonXS::DefineCoefficients ( const G4StokesVector pol0,
const G4StokesVector pol1 
)
private

Definition at line 228 of file G4PolarizedComptonXS.cc.

230{
231 polxx = pol0.x() * pol1.x();
232 polyy = pol0.y() * pol1.y();
233 polzz = pol0.z() * pol1.z();
234
235 polxz = pol0.x() * pol1.z();
236 polzx = pol0.z() * pol1.x();
237
238 polyz = pol0.y() * pol1.z();
239 polzy = pol0.z() * pol1.y();
240
241 polxy = pol0.x() * pol1.y();
242 polyx = pol0.y() * pol1.x();
243}
double z() const
double x() const
double y() const

References polxx, polxy, polxz, polyx, polyy, polyz, polzx, polzy, polzz, CLHEP::Hep3Vector::x(), CLHEP::Hep3Vector::y(), and CLHEP::Hep3Vector::z().

Referenced by Initialize().

◆ GetPol2()

G4StokesVector G4PolarizedComptonXS::GetPol2 ( )
overridevirtual

Reimplemented from G4VPolarizedXS.

Definition at line 214 of file G4PolarizedComptonXS.cc.

215{
216 // Note, mean polarization can not contain correlation effects.
217 return G4StokesVector(1. / fPhi0 * fPhi2);
218}

References fPhi0, and fPhi2.

◆ GetPol3()

G4StokesVector G4PolarizedComptonXS::GetPol3 ( )
overridevirtual

Reimplemented from G4VPolarizedXS.

Definition at line 221 of file G4PolarizedComptonXS.cc.

222{
223 // Note, mean polarization can not contain correlation effects.
224 return G4StokesVector(1. / fPhi0 * fPhi3);
225}

References fPhi0, and fPhi3.

◆ GetXmax()

G4double G4VPolarizedXS::GetXmax ( G4double  y)
virtualinherited

Reimplemented in G4PolarizedAnnihilationXS.

Definition at line 84 of file G4VPolarizedXS.cc.

84{ return fXmax; }

References G4VPolarizedXS::fXmax.

◆ GetXmin()

G4double G4VPolarizedXS::GetXmin ( G4double  y)
virtualinherited

Reimplemented in G4PolarizedAnnihilationXS.

Definition at line 81 of file G4VPolarizedXS.cc.

81{ return fXmin; }

References G4VPolarizedXS::fXmin.

◆ GetYmin()

G4double G4VPolarizedXS::GetYmin ( )
inlineinherited

Definition at line 65 of file G4VPolarizedXS.hh.

65{ return fYmin; }

References G4VPolarizedXS::fYmin.

◆ Initialize()

void G4PolarizedComptonXS::Initialize ( G4double  eps,
G4double  X,
G4double  phi,
const G4StokesVector p0,
const G4StokesVector p1,
G4int  flag = 0 
)
overridevirtual

Implements G4VPolarizedXS.

Definition at line 57 of file G4PolarizedComptonXS.cc.

61{
62 G4double cosT = 1. - (1. / eps - 1.) / X;
63 if(cosT > 1. + 1.e-8)
64 cosT = 1.;
65 else if(cosT < -1. - 1.e-8)
66 cosT = -1.;
67 G4double cosT2 = cosT * cosT;
68 G4double cosT3 = cosT2 * cosT;
69 G4double sinT2 = 1. - cosT2;
70 if(sinT2 > 1. + 1.e-8)
71 sinT2 = 1.;
72 else if(sinT2 < 0.)
73 sinT2 = 0.;
74 G4double sinT = std::sqrt(sinT2);
75 G4double cos2T = 2. * cosT2 - 1.;
76 G4double sin2T = 2. * sinT * cosT;
77 G4double eps2 = sqr(eps);
78 DefineCoefficients(pol0, pol1);
79 G4double diffXSFactor = re2 / (4. * X);
80
81 // unpolarized Cross Section
82 fUnpXS = (eps2 + 1. - eps * sinT2) / (2. * eps);
83 // initial polarization dependence
84 fPolXS = -sinT2 * pol0.x() + (1. - eps) * sinT * polzx +
85 ((eps2 - 1.) / eps) * cosT * polzz;
86 fPolXS *= 0.5;
87
89
90 if(flag == 2)
91 {
92 // polarization of outgoing photon
93 G4double phi21 = -sinT2 + 0.5 * (cos2T + 3.) * pol0.x() -
94 ((1. - eps) / eps) * sinT * polzx;
95 phi21 *= 0.5;
96 G4double phi22 = cosT * pol0.y() + ((1. - eps) / (2. * eps)) * sinT * polzy;
97 G4double phi23 = ((eps2 + 1.) / eps) * cosT * pol0.z() -
98 ((1. - eps) / eps) * (eps * cosT2 + 1.) * pol1.z();
99 phi23 += 0.5 * (1. - eps) * sin2T * pol1.x();
100 phi23 += (eps - 1.) * (-sinT2 * polxz + sinT * polyy - 0.5 * sin2T * polxx);
101 phi23 *= 0.5;
102
103 fPhi2 = G4ThreeVector(phi21, phi22, phi23);
104
105 // polarization of outgoing electron
106 G4double phi32 = -sinT2 * polxy + ((1. - eps) / eps) * sinT * polyz +
107 0.5 * (cos2T + 3.) * pol1.y();
108 phi32 *= 0.5;
109
110 G4double phi31 = 0.;
111 G4double phi31add = 0.;
112 G4double phi33 = 0.;
113 G4double phi33add = 0.;
114
115 if((1. - eps) > 1.e-12)
116 {
117 G4double helpVar = std::sqrt(eps2 - 2. * cosT * eps + 1.);
118
119 phi31 = (1. - eps) * (1. + cosT) * sinT * pol0.z();
120 phi31 +=
121 (-eps * cosT3 + eps * cosT2 + (eps - 2.) * cosT + eps) * pol1.x();
122 phi31 += -(eps * cosT2 - eps * cosT + cosT + 1.) * sinT * pol1.z();
123 phi31 /= 2. * helpVar;
124
125 phi31add = -eps * sqr(1. - cosT) * (1. + cosT) * polxx;
126 phi31add += (1. - eps) * sinT2 * polyy;
127 phi31add += -(-eps2 + cosT * (cosT * eps - eps + 1.) * eps + eps - 1.) *
128 sinT * polxz / eps;
129 phi31add /= 2. * helpVar;
130
131 phi33 = ((1. - eps) / eps) *
132 (-eps * cosT2 + eps * (eps + 1.) * cosT - 1.) * pol0.z();
133 phi33 += -(eps * cosT2 + (1. - eps) * eps * cosT + 1.) * sinT * pol1.x();
134 phi33 +=
135 -(-eps2 * cosT3 + eps * (eps2 - eps + 1.) * cosT2 - cosT + eps2) *
136 pol1.z() / eps;
137 phi33 /= -2. * helpVar;
138
139 phi33add = (eps * (eps - cosT - 1.) * cosT + 1.) * sinT * polxx;
140 phi33add += -(-eps2 + cosT * eps + eps - 1.) * sinT2 * polxz;
141 phi33add += (eps - 1.) * (cosT - eps) * sinT * polyy;
142 phi33add /= -2. * helpVar;
143 }
144 else
145 {
146 phi31 = -pol1.z() -
147 (X - 1.) * std::sqrt(1. - eps) * pol1.x() / std::sqrt(2. * X);
148 phi31add = -(-X * X * pol1.z() - 2. * X * (2. * pol0.z() - pol1.z()) -
149 (4. * pol0.x() + 5.) * pol1.z()) *
150 (1. - eps) / (4. * X);
151
152 phi33 = pol1.x() -
153 (X - 1.) * std::sqrt(1. - eps) * pol1.z() / std::sqrt(2. * X);
154 phi33add = -(X * X - 2. * X + 4. * pol0.x() + 5.) * (1. - eps) *
155 pol1.x() / (4. * X);
156 }
157 fPhi3 = G4ThreeVector(phi31 + phi31add, phi32, phi33 + phi33add);
158 }
159 fUnpXS *= diffXSFactor;
160 fPolXS *= diffXSFactor;
161 fPhi0 *= diffXSFactor;
162 fPhi2 *= diffXSFactor;
163 fPhi3 *= diffXSFactor;
164}
static const G4double eps
double G4double
Definition: G4Types.hh:83
static constexpr G4double re2
void DefineCoefficients(const G4StokesVector &pol0, const G4StokesVector &pol1)
T sqr(const T &x)
Definition: templates.hh:128

References DefineCoefficients(), eps, fPhi0, fPhi2, fPhi3, fPolXS, fUnpXS, polxx, polxy, polxz, polyy, polyz, polzx, polzy, polzz, re2, sqr(), CLHEP::Hep3Vector::x(), CLHEP::Hep3Vector::y(), and CLHEP::Hep3Vector::z().

◆ operator=()

G4PolarizedComptonXS & G4PolarizedComptonXS::operator= ( const G4PolarizedComptonXS right)
delete

◆ SetMaterial()

void G4VPolarizedXS::SetMaterial ( G4double  A,
G4double  Z,
G4double  coul 
)
inlineinherited

◆ SetXmax()

void G4VPolarizedXS::SetXmax ( G4double  xmax)
inlineprotectedinherited

◆ SetXmin()

void G4VPolarizedXS::SetXmin ( G4double  xmin)
inlineprotectedinherited

Definition at line 83 of file G4VPolarizedXS.hh.

83{ fXmin = xmin; }

References G4VPolarizedXS::fXmin.

◆ SetYmin()

void G4VPolarizedXS::SetYmin ( G4double  ymin)
inlineprotectedinherited

Definition at line 85 of file G4VPolarizedXS.hh.

85{ fYmin = ymin; }

References G4VPolarizedXS::fYmin.

Referenced by G4PolarizedComptonXS().

◆ TotalXSection()

G4double G4PolarizedComptonXS::TotalXSection ( G4double  xmin,
G4double  xmax,
G4double  y,
const G4StokesVector pol0,
const G4StokesVector pol1 
)
overridevirtual

Reimplemented from G4VPolarizedXS.

Definition at line 193 of file G4PolarizedComptonXS.cc.

197{
198 G4double k1 = 1. + 2. * k0;
199
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) +
206 2. * k0 * (k0 * (k0 + 1.) * (k0 + 8.) + 2.);
207 G4double xs_pol = (k0 + 1.) * sqr(k1) * std::log(k1) -
208 2. * k0 * (5. * sqr(k0) + 4. * k0 + 1.);
209
210 return pre * (xs_0 / k0 + pol0.p3() * pol1.z() * xs_pol);
211}
G4double p3() const
static constexpr double classic_electr_radius
static constexpr double pi
Definition: SystemOfUnits.h:55

References CLHEP::classic_electr_radius, G4VPolarizedXS::fZ, G4InuclParticleNames::k0, G4StokesVector::p3(), CLHEP::pi, sqr(), and CLHEP::Hep3Vector::z().

◆ XSection()

G4double G4PolarizedComptonXS::XSection ( const G4StokesVector pol2,
const G4StokesVector pol3 
)
overridevirtual

Implements G4VPolarizedXS.

Definition at line 166 of file G4PolarizedComptonXS.cc.

168{
169 G4bool gammaPol2 = !(pol2 == G4StokesVector::ZERO);
170 G4bool electronPol3 = !(pol3 == G4StokesVector::ZERO);
171
172 G4double phi = 0.;
173 // polarization independent part
174 phi += fPhi0;
175
176 if(gammaPol2)
177 {
178 // part depending on the polarization of the final photon
179 phi += fPhi2 * pol2;
180 }
181
182 if(electronPol3)
183 {
184 // part depending on the polarization of the final electron
185 phi += fPhi3 * pol3;
186 }
187
188 // return cross section.
189 return phi;
190}
bool G4bool
Definition: G4Types.hh:86
static const G4StokesVector ZERO

References fPhi0, fPhi2, fPhi3, and G4StokesVector::ZERO.

Field Documentation

◆ fA

G4double G4VPolarizedXS::fA
protectedinherited

Definition at line 90 of file G4VPolarizedXS.hh.

Referenced by G4VPolarizedXS::SetMaterial().

◆ fCoul

G4double G4VPolarizedXS::fCoul
protectedinherited

◆ fPhi0

G4double G4PolarizedComptonXS::fPhi0
private

Definition at line 95 of file G4PolarizedComptonXS.hh.

Referenced by G4PolarizedComptonXS(), GetPol2(), GetPol3(), Initialize(), and XSection().

◆ fPhi2

G4ThreeVector G4PolarizedComptonXS::fPhi2
private

Definition at line 91 of file G4PolarizedComptonXS.hh.

Referenced by G4PolarizedComptonXS(), GetPol2(), Initialize(), and XSection().

◆ fPhi3

G4ThreeVector G4PolarizedComptonXS::fPhi3
private

Definition at line 93 of file G4PolarizedComptonXS.hh.

Referenced by G4PolarizedComptonXS(), GetPol3(), Initialize(), and XSection().

◆ fPolXS

G4double G4PolarizedComptonXS::fPolXS
private

Definition at line 100 of file G4PolarizedComptonXS.hh.

Referenced by G4PolarizedComptonXS(), and Initialize().

◆ fUnpXS

G4double G4PolarizedComptonXS::fUnpXS
private

Definition at line 100 of file G4PolarizedComptonXS.hh.

Referenced by G4PolarizedComptonXS(), and Initialize().

◆ fXmax

G4double G4VPolarizedXS::fXmax
protectedinherited

Definition at line 88 of file G4VPolarizedXS.hh.

Referenced by G4VPolarizedXS::GetXmax(), and G4VPolarizedXS::SetXmax().

◆ fXmin

G4double G4VPolarizedXS::fXmin
protectedinherited

Definition at line 88 of file G4VPolarizedXS.hh.

Referenced by G4VPolarizedXS::GetXmin(), and G4VPolarizedXS::SetXmin().

◆ fYmin

G4double G4VPolarizedXS::fYmin
protectedinherited

Definition at line 88 of file G4VPolarizedXS.hh.

Referenced by G4VPolarizedXS::GetYmin(), and G4VPolarizedXS::SetYmin().

◆ fZ

G4double G4VPolarizedXS::fZ
protectedinherited

◆ polxx

G4double G4PolarizedComptonXS::polxx
private

Definition at line 97 of file G4PolarizedComptonXS.hh.

Referenced by DefineCoefficients(), G4PolarizedComptonXS(), and Initialize().

◆ polxy

G4double G4PolarizedComptonXS::polxy
private

Definition at line 97 of file G4PolarizedComptonXS.hh.

Referenced by DefineCoefficients(), G4PolarizedComptonXS(), and Initialize().

◆ polxz

G4double G4PolarizedComptonXS::polxz
private

Definition at line 97 of file G4PolarizedComptonXS.hh.

Referenced by DefineCoefficients(), G4PolarizedComptonXS(), and Initialize().

◆ polyx

G4double G4PolarizedComptonXS::polyx
private

Definition at line 97 of file G4PolarizedComptonXS.hh.

Referenced by DefineCoefficients(), and G4PolarizedComptonXS().

◆ polyy

G4double G4PolarizedComptonXS::polyy
private

Definition at line 97 of file G4PolarizedComptonXS.hh.

Referenced by DefineCoefficients(), G4PolarizedComptonXS(), and Initialize().

◆ polyz

G4double G4PolarizedComptonXS::polyz
private

Definition at line 97 of file G4PolarizedComptonXS.hh.

Referenced by DefineCoefficients(), G4PolarizedComptonXS(), and Initialize().

◆ polzx

G4double G4PolarizedComptonXS::polzx
private

Definition at line 97 of file G4PolarizedComptonXS.hh.

Referenced by DefineCoefficients(), G4PolarizedComptonXS(), and Initialize().

◆ polzy

G4double G4PolarizedComptonXS::polzy
private

Definition at line 97 of file G4PolarizedComptonXS.hh.

Referenced by DefineCoefficients(), G4PolarizedComptonXS(), and Initialize().

◆ polzz

G4double G4PolarizedComptonXS::polzz
private

Definition at line 97 of file G4PolarizedComptonXS.hh.

Referenced by DefineCoefficients(), G4PolarizedComptonXS(), and Initialize().

◆ re2

constexpr G4double G4PolarizedComptonXS::re2
staticconstexprprivate
Initial value:

Definition at line 83 of file G4PolarizedComptonXS.hh.

Referenced by Initialize().


The documentation for this class was generated from the following files: