Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4PolarizedPEEffectCrossSection.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// GEANT4 Class file
28//
29//
30// File name: G4PolarizedPEEffectCrossSection
31//
32// Author: Karim Laihem
33//
34// Creation date: 15.03.2007
35//
36// Modifications:
37// 19-03-07 Modified to fit in g4.8.2 framework (A.Schaelicke)
38//
39// Class Description:
40//
43
44using namespace std;
45
46//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
48 {
49 cout<<"G4PolarizedPEEffectCrossSection() init\n";
50
51}
52
53//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
55{}
56
57//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
59 G4double aLept0E,
60 G4double sinTheta,
61 const G4StokesVector & beamPol,
62 const G4StokesVector & /*p1*/,
63 G4int /*flag*/)
64{
65 // cout<<"G4PolarizedPEEffectCrossSection::Initialize()\n";
66
67// G4StokesVector PolarizedPhotoElectricEffect::Transfer_G4StokesVector(
68// G4double aGammaE, // Incoming Primary Gamma Energy.
69// G4ThreeVector aGammaDir, // Incoming Primary Gamma Direction.
70// G4StokesVector beamPol, // Incoming Primary Gamma polarization.
71// G4double aLept0E, // The Lepton e- of interest Total energy.
72// G4ThreeVector aParticl_01_Dir, // The Lepton e- of interest direction.
73// G4double cos_aTetha_Angle // The lepton of interest Scattering angle.
74// )
75
76// ***********************************************************
77// ************ added by Karim Polarization transfer to e- in PhotoelectricEffect.
78// ************
79// ***********************************************************
80 G4double Gfactor = aLept0E/electron_mass_c2+1.;
81 G4double Gfactor_2 = Gfactor * Gfactor;
82
83 G4double BETA = sqrt(1. - 1./(Gfactor_2));
84
85 G4double Stokes_P3 = beamPol.z() ;
86
87 G4double m0_c2 = electron_mass_c2;
88 G4double Lept0E = aLept0E/m0_c2+1., Lept0E2 = Lept0E * Lept0E ;
89 G4double GammaE = aGammaE/m0_c2;
90
91
92// G4double cosTheta = cos_aTetha_Angle;
93// G4double sinTheta = sqrt(1- cos_aTetha_Angle * cos_aTetha_Angle);
94 G4double cosTheta = std::sqrt(1. - sinTheta*sinTheta);
95
96 G4double D_Lepton0 = (1./GammaE) * ((2./(GammaE*Lept0E*(1-BETA*cosTheta)))-1.);
97
98 G4double I_Lepton0 = 1.0+D_Lepton0;
99
100 G4double A_Lepton0 = (Lept0E/(Lept0E+1))*(2.0/(GammaE*Lept0E)
101 + BETA*cosTheta
102 +(2.0/((GammaE*Lept0E2)*(1.0-BETA*cosTheta)))) / I_Lepton0 ;
103
104 G4double B_Lepton0 = (Lept0E/(Lept0E+1.0)) * BETA * sinTheta * (2.0/(GammaE*Lept0E*(1-BETA*cosTheta))-1.0)/I_Lepton0;
105
106 G4double Stokes_S1 = (Stokes_P3 * B_Lepton0) ;
107 G4double Stokes_S2 = 0.;
108 G4double Stokes_S3 = (Stokes_P3 * A_Lepton0) ;
109
110
111 theFinalElectronPolarization.setX(Stokes_S1);
112 theFinalElectronPolarization.setY(Stokes_S2);
113 theFinalElectronPolarization.setZ(Stokes_S3);
114
115 if((theFinalElectronPolarization.x()*theFinalElectronPolarization.x()
116 + theFinalElectronPolarization.y()* theFinalElectronPolarization.y()
117 + theFinalElectronPolarization.z()* theFinalElectronPolarization.z())>1)
118
119 {
120 cout<<"Warning: PhotoelectricEffect Problem in pol-transfer photon to lepton:Px2 + Py2 + Pz2 > 1"<<endl;
121 cout<<"Polarization transfer forced to be total and similar as incoming Photo"<<endl;
122 // *KL* Surprising if it arrives (never seen it up to now)
123 theFinalElectronPolarization = beamPol; // suplement de securite
124// cout<<"PhotoEffect okay :"
125// <<"\t"<<(aLept0E-m0_c2)/aGammaE
126// <<"\t"<<aGammaE
127// <<"\t"<<aLept0E
128// <<"\t"<<cos_aTetha_Angle
129// <<"\t"<<beamPol
130// <<"\t"<<theFinalElectronPolarization
131// <<"\t"<<A_Lepton0
132// <<endl;
133 }
134}
135
136
137
138//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
139
141 const G4StokesVector & /*pol3*/)
142{
143 cout<<"ERROR dummy routine G4PolarizedPEEffectCrossSection::XSection() called\n";
144 return 0.;
145}
146
147//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
148
150{
151 return theFinalElectronPolarization;
152}
153
154//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
155
157{
158 return G4StokesVector();
159}
160
161
162
163//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
164
165
166//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
double z() const
double x() const
void setY(double)
double y() const
void setZ(double)
void setX(double)
G4double XSection(const G4StokesVector &pol2, const G4StokesVector &pol3) override
virtual void Initialize(G4double aGammaE, G4double aLept0E, G4double sintheta, const G4StokesVector &beamPol, const G4StokesVector &, G4int flag=0) override