Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4PolarizedBremsstrahlungCrossSection.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//
29// GEANT4 Class file
30//
31//
32// File name: G4PolarizedBremsstrahlungCrossSection
33//
34// Author: Andreas Schaelicke on the base of Karim Laihems code
35//
36// Creation date: 16.08.2006
37//
38
41
42G4bool G4PolarizedBremsstrahlungCrossSection::scrnInitialized=false;
43G4double G4PolarizedBremsstrahlungCrossSection::SCRN [3][20];
44// screening function lookup table;
45
46
47void G4PolarizedBremsstrahlungCrossSection::InitializeMe()
48{
49 if (!scrnInitialized) {
50 SCRN [1][1]= 0.5 ; SCRN [2][1] = 0.0145;
51 SCRN [1][2]= 1.0 ; SCRN [2][2] = 0.0490;
52 SCRN [1][3]= 2.0 ; SCRN [2][3] = 0.1400;
53 SCRN [1][4]= 4.0 ; SCRN [2][4] = 0.3312;
54 SCRN [1][5]= 8.0 ; SCRN [2][5] = 0.6758;
55 SCRN [1][6]= 15.0 ; SCRN [2][6] = 1.126;
56 SCRN [1][7]= 20.0 ; SCRN [2][7] = 1.367;
57 SCRN [1][8]= 25.0 ; SCRN [2][8] = 1.564;
58 SCRN [1][9]= 30.0 ; SCRN [2][9] = 1.731;
59 SCRN [1][10]= 35.0 ; SCRN [2][10]= 1.875;
60 SCRN [1][11]= 40.0 ; SCRN [2][11]= 2.001;
61 SCRN [1][12]= 45.0 ; SCRN [2][12]= 2.114;
62 SCRN [1][13]= 50.0 ; SCRN [2][13]= 2.216;
63 SCRN [1][14]= 60.0 ; SCRN [2][14]= 2.393;
64 SCRN [1][15]= 70.0 ; SCRN [2][15]= 2.545;
65 SCRN [1][16]= 80.0 ; SCRN [2][16]= 2.676;
66 SCRN [1][17]= 90.0 ; SCRN [2][17]= 2.793;
67 SCRN [1][18]= 100.0 ; SCRN [2][18]= 2.897;
68 SCRN [1][19]= 120.0 ; SCRN [2][19]= 3.078;
69
70 scrnInitialized=true;
71 }
72}
73
75{
76 InitializeMe();
77}
78
79
81 G4double aLept0E, G4double aGammaE, G4double sintheta,
82 const G4StokesVector & beamPol,
83 const G4StokesVector & /*p1*/,
84 G4int /*flag*/)
85{
86// G4cout<<"G4PolarizedBremsstrahlungCrossSection::Initialize \n"
87// <<"lepE = "<<aLept0E
88// <<"gamE = "<<aGammaE
89// <<"sint = "<<sintheta<<"\n"
90// <<"beamPol="<<beamPol<<"\n";
91
92 G4double aLept1E = aLept0E - aGammaE;
93
94 G4double Stokes_S1 = beamPol.x() ;
95 G4double Stokes_S2 = beamPol.y() ;
96 G4double Stokes_S3 = beamPol.z() ;
97 // **************************************************************************
98
99 G4double m0_c2 = electron_mass_c2;
100 G4double Lept0E = aLept0E/m0_c2+1., Lept0E2 = Lept0E * Lept0E ;
101 G4double GammaE = aGammaE/m0_c2, GammaE2 = GammaE * GammaE ;
102 G4double Lept1E = aLept1E/m0_c2+1., Lept1E2 = Lept1E * Lept1E ;
103
104
105 // const G4Element* theSelectedElement = theModel->SelectedAtom();
106
107 // ******* Gamma Transvers Momentum
108
109 G4double TMom = std::sqrt(Lept0E2 -1.)* sintheta;
110 G4double u = TMom , u2 =u * u ;
111 G4double Xsi = 1./(1.+u2) , Xsi2 = Xsi * Xsi ;
112
113 // G4double theZ = theSelectedElement->GetZ();
114
115 // G4double fCoul = theSelectedElement->GetfCoulomb();
116 G4double delta = 12. * std::pow(theZ, 1./3.) *
117 Lept0E * Lept1E * Xsi / (121. * GammaE);
118 G4double GG=0.;
119
120 if(delta < 0.5) {
121 GG = std::log(2.* Lept0E * Lept1E / GammaE) - 2. - fCoul;
122 }
123 else if ( delta < 120) {
124 for (G4int j=2; j<=19; j++) {
125 if(SCRN[1][j] >= delta) {
126 GG =std::log(2 * Lept0E * Lept1E / GammaE) - 2 - fCoul
127 -(SCRN[2][j-1]+(delta-SCRN[1][j-1])*(SCRN[2][j]-SCRN[2][j-1])
128 /(SCRN[1][j]-SCRN[1][j-1]));
129 break;
130 }
131 }
132 }
133 else {
134 G4double alpha_sc = (111. * std::pow(theZ, -1./3.)) / Xsi;
135 GG = std::log(alpha_sc)- 2. - fCoul;
136 }
137
138 if(GG<-1.) GG=-1.; // *KL* do we need this ?!
139
140 G4double I_Lept = (Lept0E2 + Lept1E2) * (3.+2.*GG) - 2 * Lept0E * Lept1E * (1. + 4. * u2 * Xsi2 * GG);
141 G4double F_Lept = Lept1E * 4. * GammaE * u * Xsi * (1. - 2 * Xsi) * GG / I_Lept;
142 G4double E_Lept = Lept0E * 4. * GammaE * u * Xsi * (2. * Xsi - 1.) * GG / I_Lept;
143 G4double M_Lept = 4. * Lept0E * Lept1E * (1. + GG - 2. * Xsi2 * u2 * GG) / I_Lept ;
144 G4double P_Lept = GammaE2 * (1. + 8. * GG * (Xsi - 0.5)*(Xsi - 0.5)) / I_Lept ;
145
146 G4double Stokes_SS1 = M_Lept * Stokes_S1 + E_Lept * Stokes_S3;
147 G4double Stokes_SS2 = M_Lept * Stokes_S2 ;
148 G4double Stokes_SS3 = (M_Lept + P_Lept) * Stokes_S3 + F_Lept * Stokes_S1;
149
150 theFinalLeptonPolarization.setX(Stokes_SS1);
151 theFinalLeptonPolarization.setY(Stokes_SS2);
152 theFinalLeptonPolarization.setZ(Stokes_SS3);
153
154 if(theFinalLeptonPolarization.mag2()>1) {
155 G4cout<<" WARNING in pol-brem theFinalLeptonPolarization \n";
156 G4cout
157 <<"\t"<<theFinalLeptonPolarization
158 <<"\t GG\t"<<GG
159 <<"\t delta\t"<<delta
160 <<G4endl;
161 theFinalLeptonPolarization.setX(0);
162 theFinalLeptonPolarization.setY(0);
163 theFinalLeptonPolarization.setZ(Stokes_SS3);
164 if(Stokes_SS3>1) theFinalLeptonPolarization.setZ(1);
165 }
166
167
168 G4double I_Gamma = (Lept0E2 + Lept1E2)*(3.+2.*GG) - 2. * Lept0E * Lept1E * (1. + 4. * u2 * Xsi2 * GG);
169 G4double D_Gamma = 8. * Lept0E * Lept1E * u2 * Xsi2 * GG / I_Gamma;
170 G4double L_Gamma = GammaE * ((Lept0E + Lept1E) * (3. + 2. * GG)
171 - 2. * Lept1E * (1. + 4. * u2 * Xsi2 * GG))/I_Gamma;
172 G4double T_Gamma = 4. * GammaE * Lept1E * Xsi * u * (2. * Xsi - 1.) * GG / I_Gamma ;
173
174 G4double Stokes_P1 = D_Gamma ;
175 G4double Stokes_P2 = 0.;
176 G4double Stokes_P3 = (Stokes_S3*L_Gamma + Stokes_S1*T_Gamma) ;
177
178 theFinalGammaPolarization.SetPhoton();
179
180 theFinalGammaPolarization.setX(Stokes_P1);
181 theFinalGammaPolarization.setY(Stokes_P2);
182 theFinalGammaPolarization.setZ(Stokes_P3);
183
184 if(theFinalGammaPolarization.mag2()>1) {
185 G4cout<<" WARNING in pol-brem theFinalGammaPolarization \n";
186 G4cout
187 <<"\t"<<theFinalGammaPolarization
188 <<"\t GG\t"<<GG
189 <<"\t delta\t"<<delta
190 <<G4endl;
191 }
192}
193
195 const G4StokesVector & /*pol3*/)
196{
197 G4cout<<"ERROR dummy routine G4PolarizedBremsstrahlungCrossSection::XSection called \n";
198 return 0.;
199}
200
201 // return expected mean polarisation
203{
204 // electron/positron
205 return theFinalLeptonPolarization;
206}
208{
209 // photon
210 return theFinalGammaPolarization;;
211}
212
213
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
double z() const
double x() const
void setY(double)
double mag2() const
double y() const
void setZ(double)
void setX(double)
virtual void Initialize(G4double eps, G4double X, G4double phi, const G4StokesVector &p0, const G4StokesVector &p1, G4int flag=0) override
virtual G4double XSection(const G4StokesVector &pol2, const G4StokesVector &pol3) override