Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4PolarizedPairProductionCrossSection Class Reference

#include <G4PolarizedPairProductionCrossSection.hh>

+ Inheritance diagram for G4PolarizedPairProductionCrossSection:

Public Member Functions

 G4PolarizedPairProductionCrossSection ()
 
virtual void Initialize (G4double eps, G4double X, G4double phi, const G4StokesVector &p0, const G4StokesVector &p1, G4int flag=0)
 
virtual G4double XSection (const G4StokesVector &pol2, const G4StokesVector &pol3)
 
G4StokesVector GetPol2 ()
 
G4StokesVector GetPol3 ()
 
- Public Member Functions inherited from G4VPolarizedCrossSection
 G4VPolarizedCrossSection ()
 
virtual ~G4VPolarizedCrossSection ()
 
virtual void Initialize (G4double, G4double, G4double, const G4StokesVector &p0, const G4StokesVector &p1, G4int flag=0)
 
virtual G4double XSection (const G4StokesVector &pol2, const G4StokesVector &pol3)=0
 
virtual G4double TotalXSection (G4double xmin, G4double xmax, G4double y, const G4StokesVector &pol0, const G4StokesVector &pol1)
 
virtual G4StokesVector GetPol2 ()
 
virtual G4StokesVector GetPol3 ()
 
G4double GetYmin ()
 
virtual G4double GetXmin (G4double y)
 
virtual G4double GetXmax (G4double y)
 
void SetMaterial (G4double A, G4double Z, G4double coul)
 

Additional Inherited Members

- Protected Member Functions inherited from G4VPolarizedCrossSection
void SetXmin (G4double xmin)
 
void SetXmax (G4double xmax)
 
void SetYmin (G4double ymin)
 
- Protected Attributes inherited from G4VPolarizedCrossSection
G4double fXmin
 
G4double fXmax
 
G4double fYmin
 
G4double theA
 
G4double theZ
 
G4double fCoul
 

Detailed Description

Definition at line 53 of file G4PolarizedPairProductionCrossSection.hh.

Constructor & Destructor Documentation

◆ G4PolarizedPairProductionCrossSection()

G4PolarizedPairProductionCrossSection::G4PolarizedPairProductionCrossSection ( )

Definition at line 74 of file G4PolarizedPairProductionCrossSection.cc.

75{
76 InitializeMe();
77}

Member Function Documentation

◆ GetPol2()

G4StokesVector G4PolarizedPairProductionCrossSection::GetPol2 ( )
virtual

Reimplemented from G4VPolarizedCrossSection.

Definition at line 195 of file G4PolarizedPairProductionCrossSection.cc.

196{
197 // electron/positron
198 return theFinalElectronPolarization;
199}

◆ GetPol3()

G4StokesVector G4PolarizedPairProductionCrossSection::GetPol3 ( )
virtual

Reimplemented from G4VPolarizedCrossSection.

Definition at line 200 of file G4PolarizedPairProductionCrossSection.cc.

201{
202 // photon
203 return theFinalPositronPolarization;;
204}

◆ Initialize()

void G4PolarizedPairProductionCrossSection::Initialize ( G4double  eps,
G4double  X,
G4double  phi,
const G4StokesVector p0,
const G4StokesVector p1,
G4int  flag = 0 
)
virtual

Reimplemented from G4VPolarizedCrossSection.

Definition at line 80 of file G4PolarizedPairProductionCrossSection.cc.

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
145 theFinalElectronPolarization.setX(Stokes_S1);
146 theFinalElectronPolarization.setY(Stokes_S2);
147 theFinalElectronPolarization.setZ(Stokes_S3);
148
149 if(theFinalElectronPolarization.mag2()>1) {
150 G4cout<<" WARNING in pol-conv theFinalElectronPolarization \n";
151 G4cout
152 <<"\t"<<theFinalElectronPolarization
153 <<"\t GG\t"<<GG
154 <<"\t delta\t"<<delta
155 <<G4endl;
156 theFinalElectronPolarization.setX(0);
157 theFinalElectronPolarization.setY(0);
158 theFinalElectronPolarization.setZ(Stokes_S3);
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
171 theFinalPositronPolarization.SetPhoton();
172
173 theFinalPositronPolarization.setX(Stokes_SS1);
174 theFinalPositronPolarization.setY(Stokes_SS2);
175 theFinalPositronPolarization.setZ(Stokes_SS3);
176
177 if(theFinalPositronPolarization.mag2()>1) {
178 G4cout<<" WARNING in pol-conv theFinalPositronPolarization \n";
179 G4cout
180 <<"\t"<<theFinalPositronPolarization
181 <<"\t GG\t"<<GG
182 <<"\t delta\t"<<delta
183 <<G4endl;
184 }
185}
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
void setY(double)
double mag2() const
void setZ(double)
void setX(double)

◆ XSection()

G4double G4PolarizedPairProductionCrossSection::XSection ( const G4StokesVector pol2,
const G4StokesVector pol3 
)
virtual

Implements G4VPolarizedCrossSection.

Definition at line 187 of file G4PolarizedPairProductionCrossSection.cc.

189{
190 G4cout<<"ERROR dummy routine G4PolarizedPairProductionCrossSection::XSection called \n";
191 return 0;
192}

The documentation for this class was generated from the following files: