Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4PolarizedPEEffectCrossSection.cc
이 파일의 문서화 페이지로 가기
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 // $Id: G4PolarizedPEEffectCrossSection.cc 96114 2016-03-16 18:51:33Z gcosmo $
27 //
28 // GEANT4 Class file
29 //
30 //
31 // File name: G4PolarizedPEEffectCrossSection
32 //
33 // Author: Karim Laihem
34 //
35 // Creation date: 15.03.2007
36 //
37 // Modifications:
38 // 19-03-07 Modified to fit in g4.8.2 framework (A.Schaelicke)
39 //
40 // Class Description:
41 //
43 #include "G4PhysicalConstants.hh"
44 
45 using namespace std;
46 
47 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
49  {
50  cout<<"G4PolarizedPEEffectCrossSection() init\n";
51 
52 }
53 
54 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
56 {}
57 
58 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
60  G4double aLept0E,
61  G4double sinTheta,
62  const G4StokesVector & beamPol,
63  const G4StokesVector & /*p1*/,
64  G4int /*flag*/)
65 {
66  // cout<<"G4PolarizedPEEffectCrossSection::Initialize()\n";
67 
68 // G4StokesVector PolarizedPhotoElectricEffect::Transfer_G4StokesVector(
69 // G4double aGammaE, // Incoming Primary Gamma Energy.
70 // G4ThreeVector aGammaDir, // Incoming Primary Gamma Direction.
71 // G4StokesVector beamPol, // Incoming Primary Gamma polarization.
72 // G4double aLept0E, // The Lepton e- of interest Total energy.
73 // G4ThreeVector aParticl_01_Dir, // The Lepton e- of interest direction.
74 // G4double cos_aTetha_Angle // The lepton of interest Scattering angle.
75 // )
76 
77 // ***********************************************************
78 // ************ added by Karim Polarization transfer to e- in PhotoelectricEffect.
79 // ************
80 // ***********************************************************
81  G4double Gfactor = aLept0E/electron_mass_c2+1.;
82  G4double Gfactor_2 = Gfactor * Gfactor;
83 
84  G4double BETA = sqrt(1. - 1./(Gfactor_2));
85 
86  G4double Stokes_P3 = beamPol.z() ;
87 
88  G4double m0_c2 = electron_mass_c2;
89  G4double Lept0E = aLept0E/m0_c2+1., Lept0E2 = Lept0E * Lept0E ;
90  G4double GammaE = aGammaE/m0_c2;
91 
92 
93 // G4double cosTheta = cos_aTetha_Angle;
94 // G4double sinTheta = sqrt(1- cos_aTetha_Angle * cos_aTetha_Angle);
95  G4double cosTheta = std::sqrt(1. - sinTheta*sinTheta);
96 
97  G4double D_Lepton0 = (1./GammaE) * ((2./(GammaE*Lept0E*(1-BETA*cosTheta)))-1.);
98 
99  G4double I_Lepton0 = 1.0+D_Lepton0;
100 
101  G4double A_Lepton0 = (Lept0E/(Lept0E+1))*(2.0/(GammaE*Lept0E)
102  + BETA*cosTheta
103  +(2.0/((GammaE*Lept0E2)*(1.0-BETA*cosTheta)))) / I_Lepton0 ;
104 
105  G4double B_Lepton0 = (Lept0E/(Lept0E+1.0)) * BETA * sinTheta * (2.0/(GammaE*Lept0E*(1-BETA*cosTheta))-1.0)/I_Lepton0;
106 
107  G4double Stokes_S1 = (Stokes_P3 * B_Lepton0) ;
108  G4double Stokes_S2 = 0.;
109  G4double Stokes_S3 = (Stokes_P3 * A_Lepton0) ;
110 
111 
112  theFinalElectronPolarization.setX(Stokes_S1);
113  theFinalElectronPolarization.setY(Stokes_S2);
114  theFinalElectronPolarization.setZ(Stokes_S3);
115 
116  if((theFinalElectronPolarization.x()*theFinalElectronPolarization.x()
117  + theFinalElectronPolarization.y()* theFinalElectronPolarization.y()
118  + theFinalElectronPolarization.z()* theFinalElectronPolarization.z())>1)
119 
120  {
121  cout<<"Warning: PhotoelectricEffect Problem in pol-transfer photon to lepton:Px2 + Py2 + Pz2 > 1"<<endl;
122  cout<<"Polarization transfer forced to be total and similar as incoming Photo"<<endl;
123  // *KL* Surprising if it arrives (never seen it up to now)
124  theFinalElectronPolarization = beamPol; // suplement de securite
125 // cout<<"PhotoEffect okay :"
126 // <<"\t"<<(aLept0E-m0_c2)/aGammaE
127 // <<"\t"<<aGammaE
128 // <<"\t"<<aLept0E
129 // <<"\t"<<cos_aTetha_Angle
130 // <<"\t"<<beamPol
131 // <<"\t"<<theFinalElectronPolarization
132 // <<"\t"<<A_Lepton0
133 // <<endl;
134  }
135 }
136 
137 
138 
139 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
140 
142  const G4StokesVector & /*pol3*/)
143 {
144  cout<<"ERROR dummy routine G4PolarizedPEEffectCrossSection::XSection() called\n";
145  return 0.;
146 }
147 
148 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
149 
151 {
152  return theFinalElectronPolarization;
153 }
154 
155 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
156 
158 {
159  return G4StokesVector();
160 }
161 
162 
163 
164 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
165 
166 
167 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
double z() const
virtual void Initialize(G4double aGammaE, G4double aLept0E, G4double sintheta, const G4StokesVector &beamPol, const G4StokesVector &, G4int flag=0) override
double G4double
Definition: G4Types.hh:76
static constexpr double electron_mass_c2
int G4int
Definition: G4Types.hh:78
G4double XSection(const G4StokesVector &pol2, const G4StokesVector &pol3) override