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

#include <G4PolarizedAnnihilationXS.hh>

Inheritance diagram for G4PolarizedAnnihilationXS:
G4VPolarizedXS

Public Member Functions

G4double DiceEpsilon ()
 
 G4PolarizedAnnihilationXS ()
 
 G4PolarizedAnnihilationXS (const G4PolarizedAnnihilationXS &)=delete
 
void getCoeff ()
 
virtual G4StokesVector GetPol2 () override
 
virtual G4StokesVector GetPol3 () override
 
G4double getVar (G4int)
 
virtual G4double GetXmax (G4double y) override
 
virtual G4double GetXmin (G4double y) override
 
G4double GetYmin ()
 
virtual void Initialize (G4double eps, G4double gamma, G4double phi, const G4StokesVector &p0, const G4StokesVector &p1, G4int flag=0) override
 
G4PolarizedAnnihilationXSoperator= (const G4PolarizedAnnihilationXS &right)=delete
 
void SetMaterial (G4double A, G4double Z, G4double coul)
 
virtual G4double TotalXSection (G4double xmin, G4double xmax, G4double y, const G4StokesVector &pol0, const G4StokesVector &pol1) override
 
virtual G4double XSection (const G4StokesVector &pol2, const G4StokesVector &pol3) override
 
virtual ~G4PolarizedAnnihilationXS () 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)
 
void TotalXS ()
 

Private Attributes

G4double fDice
 
G4double fPhi0
 
G4ThreeVector fPhi2
 
G4ThreeVector fPhi3
 
G4double fPolXS
 
G4double fUnpXS
 
G4double ISPnd
 
G4double ISPxx
 
G4double ISPyy
 
G4double ISPzz
 
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 47 of file G4PolarizedAnnihilationXS.hh.

Constructor & Destructor Documentation

◆ G4PolarizedAnnihilationXS() [1/2]

G4PolarizedAnnihilationXS::G4PolarizedAnnihilationXS ( )

Definition at line 39 of file G4PolarizedAnnihilationXS.cc.

40 : polxx(0.)
41 , polyy(0.)
42 , polzz(0.)
43 , polxz(0.)
44 , polzx(0.)
45 , polxy(0.)
46 , polyx(0.)
47 , polyz(0.)
48 , polzy(0.)
49 , fPhi0(0.)
50{
51 fPhi2 = G4ThreeVector(0., 0., 0.);
52 fPhi3 = G4ThreeVector(0., 0., 0.);
53 fDice = 0.;
54 fPolXS = 0.;
55 fUnpXS = 0.;
56 ISPxx = ISPyy = ISPzz = ISPnd = 0.;
57}
CLHEP::Hep3Vector G4ThreeVector

References fDice, fPhi2, fPhi3, fPolXS, fUnpXS, ISPnd, ISPxx, ISPyy, and ISPzz.

◆ ~G4PolarizedAnnihilationXS()

G4PolarizedAnnihilationXS::~G4PolarizedAnnihilationXS ( )
overridevirtual

Definition at line 60 of file G4PolarizedAnnihilationXS.cc.

60{}

◆ G4PolarizedAnnihilationXS() [2/2]

G4PolarizedAnnihilationXS::G4PolarizedAnnihilationXS ( const G4PolarizedAnnihilationXS )
delete

Member Function Documentation

◆ DefineCoefficients()

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

Definition at line 264 of file G4PolarizedAnnihilationXS.cc.

266{
267 polxx = pol0.x() * pol1.x();
268 polyy = pol0.y() * pol1.y();
269 polzz = pol0.z() * pol1.z();
270
271 polxz = pol0.x() * pol1.z();
272 polzx = pol0.z() * pol1.x();
273
274 polyz = pol0.y() * pol1.z();
275 polzy = pol0.z() * pol1.y();
276
277 polxy = pol0.x() * pol1.y();
278 polyx = pol0.y() * pol1.x();
279}
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(), and TotalXSection().

◆ DiceEpsilon()

G4double G4PolarizedAnnihilationXS::DiceEpsilon ( )

Definition at line 294 of file G4PolarizedAnnihilationXS.cc.

294{ return fDice; }

References fDice.

Referenced by G4PolarizedAnnihilationModel::SampleSecondaries().

◆ getCoeff()

void G4PolarizedAnnihilationXS::getCoeff ( )

◆ GetPol2()

G4StokesVector G4PolarizedAnnihilationXS::GetPol2 ( )
overridevirtual

Reimplemented from G4VPolarizedXS.

Definition at line 251 of file G4PolarizedAnnihilationXS.cc.

252{
253 // Note, mean polarization can not contain correlation effects.
254 return G4StokesVector(1. / fPhi0 * fPhi2);
255}

References fPhi0, and fPhi2.

Referenced by G4PolarizedAnnihilationModel::SampleSecondaries().

◆ GetPol3()

G4StokesVector G4PolarizedAnnihilationXS::GetPol3 ( )
overridevirtual

Reimplemented from G4VPolarizedXS.

Definition at line 258 of file G4PolarizedAnnihilationXS.cc.

259{
260 // Note, mean polarization can not contain correlation effects.
261 return G4StokesVector(1. / fPhi0 * fPhi3);
262}

References fPhi0, and fPhi3.

Referenced by G4PolarizedAnnihilationModel::SampleSecondaries().

◆ getVar()

G4double G4PolarizedAnnihilationXS::getVar ( G4int  choice)

Definition at line 297 of file G4PolarizedAnnihilationXS.cc.

298{
299 if(choice == -1)
300 return fPolXS / fUnpXS;
301 if(choice == 0)
302 return fUnpXS;
303 if(choice == 1)
304 return ISPxx;
305 if(choice == 2)
306 return ISPyy;
307 if(choice == 3)
308 return ISPzz;
309 if(choice == 4)
310 return ISPnd;
311 return 0;
312}

References fPolXS, fUnpXS, ISPnd, ISPxx, ISPyy, and ISPzz.

Referenced by G4PolarizedAnnihilationModel::SampleSecondaries().

◆ GetXmax()

G4double G4PolarizedAnnihilationXS::GetXmax ( G4double  y)
overridevirtual

Reimplemented from G4VPolarizedXS.

Definition at line 288 of file G4PolarizedAnnihilationXS.cc.

289{
290 return 0.5 * (1. + std::sqrt((y - 1.) / (y + 1.)));
291}

◆ GetXmin()

G4double G4PolarizedAnnihilationXS::GetXmin ( G4double  y)
overridevirtual

Reimplemented from G4VPolarizedXS.

Definition at line 282 of file G4PolarizedAnnihilationXS.cc.

283{
284 return 0.5 * (1. - std::sqrt((y - 1.) / (y + 1.)));
285}

◆ GetYmin()

G4double G4VPolarizedXS::GetYmin ( )
inlineinherited

Definition at line 65 of file G4VPolarizedXS.hh.

65{ return fYmin; }

References G4VPolarizedXS::fYmin.

◆ Initialize()

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

Implements G4VPolarizedXS.

Definition at line 77 of file G4PolarizedAnnihilationXS.cc.

83{
84 G4double diffXSFactor = re2 / (gam - 1.);
85 DefineCoefficients(pol0, pol1);
86
87 // prepare dicing
88 fDice = 0.;
89 G4double symmXS =
90 0.125 * ((-1. / sqr(gam + 1.)) / sqr(eps) +
91 ((sqr(gam) + 4. * gam - 1.) / sqr(gam + 1.)) / eps - 1.);
92
93 G4ThreeVector epsVector(1. / sqr(eps), 1. / eps, 1.);
94 G4ThreeVector oneEpsVector(1. / sqr(1. - eps), 1. / (1. - eps), 1.);
95 G4ThreeVector sumEpsVector(epsVector + oneEpsVector);
96 G4ThreeVector difEpsVector(epsVector - oneEpsVector);
97 G4ThreeVector calcVector(0., 0., 0.);
98
99 // temporary variables
100 G4double helpVar2 = 0., helpVar1 = 0.;
101
102 // unpolarised contribution
103 helpVar1 = (gam * gam + 4. * gam + 1.) / sqr(gam + 1.);
104 helpVar2 = -1. / sqr(gam + 1.);
105 calcVector = G4ThreeVector(helpVar2, helpVar1, -1.);
106 fUnpXS = 0.125 * calcVector * sumEpsVector;
107
108 // initial particles polarised contribution
109 helpVar2 = 1. / sqr(gam + 1.);
110 helpVar1 = -(gam * gam + 4. * gam + 1.) / sqr(gam + 1.);
111 calcVector = G4ThreeVector(helpVar2, helpVar1, 0.5 * (gam + 3.));
112 ISPxx = 0.25 * (calcVector * sumEpsVector) / (gam - 1.);
113
114 helpVar1 = 1. / sqr(gam + 1.);
115 calcVector = G4ThreeVector(-helpVar1, 2. * gam * helpVar1, -1.);
116 ISPyy = 0.125 * calcVector * sumEpsVector;
117
118 helpVar1 = 1. / (gam - 1.);
119 helpVar2 = 1. / sqr(gam + 1.);
120 calcVector = G4ThreeVector(
121 -(gam * gam + 1.) * helpVar2,
122 (gam * gam * (gam + 1.) + 7. * gam + 3.) * helpVar2, -(gam + 3.));
123 ISPzz = 0.125 * helpVar1 * (calcVector * sumEpsVector);
124
125 helpVar1 = std::sqrt(std::fabs(eps * (1. - eps) * 2. * (gam + 1.) - 1.));
126 calcVector = G4ThreeVector(-1. / (gam * gam - 1.), 2. / (gam - 1.), 0.);
127 ISPnd = 0.125 * (calcVector * difEpsVector) * helpVar1;
128
129 fPolXS = ISPxx * polxx + ISPyy * polyy + ISPzz * polzz;
130 fPolXS += ISPnd * (polzx + polxz);
131 fPhi0 = fUnpXS + fPolXS;
132 fDice = symmXS;
133
134 if(polzz != 0.)
135 {
136 fDice *= (1. + (polzz * ISPzz / fUnpXS));
137 if(fDice < 0.)
138 fDice = 0.;
139 }
140 // prepare final state coefficients
141 if(flag == 2)
142 {
143 // circular polarisation
144 G4double circ1 = 0., circ2 = 0., circ3 = 0.;
145 helpVar1 = 8. * sqr(1. - eps) * sqr(eps) * (gam - 1.) * sqr(gam + 1.) /
146 std::sqrt(gam * gam - 1.);
147 helpVar2 = sqr(gam + 1.) * sqr(eps) * (-2. * eps + 3.) -
148 (gam * gam + 3. * gam + 2.) * eps;
149 circ1 = helpVar2 + gam;
150 circ1 /= helpVar1;
151 circ2 = helpVar2 + 1.;
152 circ2 /= helpVar1;
153 helpVar1 = std::sqrt(std::fabs(eps * (1. - eps) * 2. * (gam + 1.) - 1.));
154 helpVar1 /= std::sqrt(gam * gam - 1.);
155 calcVector = G4ThreeVector(1., -2. * gam, 0.);
156 circ3 = 0.125 * (calcVector * sumEpsVector) / (gam + 1.);
157 circ3 *= helpVar1;
158
159 fPhi2.setZ(circ2 * pol1.z() + circ1 * pol0.z() +
160 circ3 * (pol1.x() + pol0.x()));
161 fPhi3.setZ(-circ1 * pol1.z() - circ2 * pol0.z() -
162 circ3 * (pol1.x() + pol0.x()));
163
164 // common to both linear polarisation
165 calcVector = G4ThreeVector(-1., 2. * gam, 0.);
166 G4double linearZero = 0.125 * (calcVector * sumEpsVector) / sqr(gam + 1.);
167
168 // Linear Polarisation #1
169 helpVar1 = std::sqrt(std::fabs(2. * (gam + 1.) * (1. - eps) * eps - 1.)) /
170 ((gam + 1.) * eps * (1. - eps));
171 helpVar2 = helpVar1 * helpVar1;
172
173 // photon 1
174 G4double diagContrib = 0.125 * helpVar2 * (polxx + polyy - polzz);
175 G4double nonDiagContrib =
176 0.125 * helpVar1 * (-polxz / (1. - eps) + polzx / eps);
177
178 fPhi2.setX(linearZero + diagContrib + nonDiagContrib);
179
180 // photon 2
181 nonDiagContrib = 0.125 * helpVar1 * (polxz / eps - polzx / (1. - eps));
182
183 fPhi3.setX(linearZero + diagContrib + nonDiagContrib);
184
185 // Linear Polarisation #2
186 helpVar1 =
187 std::sqrt(gam * gam - 1.) * (2. * (gam + 1.) * eps * (1. - eps) - 1.);
188 helpVar1 /= 8. * sqr(1. - eps) * sqr(eps) * sqr(gam + 1.) * (gam - 1.);
189 helpVar2 = std::sqrt((gam * gam - 1.) *
190 std::fabs(2. * (gam + 1.) * eps * (1. - eps) - 1.));
191 helpVar2 /= 8. * sqr(1. - eps) * sqr(eps) * sqr(gam + 1.) * (gam - 1.);
192
193 G4double contrib21 = (-polxy + polyx) * helpVar1;
194 G4double contrib32 =
195 -(eps * (gam + 1.) - 1.) * polyz + (eps * (gam + 1.) - gam) * polzy;
196
197 contrib32 *= helpVar2;
198 fPhi2.setY(contrib21 + contrib32);
199
200 contrib32 =
201 -(eps * (gam + 1.) - gam) * polyz + (eps * (gam + 1.) - 1.) * polzy;
202 contrib32 *= helpVar2;
203 fPhi3.setY(contrib21 + contrib32);
204 }
205 fPhi0 *= diffXSFactor;
206 fPhi2 *= diffXSFactor;
207 fPhi3 *= diffXSFactor;
208}
static const G4double eps
double G4double
Definition: G4Types.hh:83
void setY(double)
void setZ(double)
void setX(double)
void DefineCoefficients(const G4StokesVector &pol0, const G4StokesVector &pol1)
static constexpr G4double re2
T sqr(const T &x)
Definition: templates.hh:128

References DefineCoefficients(), eps, fDice, fPhi0, fPhi2, fPhi3, fPolXS, fUnpXS, G4InuclParticleNames::gam, ISPnd, ISPxx, ISPyy, ISPzz, polxx, polxy, polxz, polyx, polyy, polyz, polzx, polzy, polzz, re2, CLHEP::Hep3Vector::setX(), CLHEP::Hep3Vector::setY(), CLHEP::Hep3Vector::setZ(), sqr(), CLHEP::Hep3Vector::x(), and CLHEP::Hep3Vector::z().

Referenced by G4PolarizedAnnihilationModel::SampleSecondaries().

◆ operator=()

G4PolarizedAnnihilationXS & G4PolarizedAnnihilationXS::operator= ( const G4PolarizedAnnihilationXS 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::G4PolarizedComptonXS().

◆ TotalXS()

void G4PolarizedAnnihilationXS::TotalXS ( )
private

Definition at line 63 of file G4PolarizedAnnihilationXS.cc.

64{
65 // total cross section is sum of
66 // + unpol xsec "sigma0"
67 // + longitudinal polarised cross section "sigma_zz" times
68 // pol_3(e-)*pol_3(e+)
69 // + transverse contribution "(sigma_xx+sigma_yy)/2" times
70 // pol_T(e-)*pol_T(e+)
71 // (Note: if both beams are transversely polarised, i.e. pol_T(e-)!=0 and
72 // pol_T(e+)!=0, and sigma_xx!=sigma_yy, then the diff. cross section
73 // will exhibit a azimuthal asymmetry even if pol_T(e-)*pol_T(e+)==0)
74}

◆ TotalXSection()

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

Reimplemented from G4VPolarizedXS.

Definition at line 220 of file G4PolarizedAnnihilationXS.cc.

224{
225 G4double totalXSFactor = pi * re2 / (gam + 1.); // atomic number ignored
226 DefineCoefficients(pol0, pol1);
227
228 G4double xs = 0.;
229
230 G4double gam2 = gam * gam;
231 G4double sqrtgam1 = std::sqrt(gam2 - 1.);
232 G4double logMEM = std::log(gam + sqrtgam1);
233 G4double unpME = (gam * (gam + 4.) + 1.) * logMEM;
234 unpME += -(gam + 3.) * sqrtgam1;
235 unpME /= 4. * (gam2 - 1.);
236 G4double longPart = (3 + gam * (gam * (gam + 1.) + 7.)) * logMEM;
237 longPart += -(5. + gam * (3 * gam + 4.)) * sqrtgam1;
238 longPart /= 4. * sqr(gam - 1.) * (gam + 1.);
239 G4double tranPart = -(5 * gam + 1.) * logMEM;
240 tranPart += (gam + 5.) * sqrtgam1;
241 tranPart /= 4. * sqr(gam - 1.) * (gam + 1.);
242
243 xs += unpME;
244 xs += polzz * longPart;
245 xs += (polxx + polyy) * tranPart;
246
247 return xs * totalXSFactor;
248}
static constexpr double pi
Definition: G4SIunits.hh:55

References DefineCoefficients(), G4InuclParticleNames::gam, pi, polxx, polyy, polzz, re2, and sqr().

Referenced by G4PolarizedAnnihilationModel::ComputeAsymmetriesPerElectron().

◆ XSection()

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

Implements G4VPolarizedXS.

Definition at line 211 of file G4PolarizedAnnihilationXS.cc.

213{
214 G4double xs = fPhi0 + pol2 * fPhi2 + pol3 * fPhi3;
215 return xs;
216}

References fPhi0, fPhi2, and fPhi3.

Field Documentation

◆ fA

G4double G4VPolarizedXS::fA
protectedinherited

Definition at line 90 of file G4VPolarizedXS.hh.

Referenced by G4VPolarizedXS::SetMaterial().

◆ fCoul

G4double G4VPolarizedXS::fCoul
protectedinherited

◆ fDice

G4double G4PolarizedAnnihilationXS::fDice
private

◆ fPhi0

G4double G4PolarizedAnnihilationXS::fPhi0
private

Definition at line 100 of file G4PolarizedAnnihilationXS.hh.

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

◆ fPhi2

G4ThreeVector G4PolarizedAnnihilationXS::fPhi2
private

◆ fPhi3

G4ThreeVector G4PolarizedAnnihilationXS::fPhi3
private

◆ fPolXS

G4double G4PolarizedAnnihilationXS::fPolXS
private

Definition at line 103 of file G4PolarizedAnnihilationXS.hh.

Referenced by G4PolarizedAnnihilationXS(), getVar(), and Initialize().

◆ fUnpXS

G4double G4PolarizedAnnihilationXS::fUnpXS
private

Definition at line 103 of file G4PolarizedAnnihilationXS.hh.

Referenced by G4PolarizedAnnihilationXS(), getVar(), 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

◆ ISPnd

G4double G4PolarizedAnnihilationXS::ISPnd
private

Definition at line 104 of file G4PolarizedAnnihilationXS.hh.

Referenced by G4PolarizedAnnihilationXS(), getVar(), and Initialize().

◆ ISPxx

G4double G4PolarizedAnnihilationXS::ISPxx
private

Definition at line 104 of file G4PolarizedAnnihilationXS.hh.

Referenced by G4PolarizedAnnihilationXS(), getVar(), and Initialize().

◆ ISPyy

G4double G4PolarizedAnnihilationXS::ISPyy
private

Definition at line 104 of file G4PolarizedAnnihilationXS.hh.

Referenced by G4PolarizedAnnihilationXS(), getVar(), and Initialize().

◆ ISPzz

G4double G4PolarizedAnnihilationXS::ISPzz
private

Definition at line 104 of file G4PolarizedAnnihilationXS.hh.

Referenced by G4PolarizedAnnihilationXS(), getVar(), and Initialize().

◆ polxx

G4double G4PolarizedAnnihilationXS::polxx
private

Definition at line 97 of file G4PolarizedAnnihilationXS.hh.

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

◆ polxy

G4double G4PolarizedAnnihilationXS::polxy
private

Definition at line 97 of file G4PolarizedAnnihilationXS.hh.

Referenced by DefineCoefficients(), and Initialize().

◆ polxz

G4double G4PolarizedAnnihilationXS::polxz
private

Definition at line 97 of file G4PolarizedAnnihilationXS.hh.

Referenced by DefineCoefficients(), and Initialize().

◆ polyx

G4double G4PolarizedAnnihilationXS::polyx
private

Definition at line 97 of file G4PolarizedAnnihilationXS.hh.

Referenced by DefineCoefficients(), and Initialize().

◆ polyy

G4double G4PolarizedAnnihilationXS::polyy
private

Definition at line 97 of file G4PolarizedAnnihilationXS.hh.

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

◆ polyz

G4double G4PolarizedAnnihilationXS::polyz
private

Definition at line 97 of file G4PolarizedAnnihilationXS.hh.

Referenced by DefineCoefficients(), and Initialize().

◆ polzx

G4double G4PolarizedAnnihilationXS::polzx
private

Definition at line 97 of file G4PolarizedAnnihilationXS.hh.

Referenced by DefineCoefficients(), and Initialize().

◆ polzy

G4double G4PolarizedAnnihilationXS::polzy
private

Definition at line 97 of file G4PolarizedAnnihilationXS.hh.

Referenced by DefineCoefficients(), and Initialize().

◆ polzz

G4double G4PolarizedAnnihilationXS::polzz
private

Definition at line 97 of file G4PolarizedAnnihilationXS.hh.

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

◆ re2

constexpr G4double G4PolarizedAnnihilationXS::re2
staticconstexprprivate
Initial value:

Definition at line 88 of file G4PolarizedAnnihilationXS.hh.

Referenced by Initialize(), and TotalXSection().


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