Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4PolarizedGammaConversionXS Class Reference

#include <G4PolarizedGammaConversionXS.hh>

+ Inheritance diagram for G4PolarizedGammaConversionXS:

Public Member Functions

 G4PolarizedGammaConversionXS ()
 
 ~G4PolarizedGammaConversionXS () override
 
void Initialize (G4double eps, G4double X, G4double phi, const G4StokesVector &p0, const G4StokesVector &p1, G4int flag=0) override
 
G4double XSection (const G4StokesVector &pol2, const G4StokesVector &pol3) override
 
G4StokesVector GetPol2 () override
 
G4StokesVector GetPol3 () override
 
G4PolarizedGammaConversionXSoperator= (const G4PolarizedGammaConversionXS &right)=delete
 
 G4PolarizedGammaConversionXS (const G4PolarizedGammaConversionXS &)=delete
 
- Public Member Functions inherited from G4VPolarizedXS
 G4VPolarizedXS ()
 
virtual ~G4VPolarizedXS ()
 
virtual G4double TotalXSection (G4double xmin, G4double xmax, G4double y, const G4StokesVector &pol0, const G4StokesVector &pol1)
 
G4double GetYmin ()
 
virtual G4double GetXmin (G4double y)
 
virtual G4double GetXmax (G4double y)
 
void SetMaterial (G4double A, G4double Z, G4double coul)
 
G4VPolarizedXSoperator= (const G4VPolarizedXS &right)=delete
 
 G4VPolarizedXS (const G4VPolarizedXS &)=delete
 

Additional Inherited Members

- Protected Member Functions inherited from G4VPolarizedXS
void SetXmin (G4double xmin)
 
void SetXmax (G4double xmax)
 
void SetYmin (G4double ymin)
 
- Protected Attributes inherited from G4VPolarizedXS
G4double fXmin
 
G4double fXmax
 
G4double fYmin
 
G4double fA
 
G4double fZ
 
G4double fCoul
 

Detailed Description

Definition at line 38 of file G4PolarizedGammaConversionXS.hh.

Constructor & Destructor Documentation

◆ G4PolarizedGammaConversionXS() [1/2]

G4PolarizedGammaConversionXS::G4PolarizedGammaConversionXS ( )

Definition at line 43 of file G4PolarizedGammaConversionXS.cc.

44{
45 fFinalElectronPolarization = G4StokesVector::ZERO;
46 fFinalPositronPolarization = G4StokesVector::ZERO;
47}
static const G4StokesVector ZERO

◆ ~G4PolarizedGammaConversionXS()

G4PolarizedGammaConversionXS::~G4PolarizedGammaConversionXS ( )
override

Definition at line 49 of file G4PolarizedGammaConversionXS.cc.

49{}

◆ G4PolarizedGammaConversionXS() [2/2]

G4PolarizedGammaConversionXS::G4PolarizedGammaConversionXS ( const G4PolarizedGammaConversionXS & )
delete

Member Function Documentation

◆ GetPol2()

G4StokesVector G4PolarizedGammaConversionXS::GetPol2 ( )
overridevirtual

Reimplemented from G4VPolarizedXS.

Definition at line 176 of file G4PolarizedGammaConversionXS.cc.

177{
178 // electron/positron
179 return fFinalElectronPolarization;
180}

Referenced by G4PolarizedGammaConversionModel::SampleSecondaries().

◆ GetPol3()

G4StokesVector G4PolarizedGammaConversionXS::GetPol3 ( )
overridevirtual

Reimplemented from G4VPolarizedXS.

Definition at line 182 of file G4PolarizedGammaConversionXS.cc.

183{
184 // photon
185 return fFinalPositronPolarization;
186}

Referenced by G4PolarizedGammaConversionModel::SampleSecondaries().

◆ Initialize()

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

Implements G4VPolarizedXS.

Definition at line 51 of file G4PolarizedGammaConversionXS.cc.

54{
55 G4double aLept1E = aGammaE - aLept0E;
56
57 G4double Stokes_P3 = beamPol.z();
58
59 G4double Lept0E = aLept0E / CLHEP::electron_mass_c2 + 1.;
60 G4double Lept0E2 = Lept0E * Lept0E;
61 G4double GammaE = aGammaE / CLHEP::electron_mass_c2;
62 G4double Lept1E = aLept1E / CLHEP::electron_mass_c2 - 1.;
63 G4double Lept1E2 = Lept1E * Lept1E;
64
65 // ******* Gamma Transvers Momentum
66 G4double TMom = std::sqrt(Lept0E2 - 1.) * sintheta;
67 G4double u = TMom;
68 G4double u2 = u * u;
69 G4double Xsi = 1. / (1. + u2);
70 G4double Xsi2 = Xsi * Xsi;
71
72 G4double delta =
73 12. * std::pow(fZ, 1. / 3.) * Lept0E * Lept1E * Xsi / (121. * GammaE);
74 G4double GG = 0.;
75
76 if(delta < 0.5)
77 {
78 GG = std::log(2. * Lept0E * Lept1E / GammaE) - 2. - fCoul;
79 }
80 else if(delta < 120.)
81 {
82 for(G4int j = 1; j < 19; ++j)
83 {
84 if(SCRN[0][j] >= delta)
85 {
86 GG = std::log(2. * Lept0E * Lept1E / GammaE) - 2. - fCoul -
87 (SCRN[1][j - 1] + (delta - SCRN[0][j - 1]) *
88 (SCRN[1][j] - SCRN[1][j - 1]) /
89 (SCRN[0][j] - SCRN[0][j - 1]));
90 break;
91 }
92 }
93 }
94 else
95 {
96 G4double alpha_sc = (111. * std::pow(fZ, -1. / 3.)) / Xsi;
97 GG = std::log(alpha_sc) - 2. - fCoul;
98 }
99
100 if(GG < -1.)
101 GG = -1.;
102
103 G4double I_Lepton = (Lept0E2 + Lept1E2) * (3 + 2 * GG) +
104 2. * Lept0E * Lept1E * (1. + 4. * u2 * Xsi2 * GG);
105
106 G4double L_Lepton1 = GammaE *
107 ((Lept0E - Lept1E) * (3. + 2. * GG) +
108 2 * Lept1E * (1. + 4. * u2 * Xsi2 * GG)) /
109 I_Lepton;
110
111 G4double T_Lepton1 =
112 4. * GammaE * Lept1E * Xsi * u * (1. - 2. * Xsi) * GG / I_Lepton;
113
114 G4double Stokes_S1 = (Stokes_P3 * T_Lepton1);
115 G4double Stokes_S2 = 0.;
116 G4double Stokes_S3 = (Stokes_P3 * L_Lepton1);
117
118 fFinalElectronPolarization.setX(Stokes_S1);
119 fFinalElectronPolarization.setY(Stokes_S2);
120 fFinalElectronPolarization.setZ(Stokes_S3);
121
122 if(fFinalElectronPolarization.mag2() > 1.)
123 {
125 ed << "\t" << fFinalElectronPolarization << "\t GG\t" << GG << "\t delta\t"
126 << delta << "\n";
127 G4Exception("G4PolarizedGammaConversionXS::Initialize", "pol022",
128 JustWarning, ed);
129 fFinalElectronPolarization.setX(0.);
130 fFinalElectronPolarization.setY(0.);
131 fFinalElectronPolarization.setZ(Stokes_S3);
132 if(Stokes_S3 > 1.)
133 fFinalElectronPolarization.setZ(1.);
134 }
135
136 G4double L_Lepton2 = GammaE *
137 ((Lept1E - Lept0E) * (3. + 2. * GG) +
138 2 * Lept0E * (1. + 4. * u2 * Xsi2 * GG)) /
139 I_Lepton;
140
141 G4double T_Lepton2 =
142 4. * GammaE * Lept0E * Xsi * u * (1. - 2. * Xsi) * GG / I_Lepton;
143
144 G4double Stokes_SS1 = (Stokes_P3 * T_Lepton2);
145 G4double Stokes_SS2 = 0.;
146 G4double Stokes_SS3 = (Stokes_P3 * L_Lepton2);
147
148 fFinalPositronPolarization.SetPhoton();
149
150 fFinalPositronPolarization.setX(Stokes_SS1);
151 fFinalPositronPolarization.setY(Stokes_SS2);
152 fFinalPositronPolarization.setZ(Stokes_SS3);
153
154 if(fFinalPositronPolarization.mag2() > 1.)
155 {
157 ed << "\t" << fFinalPositronPolarization << "\t GG\t" << GG << "\t delta\t"
158 << delta << "\n";
159 G4Exception("G4PolarizedGammaConversionXS::Initialize", "pol023",
160 JustWarning, ed);
161 }
162}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
void setY(double)
double mag2() const
void setZ(double)
void setX(double)

Referenced by G4PolarizedGammaConversionModel::SampleSecondaries().

◆ operator=()

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

◆ XSection()

G4double G4PolarizedGammaConversionXS::XSection ( const G4StokesVector & pol2,
const G4StokesVector & pol3 )
overridevirtual

Implements G4VPolarizedXS.

Definition at line 164 of file G4PolarizedGammaConversionXS.cc.

166{
168 ed << "ERROR dummy routine G4PolarizedGammaConversionXS::XSection "
169 "called \n";
170 G4Exception("G4PolarizedGammaConversionXS::Initialize", "pol024",
171 FatalException, ed);
172 return 0.;
173}
@ FatalException

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