Geant4 10.7.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) override
 
virtual G4double XSection (const G4StokesVector &pol2, const G4StokesVector &pol3) override
 
G4StokesVector GetPol2 () override
 
G4StokesVector GetPol3 () override
 
- 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 52 of file G4PolarizedPairProductionCrossSection.hh.

Constructor & Destructor Documentation

◆ G4PolarizedPairProductionCrossSection()

G4PolarizedPairProductionCrossSection::G4PolarizedPairProductionCrossSection ( )

Definition at line 73 of file G4PolarizedPairProductionCrossSection.cc.

74{
75 InitializeMe();
76}

Member Function Documentation

◆ GetPol2()

G4StokesVector G4PolarizedPairProductionCrossSection::GetPol2 ( )
overridevirtual

Reimplemented from G4VPolarizedCrossSection.

Definition at line 194 of file G4PolarizedPairProductionCrossSection.cc.

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

◆ GetPol3()

G4StokesVector G4PolarizedPairProductionCrossSection::GetPol3 ( )
overridevirtual

Reimplemented from G4VPolarizedCrossSection.

Definition at line 199 of file G4PolarizedPairProductionCrossSection.cc.

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

◆ Initialize()

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

Reimplemented from G4VPolarizedCrossSection.

Definition at line 79 of file G4PolarizedPairProductionCrossSection.cc.

85{
86 G4double aLept1E = aGammaE - aLept0E;
87
88 G4double Stokes_P3 = beamPol.z() ;
89 // **************************************************************************
90
91 G4double m0_c2 = electron_mass_c2;
92 G4double Lept0E = aLept0E/m0_c2+1., Lept0E2 = Lept0E * Lept0E ;
93 G4double GammaE = aGammaE/m0_c2;
94 G4double Lept1E = aLept1E/m0_c2-1., Lept1E2 = Lept1E * Lept1E ;
95
96
97 // const G4Element* theSelectedElement = theModel->SelectedAtom();
98
99 // ******* Gamma Transvers Momentum
100
101 G4double TMom = std::sqrt(Lept0E2 -1.)* sintheta, u = TMom , u2 =u * u ;
102 G4double Xsi = 1./(1.+u2) , Xsi2 = Xsi * Xsi ;
103
104 // G4double theZ = theSelectedElement->GetZ();
105
106 // G4double fCoul = theSelectedElement->GetfCoulomb();
107 G4double delta = 12. * std::pow(theZ, 1./3.) * Lept0E * Lept1E * Xsi / (121. * GammaE);
108 G4double GG=0.;
109
110 if(delta < 0.5) {
111 GG = std::log(2.* Lept0E * Lept1E / GammaE) - 2. - fCoul;
112 }
113 else if ( delta < 120.) {
114 for (G4int j=2; j<=19; j++) {
115 if(SCRN[1][j] >= delta) {
116 GG =std::log(2. * Lept0E * Lept1E / GammaE) - 2. - fCoul
117 -(SCRN[2][j-1]+(delta-SCRN[1][j-1])*(SCRN[2][j]-SCRN[2][j-1])/(SCRN[1][j]-SCRN[1][j-1]));
118 break;
119 }
120 }
121 }
122 else {
123 G4double alpha_sc = (111. * std::pow(theZ, -1./3.)) / Xsi;
124 GG = std::log(alpha_sc)- 2. - fCoul;
125 }
126
127 if(GG<-1.) GG=-1.; // *KL* do we need this ?!
128
129
130 G4double I_Lepton = (Lept0E2 + Lept1E2)*(3+2*GG) + 2. * Lept0E * Lept1E * (1. + 4. * u2 * Xsi2 * GG);
131
132 // G4double D_Lepton1 = -8 * Lept0E * Lept1E * u2 * Xsi2 * GG / I_Lepton;
133
134 G4double L_Lepton1 = GammaE * ((Lept0E - Lept1E) * (3. + 2. * GG)+2 * Lept1E * (1. + 4. * u2 * Xsi2 * GG))/I_Lepton;
135
136 G4double T_Lepton1 = 4. * GammaE * Lept1E * Xsi * u * (1. - 2. * Xsi) * GG / I_Lepton ;
137
138
139 G4double Stokes_S1 = (Stokes_P3 * T_Lepton1) ;
140 G4double Stokes_S2 = 0.;
141 G4double Stokes_S3 = (Stokes_P3 * L_Lepton1) ;
142
143
144 theFinalElectronPolarization.setX(Stokes_S1);
145 theFinalElectronPolarization.setY(Stokes_S2);
146 theFinalElectronPolarization.setZ(Stokes_S3);
147
148 if(theFinalElectronPolarization.mag2()>1.) {
149 G4cout<<" WARNING in pol-conv theFinalElectronPolarization \n";
150 G4cout
151 <<"\t"<<theFinalElectronPolarization
152 <<"\t GG\t"<<GG
153 <<"\t delta\t"<<delta
154 <<G4endl;
155 theFinalElectronPolarization.setX(0.);
156 theFinalElectronPolarization.setY(0.);
157 theFinalElectronPolarization.setZ(Stokes_S3);
158 if(Stokes_S3>1.) theFinalElectronPolarization.setZ(1.);
159 }
160
161
162 G4double L_Lepton2 = GammaE * ((Lept1E - Lept0E) * (3. + 2. * GG)+2 * Lept0E * (1. + 4. * u2 * Xsi2 * GG))/I_Lepton;
163
164 G4double T_Lepton2 = 4. * GammaE * Lept0E * Xsi * u * (1. - 2. * Xsi) * GG / I_Lepton ;
165
166 G4double Stokes_SS1 = (Stokes_P3 * T_Lepton2) ;
167 G4double Stokes_SS2 = 0.;
168 G4double Stokes_SS3 = (Stokes_P3 * L_Lepton2) ;
169
170 theFinalPositronPolarization.SetPhoton();
171
172 theFinalPositronPolarization.setX(Stokes_SS1);
173 theFinalPositronPolarization.setY(Stokes_SS2);
174 theFinalPositronPolarization.setZ(Stokes_SS3);
175
176 if(theFinalPositronPolarization.mag2()>1.) {
177 G4cout<<" WARNING in pol-conv theFinalPositronPolarization \n";
178 G4cout
179 <<"\t"<<theFinalPositronPolarization
180 <<"\t GG\t"<<GG
181 <<"\t delta\t"<<delta
182 <<G4endl;
183 }
184}
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL 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 
)
overridevirtual

Implements G4VPolarizedCrossSection.

Definition at line 186 of file G4PolarizedPairProductionCrossSection.cc.

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

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