Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4PolarizedBremsstrahlungXS.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// -------------------------------------------------------------------
27//
28// Geant4 class file
29//
30// File name: G4PolarizedBremsstrahlungXS
31//
32// Author: Andreas Schaelicke on the base of Karim Laihems code
33
35
37
38G4double G4PolarizedBremsstrahlungXS::SCRN[2][19] = {
39 { 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,
40 60.0, 70.0, 80.0, 90.0, 100.0, 120.0 },
41 { 0.0145, 0.0490, 0.1400, 0.3312, 0.6758, 1.126, 1.367, 1.564, 1.731, 1.875,
42 2.001, 2.114, 2.216, 2.393, 2.545, 2.676, 2.793, 2.897, 3.078 }
43};
44
46{
47 fFinalLeptonPolarization = G4StokesVector::ZERO;
48 fFinalGammaPolarization = G4StokesVector::ZERO;
49}
50
52
54 G4double sintheta,
55 const G4StokesVector& beamPol,
56 const G4StokesVector& /*p1*/,
57 G4int /*flag*/)
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
129 if(fFinalLeptonPolarization.mag2() > 1.)
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);
137 fFinalLeptonPolarization.setX(0);
138 fFinalLeptonPolarization.setY(0);
139 fFinalLeptonPolarization.setZ(Stokes_SS3);
140 if(Stokes_SS3 > 1)
141 fFinalLeptonPolarization.setZ(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
158 fFinalGammaPolarization.SetPhoton();
159
160 fFinalGammaPolarization.setX(Stokes_P1);
161 fFinalGammaPolarization.setY(Stokes_P2);
162 fFinalGammaPolarization.setZ(Stokes_P3);
163
164 if(fFinalGammaPolarization.mag2() > 1.)
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}
174
176 const G4StokesVector& /*pol3*/)
177{
179 ed << "ERROR dummy routine G4PolarizedBremsstrahlungXS::XSection "
180 "called.\n";
181 G4Exception("G4PolarizedBremsstrahlungXS::XSection", "pol016", FatalException,
182 ed);
183
184 return 0.;
185}
186
187// return expected mean polarisation
189{
190 // electron/positron
191 return fFinalLeptonPolarization;
192}
193
195{
196 // photon
197 return fFinalGammaPolarization;
198}
@ JustWarning
@ FatalException
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
double z() const
double x() const
void setY(double)
double mag2() const
double y() const
void setZ(double)
void setX(double)
G4double XSection(const G4StokesVector &pol2, const G4StokesVector &pol3) override
void Initialize(G4double eps, G4double X, G4double phi, const G4StokesVector &p0, const G4StokesVector &p1, G4int flag=0) override
static const G4StokesVector ZERO