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

#include <G4PolarizedBremsstrahlungCrossSection.hh>

+ Inheritance diagram for G4PolarizedBremsstrahlungCrossSection:

Public Member Functions

 G4PolarizedBremsstrahlungCrossSection ()
 
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 53 of file G4PolarizedBremsstrahlungCrossSection.hh.

Constructor & Destructor Documentation

◆ G4PolarizedBremsstrahlungCrossSection()

G4PolarizedBremsstrahlungCrossSection::G4PolarizedBremsstrahlungCrossSection ( )

Definition at line 74 of file G4PolarizedBremsstrahlungCrossSection.cc.

75{
76 InitializeMe();
77}

Member Function Documentation

◆ GetPol2()

G4StokesVector G4PolarizedBremsstrahlungCrossSection::GetPol2 ( )
overridevirtual

Reimplemented from G4VPolarizedCrossSection.

Definition at line 202 of file G4PolarizedBremsstrahlungCrossSection.cc.

203{
204 // electron/positron
205 return theFinalLeptonPolarization;
206}

◆ GetPol3()

G4StokesVector G4PolarizedBremsstrahlungCrossSection::GetPol3 ( )
overridevirtual

Reimplemented from G4VPolarizedCrossSection.

Definition at line 207 of file G4PolarizedBremsstrahlungCrossSection.cc.

208{
209 // photon
210 return theFinalGammaPolarization;;
211}

◆ Initialize()

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

Reimplemented from G4VPolarizedCrossSection.

Definition at line 80 of file G4PolarizedBremsstrahlungCrossSection.cc.

85{
86// G4cout<<"G4PolarizedBremsstrahlungCrossSection::Initialize \n"
87// <<"lepE = "<<aLept0E
88// <<"gamE = "<<aGammaE
89// <<"sint = "<<sintheta<<"\n"
90// <<"beamPol="<<beamPol<<"\n";
91
92 G4double aLept1E = aLept0E - aGammaE;
93
94 G4double Stokes_S1 = beamPol.x() ;
95 G4double Stokes_S2 = beamPol.y() ;
96 G4double Stokes_S3 = beamPol.z() ;
97 // **************************************************************************
98
99 G4double m0_c2 = electron_mass_c2;
100 G4double Lept0E = aLept0E/m0_c2+1., Lept0E2 = Lept0E * Lept0E ;
101 G4double GammaE = aGammaE/m0_c2, GammaE2 = GammaE * GammaE ;
102 G4double Lept1E = aLept1E/m0_c2+1., Lept1E2 = Lept1E * Lept1E ;
103
104
105 // const G4Element* theSelectedElement = theModel->SelectedAtom();
106
107 // ******* Gamma Transvers Momentum
108
109 G4double TMom = std::sqrt(Lept0E2 -1.)* sintheta;
110 G4double u = TMom , u2 =u * u ;
111 G4double Xsi = 1./(1.+u2) , Xsi2 = Xsi * Xsi ;
112
113 // G4double theZ = theSelectedElement->GetZ();
114
115 // G4double fCoul = theSelectedElement->GetfCoulomb();
116 G4double delta = 12. * std::pow(theZ, 1./3.) *
117 Lept0E * Lept1E * Xsi / (121. * GammaE);
118 G4double GG=0.;
119
120 if(delta < 0.5) {
121 GG = std::log(2.* Lept0E * Lept1E / GammaE) - 2. - fCoul;
122 }
123 else if ( delta < 120) {
124 for (G4int j=2; j<=19; j++) {
125 if(SCRN[1][j] >= delta) {
126 GG =std::log(2 * Lept0E * Lept1E / GammaE) - 2 - fCoul
127 -(SCRN[2][j-1]+(delta-SCRN[1][j-1])*(SCRN[2][j]-SCRN[2][j-1])
128 /(SCRN[1][j]-SCRN[1][j-1]));
129 break;
130 }
131 }
132 }
133 else {
134 G4double alpha_sc = (111. * std::pow(theZ, -1./3.)) / Xsi;
135 GG = std::log(alpha_sc)- 2. - fCoul;
136 }
137
138 if(GG<-1.) GG=-1.; // *KL* do we need this ?!
139
140 G4double I_Lept = (Lept0E2 + Lept1E2) * (3.+2.*GG) - 2 * Lept0E * Lept1E * (1. + 4. * u2 * Xsi2 * GG);
141 G4double F_Lept = Lept1E * 4. * GammaE * u * Xsi * (1. - 2 * Xsi) * GG / I_Lept;
142 G4double E_Lept = Lept0E * 4. * GammaE * u * Xsi * (2. * Xsi - 1.) * GG / I_Lept;
143 G4double M_Lept = 4. * Lept0E * Lept1E * (1. + GG - 2. * Xsi2 * u2 * GG) / I_Lept ;
144 G4double P_Lept = GammaE2 * (1. + 8. * GG * (Xsi - 0.5)*(Xsi - 0.5)) / I_Lept ;
145
146 G4double Stokes_SS1 = M_Lept * Stokes_S1 + E_Lept * Stokes_S3;
147 G4double Stokes_SS2 = M_Lept * Stokes_S2 ;
148 G4double Stokes_SS3 = (M_Lept + P_Lept) * Stokes_S3 + F_Lept * Stokes_S1;
149
150 theFinalLeptonPolarization.setX(Stokes_SS1);
151 theFinalLeptonPolarization.setY(Stokes_SS2);
152 theFinalLeptonPolarization.setZ(Stokes_SS3);
153
154 if(theFinalLeptonPolarization.mag2()>1) {
155 G4cout<<" WARNING in pol-brem theFinalLeptonPolarization \n";
156 G4cout
157 <<"\t"<<theFinalLeptonPolarization
158 <<"\t GG\t"<<GG
159 <<"\t delta\t"<<delta
160 <<G4endl;
161 theFinalLeptonPolarization.setX(0);
162 theFinalLeptonPolarization.setY(0);
163 theFinalLeptonPolarization.setZ(Stokes_SS3);
164 if(Stokes_SS3>1) theFinalLeptonPolarization.setZ(1);
165 }
166
167
168 G4double I_Gamma = (Lept0E2 + Lept1E2)*(3.+2.*GG) - 2. * Lept0E * Lept1E * (1. + 4. * u2 * Xsi2 * GG);
169 G4double D_Gamma = 8. * Lept0E * Lept1E * u2 * Xsi2 * GG / I_Gamma;
170 G4double L_Gamma = GammaE * ((Lept0E + Lept1E) * (3. + 2. * GG)
171 - 2. * Lept1E * (1. + 4. * u2 * Xsi2 * GG))/I_Gamma;
172 G4double T_Gamma = 4. * GammaE * Lept1E * Xsi * u * (2. * Xsi - 1.) * GG / I_Gamma ;
173
174 G4double Stokes_P1 = D_Gamma ;
175 G4double Stokes_P2 = 0.;
176 G4double Stokes_P3 = (Stokes_S3*L_Gamma + Stokes_S1*T_Gamma) ;
177
178 theFinalGammaPolarization.SetPhoton();
179
180 theFinalGammaPolarization.setX(Stokes_P1);
181 theFinalGammaPolarization.setY(Stokes_P2);
182 theFinalGammaPolarization.setZ(Stokes_P3);
183
184 if(theFinalGammaPolarization.mag2()>1) {
185 G4cout<<" WARNING in pol-brem theFinalGammaPolarization \n";
186 G4cout
187 <<"\t"<<theFinalGammaPolarization
188 <<"\t GG\t"<<GG
189 <<"\t delta\t"<<delta
190 <<G4endl;
191 }
192}
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 G4PolarizedBremsstrahlungCrossSection::XSection ( const G4StokesVector pol2,
const G4StokesVector pol3 
)
overridevirtual

Implements G4VPolarizedCrossSection.

Definition at line 194 of file G4PolarizedBremsstrahlungCrossSection.cc.

196{
197 G4cout<<"ERROR dummy routine G4PolarizedBremsstrahlungCrossSection::XSection called \n";
198 return 0.;
199}

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