Geant4  10.02.p03
G4PolarizedBremsstrahlungCrossSection Class Reference

#include <G4PolarizedBremsstrahlungCrossSection.hh>

Inheritance diagram for G4PolarizedBremsstrahlungCrossSection:
Collaboration 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)
 
virtual G4double XSection (const G4StokesVector &pol2, const G4StokesVector &pol3)
 
G4StokesVector GetPol2 ()
 
G4StokesVector GetPol3 ()
 
- Public Member Functions inherited from G4VPolarizedCrossSection
 G4VPolarizedCrossSection ()
 
virtual ~G4VPolarizedCrossSection ()
 
virtual G4double TotalXSection (G4double xmin, G4double xmax, G4double y, const G4StokesVector &pol0, const G4StokesVector &pol1)
 
G4double GetYmin ()
 
virtual G4double GetXmin (G4double y)
 
virtual G4double GetXmax (G4double y)
 
void SetMaterial (G4double A, G4double Z, G4double coul)
 

Private Member Functions

void InitializeMe ()
 

Private Attributes

G4StokesVector theFinalLeptonPolarization
 
G4StokesVector theFinalGammaPolarization
 

Static Private Attributes

static G4bool scrnInitialized =false
 
static G4double SCRN [3][20]
 

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 54 of file G4PolarizedBremsstrahlungCrossSection.hh.

Constructor & Destructor Documentation

◆ G4PolarizedBremsstrahlungCrossSection()

G4PolarizedBremsstrahlungCrossSection::G4PolarizedBremsstrahlungCrossSection ( )

Definition at line 75 of file G4PolarizedBremsstrahlungCrossSection.cc.

Here is the call graph for this function:

Member Function Documentation

◆ GetPol2()

G4StokesVector G4PolarizedBremsstrahlungCrossSection::GetPol2 ( )
virtual

Reimplemented from G4VPolarizedCrossSection.

Definition at line 203 of file G4PolarizedBremsstrahlungCrossSection.cc.

204 {
205  // electron/positron
207 }

◆ GetPol3()

G4StokesVector G4PolarizedBremsstrahlungCrossSection::GetPol3 ( )
virtual

Reimplemented from G4VPolarizedCrossSection.

Definition at line 208 of file G4PolarizedBremsstrahlungCrossSection.cc.

209 {
210  // photon
212 }

◆ Initialize()

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

Reimplemented from G4VPolarizedCrossSection.

Definition at line 81 of file G4PolarizedBremsstrahlungCrossSection.cc.

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

◆ InitializeMe()

void G4PolarizedBremsstrahlungCrossSection::InitializeMe ( )
private

Definition at line 48 of file G4PolarizedBremsstrahlungCrossSection.cc.

49 {
50  if (!scrnInitialized) {
51  SCRN [1][1]= 0.5 ; SCRN [2][1] = 0.0145;
52  SCRN [1][2]= 1.0 ; SCRN [2][2] = 0.0490;
53  SCRN [1][3]= 2.0 ; SCRN [2][3] = 0.1400;
54  SCRN [1][4]= 4.0 ; SCRN [2][4] = 0.3312;
55  SCRN [1][5]= 8.0 ; SCRN [2][5] = 0.6758;
56  SCRN [1][6]= 15.0 ; SCRN [2][6] = 1.126;
57  SCRN [1][7]= 20.0 ; SCRN [2][7] = 1.367;
58  SCRN [1][8]= 25.0 ; SCRN [2][8] = 1.564;
59  SCRN [1][9]= 30.0 ; SCRN [2][9] = 1.731;
60  SCRN [1][10]= 35.0 ; SCRN [2][10]= 1.875;
61  SCRN [1][11]= 40.0 ; SCRN [2][11]= 2.001;
62  SCRN [1][12]= 45.0 ; SCRN [2][12]= 2.114;
63  SCRN [1][13]= 50.0 ; SCRN [2][13]= 2.216;
64  SCRN [1][14]= 60.0 ; SCRN [2][14]= 2.393;
65  SCRN [1][15]= 70.0 ; SCRN [2][15]= 2.545;
66  SCRN [1][16]= 80.0 ; SCRN [2][16]= 2.676;
67  SCRN [1][17]= 90.0 ; SCRN [2][17]= 2.793;
68  SCRN [1][18]= 100.0 ; SCRN [2][18]= 2.897;
69  SCRN [1][19]= 120.0 ; SCRN [2][19]= 3.078;
70 
71  scrnInitialized=true;
72  }
73 }
Here is the caller graph for this function:

◆ XSection()

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

Implements G4VPolarizedCrossSection.

Definition at line 195 of file G4PolarizedBremsstrahlungCrossSection.cc.

197 {
198  G4cout<<"ERROR dummy routine G4PolarizedBremsstrahlungCrossSection::XSection called \n";
199  return 0;
200 }
G4GLOB_DLL std::ostream G4cout

Member Data Documentation

◆ SCRN

G4double G4PolarizedBremsstrahlungCrossSection::SCRN
staticprivate

Definition at line 76 of file G4PolarizedBremsstrahlungCrossSection.hh.

◆ scrnInitialized

G4bool G4PolarizedBremsstrahlungCrossSection::scrnInitialized =false
staticprivate

Definition at line 75 of file G4PolarizedBremsstrahlungCrossSection.hh.

◆ theFinalGammaPolarization

G4StokesVector G4PolarizedBremsstrahlungCrossSection::theFinalGammaPolarization
private

Definition at line 71 of file G4PolarizedBremsstrahlungCrossSection.hh.

◆ theFinalLeptonPolarization

G4StokesVector G4PolarizedBremsstrahlungCrossSection::theFinalLeptonPolarization
private

Definition at line 70 of file G4PolarizedBremsstrahlungCrossSection.hh.


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