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

#include <G4PolarizedBremsstrahlungXS.hh>

Inheritance diagram for G4PolarizedBremsstrahlungXS:
G4VPolarizedXS

Public Member Functions

 G4PolarizedBremsstrahlungXS ()
 
 G4PolarizedBremsstrahlungXS (const G4PolarizedBremsstrahlungXS &)=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
 
G4PolarizedBremsstrahlungXSoperator= (const G4PolarizedBremsstrahlungXS &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)
 
G4double XSection (const G4StokesVector &pol2, const G4StokesVector &pol3) override
 
 ~G4PolarizedBremsstrahlungXS () 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 Attributes

G4StokesVector fFinalGammaPolarization
 
G4StokesVector fFinalLeptonPolarization
 

Static Private Attributes

static G4double SCRN [2][19]
 

Detailed Description

Definition at line 44 of file G4PolarizedBremsstrahlungXS.hh.

Constructor & Destructor Documentation

◆ G4PolarizedBremsstrahlungXS() [1/2]

G4PolarizedBremsstrahlungXS::G4PolarizedBremsstrahlungXS ( )

◆ ~G4PolarizedBremsstrahlungXS()

G4PolarizedBremsstrahlungXS::~G4PolarizedBremsstrahlungXS ( )
override

Definition at line 51 of file G4PolarizedBremsstrahlungXS.cc.

51{}

◆ G4PolarizedBremsstrahlungXS() [2/2]

G4PolarizedBremsstrahlungXS::G4PolarizedBremsstrahlungXS ( const G4PolarizedBremsstrahlungXS )
delete

Member Function Documentation

◆ GetPol2()

G4StokesVector G4PolarizedBremsstrahlungXS::GetPol2 ( )
overridevirtual

Reimplemented from G4VPolarizedXS.

Definition at line 188 of file G4PolarizedBremsstrahlungXS.cc.

189{
190 // electron/positron
192}

References fFinalLeptonPolarization.

◆ GetPol3()

G4StokesVector G4PolarizedBremsstrahlungXS::GetPol3 ( )
overridevirtual

Reimplemented from G4VPolarizedXS.

Definition at line 194 of file G4PolarizedBremsstrahlungXS.cc.

195{
196 // photon
198}

References fFinalGammaPolarization.

◆ 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 G4PolarizedBremsstrahlungXS::Initialize ( G4double  eps,
G4double  X,
G4double  phi,
const G4StokesVector p0,
const G4StokesVector p1,
G4int  flag = 0 
)
overridevirtual

Implements G4VPolarizedXS.

Definition at line 53 of file G4PolarizedBremsstrahlungXS.cc.

58{
59 G4double aLept1E = aLept0E - aGammaE;
60
61 G4double Stokes_S1 = beamPol.x();
62 G4double Stokes_S2 = beamPol.y();
63 G4double Stokes_S3 = beamPol.z();
64
65 G4double Lept0E = aLept0E / electron_mass_c2 + 1.;
66 G4double Lept0E2 = Lept0E * Lept0E;
67 G4double GammaE = aGammaE / electron_mass_c2;
68 G4double GammaE2 = GammaE * GammaE;
69 G4double Lept1E = aLept1E / electron_mass_c2 + 1.;
70 G4double Lept1E2 = Lept1E * Lept1E;
71
72 // ******* Gamma Transverse Momentum
73 G4double TMom = std::sqrt(Lept0E2 - 1.) * sintheta;
74 G4double u = TMom;
75 G4double u2 = u * u;
76 G4double Xsi = 1. / (1. + u2);
77 G4double Xsi2 = Xsi * Xsi;
78
79 G4double delta =
80 12. * std::pow(fZ, 1. / 3.) * Lept0E * Lept1E * Xsi / (121. * GammaE);
81
82 G4double GG = 0.;
83 if(delta < 0.5)
84 {
85 GG = std::log(2. * Lept0E * Lept1E / GammaE) - 2. - fCoul;
86 }
87 else if(delta < 120)
88 {
89 for(G4int j = 1; j < 19; ++j)
90 {
91 if(SCRN[0][j] >= delta)
92 {
93 GG = std::log(2 * Lept0E * Lept1E / GammaE) - 2 - fCoul -
94 (SCRN[1][j - 1] + (delta - SCRN[0][j - 1]) *
95 (SCRN[1][j] - SCRN[1][j - 1]) /
96 (SCRN[0][j] - SCRN[0][j - 1]));
97 break;
98 }
99 }
100 }
101 else
102 {
103 G4double alpha_sc = (111. * std::pow(fZ, -1. / 3.)) / Xsi;
104 GG = std::log(alpha_sc) - 2. - fCoul;
105 }
106
107 if(GG < -1.)
108 GG = -1.;
109
110 G4double I_Lept = (Lept0E2 + Lept1E2) * (3. + 2. * GG) -
111 2 * Lept0E * Lept1E * (1. + 4. * u2 * Xsi2 * GG);
112 G4double F_Lept =
113 Lept1E * 4. * GammaE * u * Xsi * (1. - 2 * Xsi) * GG / I_Lept;
114 G4double E_Lept =
115 Lept0E * 4. * GammaE * u * Xsi * (2. * Xsi - 1.) * GG / I_Lept;
116 G4double M_Lept =
117 4. * Lept0E * Lept1E * (1. + GG - 2. * Xsi2 * u2 * GG) / I_Lept;
118 G4double P_Lept =
119 GammaE2 * (1. + 8. * GG * (Xsi - 0.5) * (Xsi - 0.5)) / I_Lept;
120
121 G4double Stokes_SS1 = M_Lept * Stokes_S1 + E_Lept * Stokes_S3;
122 G4double Stokes_SS2 = M_Lept * Stokes_S2;
123 G4double Stokes_SS3 = (M_Lept + P_Lept) * Stokes_S3 + F_Lept * Stokes_S1;
124
125 fFinalLeptonPolarization.setX(Stokes_SS1);
126 fFinalLeptonPolarization.setY(Stokes_SS2);
127 fFinalLeptonPolarization.setZ(Stokes_SS3);
128
130 {
132 ed << " WARNING in pol-brem fFinalLeptonPolarization \n";
133 ed << "\t" << fFinalLeptonPolarization << "\t GG\t" << GG << "\t delta\t"
134 << delta;
135 G4Exception("G4PolarizedBremsstrahlungXS::Initialize", "pol014",
136 JustWarning, ed);
139 fFinalLeptonPolarization.setZ(Stokes_SS3);
140 if(Stokes_SS3 > 1)
142 }
143
144 G4double I_Gamma = (Lept0E2 + Lept1E2) * (3. + 2. * GG) -
145 2. * Lept0E * Lept1E * (1. + 4. * u2 * Xsi2 * GG);
146 G4double D_Gamma = 8. * Lept0E * Lept1E * u2 * Xsi2 * GG / I_Gamma;
147 G4double L_Gamma = GammaE *
148 ((Lept0E + Lept1E) * (3. + 2. * GG) -
149 2. * Lept1E * (1. + 4. * u2 * Xsi2 * GG)) /
150 I_Gamma;
151 G4double T_Gamma =
152 4. * GammaE * Lept1E * Xsi * u * (2. * Xsi - 1.) * GG / I_Gamma;
153
154 G4double Stokes_P1 = D_Gamma;
155 G4double Stokes_P2 = 0.;
156 G4double Stokes_P3 = (Stokes_S3 * L_Gamma + Stokes_S1 * T_Gamma);
157
159
160 fFinalGammaPolarization.setX(Stokes_P1);
161 fFinalGammaPolarization.setY(Stokes_P2);
162 fFinalGammaPolarization.setZ(Stokes_P3);
163
165 {
167 ed << " WARNING in pol-brem fFinalGammaPolarization \n";
168 ed << "\t" << fFinalGammaPolarization << "\t GG\t" << GG << "\t delta\t"
169 << delta;
170 G4Exception("G4PolarizedBremsstrahlungXS::Initialize", "pol015",
171 JustWarning, ed);
172 }
173}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
void setY(double)
double mag2() const
void setZ(double)
void setX(double)
float electron_mass_c2
Definition: hepunit.py:273

References source.hepunit::electron_mass_c2, G4VPolarizedXS::fCoul, fFinalGammaPolarization, fFinalLeptonPolarization, G4VPolarizedXS::fZ, G4Exception(), JustWarning, CLHEP::Hep3Vector::mag2(), SCRN, G4StokesVector::SetPhoton(), CLHEP::Hep3Vector::setX(), CLHEP::Hep3Vector::setY(), CLHEP::Hep3Vector::setZ(), CLHEP::Hep3Vector::x(), CLHEP::Hep3Vector::y(), and CLHEP::Hep3Vector::z().

◆ operator=()

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

◆ TotalXSection()

G4double G4VPolarizedXS::TotalXSection ( G4double  xmin,
G4double  xmax,
G4double  y,
const G4StokesVector pol0,
const G4StokesVector pol1 
)
virtualinherited

Reimplemented in G4PolarizedAnnihilationXS, G4PolarizedComptonXS, G4PolarizedIonisationBhabhaXS, and G4PolarizedIonisationMollerXS.

Definition at line 86 of file G4VPolarizedXS.cc.

89{
91 ed << "WARNING virtual function G4VPolarizedXS::TotalXSection() "
92 "called.\n";
93 G4Exception("G4VPolarizedXS::TotalXSection", "pol032", FatalException, ed);
94 return 0.;
95}
@ FatalException

References FatalException, and G4Exception().

Referenced by G4PolarizedIonisationModel::ComputeCrossSectionPerElectron().

◆ XSection()

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

Implements G4VPolarizedXS.

Definition at line 175 of file G4PolarizedBremsstrahlungXS.cc.

177{
179 ed << "ERROR dummy routine G4PolarizedBremsstrahlungXS::XSection "
180 "called.\n";
181 G4Exception("G4PolarizedBremsstrahlungXS::XSection", "pol016", FatalException,
182 ed);
183
184 return 0.;
185}

References FatalException, and G4Exception().

Field Documentation

◆ fA

G4double G4VPolarizedXS::fA
protectedinherited

Definition at line 90 of file G4VPolarizedXS.hh.

Referenced by G4VPolarizedXS::SetMaterial().

◆ fCoul

G4double G4VPolarizedXS::fCoul
protectedinherited

◆ fFinalGammaPolarization

G4StokesVector G4PolarizedBremsstrahlungXS::fFinalGammaPolarization
private

Definition at line 68 of file G4PolarizedBremsstrahlungXS.hh.

Referenced by G4PolarizedBremsstrahlungXS(), GetPol3(), and Initialize().

◆ fFinalLeptonPolarization

G4StokesVector G4PolarizedBremsstrahlungXS::fFinalLeptonPolarization
private

Definition at line 67 of file G4PolarizedBremsstrahlungXS.hh.

Referenced by G4PolarizedBremsstrahlungXS(), GetPol2(), 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

◆ SCRN

G4double G4PolarizedBremsstrahlungXS::SCRN
staticprivate
Initial value:
= {
{ 0.5, 1.0, 2.0, 4.0, 8.0, 15.0, 20.0, 25.0, 30.0, 35.0, 40.0, 45.0, 50.0,
60.0, 70.0, 80.0, 90.0, 100.0, 120.0 },
{ 0.0145, 0.0490, 0.1400, 0.3312, 0.6758, 1.126, 1.367, 1.564, 1.731, 1.875,
2.001, 2.114, 2.216, 2.393, 2.545, 2.676, 2.793, 2.897, 3.078 }
}

Definition at line 65 of file G4PolarizedBremsstrahlungXS.hh.

Referenced by Initialize().


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