Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4PolarizedPairProductionCrossSection.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: G4PolarizedPairProductionCrossSection.cc 96114 2016-03-16 18:51:33Z gcosmo $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class file
31 //
32 //
33 // File name: G4PolarizedPairProductionCrossSection
34 //
35 // Author: Andreas Schaelicke on the base of Karim Laihems code
36 //
37 // Creation date: 16.08.2006
38 //
39 
41 #include "G4PhysicalConstants.hh"
42 
45 // screening function lookup table;
46 
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,
82  G4double sintheta,
83  const G4StokesVector & beamPol,
84  const G4StokesVector & /*p1*/,
85  G4int /*flag*/)
86 {
87  G4double aLept1E = aGammaE - aLept0E;
88 
89  G4double Stokes_P3 = beamPol.z() ;
90  // **************************************************************************
91 
92  G4double m0_c2 = electron_mass_c2;
93  G4double Lept0E = aLept0E/m0_c2+1., Lept0E2 = Lept0E * Lept0E ;
94  G4double GammaE = aGammaE/m0_c2;
95  G4double Lept1E = aLept1E/m0_c2-1., Lept1E2 = Lept1E * Lept1E ;
96 
97 
98  // const G4Element* theSelectedElement = theModel->SelectedAtom();
99 
100  // ******* Gamma Transvers Momentum
101 
102  G4double TMom = std::sqrt(Lept0E2 -1.)* sintheta, u = TMom , u2 =u * u ;
103  G4double Xsi = 1./(1.+u2) , Xsi2 = Xsi * Xsi ;
104 
105  // G4double theZ = theSelectedElement->GetZ();
106 
107  // G4double fCoul = theSelectedElement->GetfCoulomb();
108  G4double delta = 12. * std::pow(theZ, 1./3.) * Lept0E * Lept1E * Xsi / (121. * GammaE);
109  G4double GG=0.;
110 
111  if(delta < 0.5) {
112  GG = std::log(2.* Lept0E * Lept1E / GammaE) - 2. - fCoul;
113  }
114  else if ( delta < 120.) {
115  for (G4int j=2; j<=19; j++) {
116  if(SCRN[1][j] >= delta) {
117  GG =std::log(2. * Lept0E * Lept1E / GammaE) - 2. - fCoul
118  -(SCRN[2][j-1]+(delta-SCRN[1][j-1])*(SCRN[2][j]-SCRN[2][j-1])/(SCRN[1][j]-SCRN[1][j-1]));
119  break;
120  }
121  }
122  }
123  else {
124  G4double alpha_sc = (111. * std::pow(theZ, -1./3.)) / Xsi;
125  GG = std::log(alpha_sc)- 2. - fCoul;
126  }
127 
128  if(GG<-1.) GG=-1.; // *KL* do we need this ?!
129 
130 
131  G4double I_Lepton = (Lept0E2 + Lept1E2)*(3+2*GG) + 2. * Lept0E * Lept1E * (1. + 4. * u2 * Xsi2 * GG);
132 
133  // G4double D_Lepton1 = -8 * Lept0E * Lept1E * u2 * Xsi2 * GG / I_Lepton;
134 
135  G4double L_Lepton1 = GammaE * ((Lept0E - Lept1E) * (3. + 2. * GG)+2 * Lept1E * (1. + 4. * u2 * Xsi2 * GG))/I_Lepton;
136 
137  G4double T_Lepton1 = 4. * GammaE * Lept1E * Xsi * u * (1. - 2. * Xsi) * GG / I_Lepton ;
138 
139 
140  G4double Stokes_S1 = (Stokes_P3 * T_Lepton1) ;
141  G4double Stokes_S2 = 0.;
142  G4double Stokes_S3 = (Stokes_P3 * L_Lepton1) ;
143 
144 
148 
150  G4cout<<" WARNING in pol-conv theFinalElectronPolarization \n";
151  G4cout
153  <<"\t GG\t"<<GG
154  <<"\t delta\t"<<delta
155  <<G4endl;
159  if(Stokes_S3>1.) theFinalElectronPolarization.setZ(1.);
160  }
161 
162 
163  G4double L_Lepton2 = GammaE * ((Lept1E - Lept0E) * (3. + 2. * GG)+2 * Lept0E * (1. + 4. * u2 * Xsi2 * GG))/I_Lepton;
164 
165  G4double T_Lepton2 = 4. * GammaE * Lept0E * Xsi * u * (1. - 2. * Xsi) * GG / I_Lepton ;
166 
167  G4double Stokes_SS1 = (Stokes_P3 * T_Lepton2) ;
168  G4double Stokes_SS2 = 0.;
169  G4double Stokes_SS3 = (Stokes_P3 * L_Lepton2) ;
170 
172 
176 
178  G4cout<<" WARNING in pol-conv theFinalPositronPolarization \n";
179  G4cout
181  <<"\t GG\t"<<GG
182  <<"\t delta\t"<<delta
183  <<G4endl;
184  }
185 }
186 
188  const G4StokesVector & /*pol3*/)
189 {
190  G4cout<<"ERROR dummy routine G4PolarizedPairProductionCrossSection::XSection called \n";
191  return 0.;
192 }
193 
194  // return expected mean polarisation
196 {
197  // electron/positron
199 }
201 {
202  // photon
204 }
205 
206 
#define G4endl
Definition: G4ios.hh:61
double z() const
void setX(double)
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
void setZ(double)
static constexpr double electron_mass_c2
virtual G4double XSection(const G4StokesVector &pol2, const G4StokesVector &pol3) override
double mag2() const
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
virtual void Initialize(G4double eps, G4double X, G4double phi, const G4StokesVector &p0, const G4StokesVector &p1, G4int flag=0) override
void setY(double)