G4PolarizedBremsstrahlungCrossSection Class Reference

#include <G4PolarizedBremsstrahlungCrossSection.hh>

Inheritance diagram for G4PolarizedBremsstrahlungCrossSection:

G4VPolarizedCrossSection

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 ()

Detailed Description

Definition at line 54 of file G4PolarizedBremsstrahlungCrossSection.hh.


Constructor & Destructor Documentation

G4PolarizedBremsstrahlungCrossSection::G4PolarizedBremsstrahlungCrossSection (  ) 

Definition at line 75 of file G4PolarizedBremsstrahlungCrossSection.cc.

00076 {
00077   InitializeMe();
00078 }


Member Function Documentation

G4StokesVector G4PolarizedBremsstrahlungCrossSection::GetPol2 (  )  [virtual]

Reimplemented from G4VPolarizedCrossSection.

Definition at line 203 of file G4PolarizedBremsstrahlungCrossSection.cc.

00204 {
00205   // electron/positron
00206   return theFinalLeptonPolarization;
00207 }

G4StokesVector G4PolarizedBremsstrahlungCrossSection::GetPol3 (  )  [virtual]

Reimplemented from G4VPolarizedCrossSection.

Definition at line 208 of file G4PolarizedBremsstrahlungCrossSection.cc.

00209 {
00210   // photon
00211   return theFinalGammaPolarization;;
00212 }

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.

References G4VPolarizedCrossSection::fCoul, G4cout, G4endl, G4StokesVector::SetPhoton(), and G4VPolarizedCrossSection::theZ.

00086 {
00087 //   G4cout<<"G4PolarizedBremsstrahlungCrossSection::Initialize \n"
00088 //      <<"lepE = "<<aLept0E
00089 //      <<"gamE = "<<aGammaE
00090 //      <<"sint = "<<sintheta<<"\n"
00091 //      <<"beamPol="<<beamPol<<"\n";
00092 
00093   G4double aLept1E = aLept0E - aGammaE;
00094 
00095   G4double Stokes_S1  = beamPol.x()   ; 
00096   G4double Stokes_S2  = beamPol.y()   ; 
00097   G4double Stokes_S3  = beamPol.z()   ; 
00098   // **************************************************************************  
00099     
00100   G4double m0_c2  = electron_mass_c2; 
00101   G4double Lept0E = aLept0E/m0_c2+1.,   Lept0E2 = Lept0E * Lept0E ;
00102   G4double GammaE = aGammaE/m0_c2,      GammaE2 = GammaE * GammaE ;
00103   G4double Lept1E = aLept1E/m0_c2+1.,   Lept1E2 = Lept1E * Lept1E ;
00104     
00105 
00106   //  const G4Element* theSelectedElement = theModel->SelectedAtom();
00107 
00108   // *******  Gamma Transvers Momentum
00109     
00110   G4double TMom = std::sqrt(Lept0E2 -1.)* sintheta;
00111   G4double u    = TMom       , u2 =u * u ;
00112   G4double Xsi  = 1./(1.+u2)                      , Xsi2 = Xsi * Xsi  ; 
00113 
00114   //  G4double theZ  = theSelectedElement->GetZ();
00115     
00116   //  G4double fCoul = theSelectedElement->GetfCoulomb();
00117   G4double delta = 12. * std::pow(theZ, 1./3.) * 
00118     Lept0E * Lept1E * Xsi / (121. * GammaE); 
00119   G4double GG=0.;
00120 
00121   if(delta < 0.5) {
00122     GG = std::log(2.* Lept0E * Lept1E / GammaE) - 2. - fCoul; 
00123   }
00124   else if ( delta < 120) {
00125     for (G4int j=2; j<=19; j++)  {
00126       if(SCRN[1][j] >= delta)    {
00127         GG =std::log(2 * Lept0E * Lept1E / GammaE) - 2 - fCoul
00128           -(SCRN[2][j-1]+(delta-SCRN[1][j-1])*(SCRN[2][j]-SCRN[2][j-1])
00129             /(SCRN[1][j]-SCRN[1][j-1]));
00130         break;
00131       }
00132     }
00133   }
00134   else  {
00135     G4double alpha_sc  = (111 * std::pow(theZ, -1./3.)) / Xsi;
00136     GG = std::log(alpha_sc)- 2 - fCoul;
00137   }
00138 
00139   if(GG<-1) GG=-1;     // *KL* do we need this ?!
00140 
00141   G4double I_Lept   = (Lept0E2 + Lept1E2) * (3.+2.*GG) - 2 * Lept0E * Lept1E * (1. + 4. * u2 * Xsi2 * GG);
00142   G4double F_Lept   = Lept1E * 4. * GammaE *  u * Xsi * (1. - 2 * Xsi) * GG / I_Lept;
00143   G4double E_Lept   = Lept0E * 4. * GammaE *  u * Xsi * (2. * Xsi - 1.) * GG / I_Lept; 
00144   G4double M_Lept   = 4. * Lept0E * Lept1E * (1. + GG - 2. * Xsi2 * u2 * GG) / I_Lept ;
00145   G4double P_Lept   = GammaE2 * (1. + 8. * GG * (Xsi - 0.5)*(Xsi - 0.5)) / I_Lept ;
00146     
00147   G4double Stokes_SS1 = M_Lept * Stokes_S1 + E_Lept * Stokes_S3;
00148   G4double Stokes_SS2 = M_Lept * Stokes_S2 ;
00149   G4double Stokes_SS3 = (M_Lept + P_Lept) * Stokes_S3 + F_Lept * Stokes_S1; 
00150 
00151   theFinalLeptonPolarization.setX(Stokes_SS1);
00152   theFinalLeptonPolarization.setY(Stokes_SS2);
00153   theFinalLeptonPolarization.setZ(Stokes_SS3);
00154 
00155   if(theFinalLeptonPolarization.mag2()>1) { 
00156     G4cout<<" WARNING in pol-brem theFinalLeptonPolarization \n";
00157     G4cout
00158       <<"\t"<<theFinalLeptonPolarization
00159       <<"\t GG\t"<<GG
00160       <<"\t delta\t"<<delta
00161       <<G4endl;
00162     theFinalLeptonPolarization.setX(0);
00163     theFinalLeptonPolarization.setY(0);
00164     theFinalLeptonPolarization.setZ(Stokes_SS3);
00165     if(Stokes_SS3>1) theFinalLeptonPolarization.setZ(1);
00166   }
00167 
00168 
00169   G4double I_Gamma   = (Lept0E2 + Lept1E2)*(3+2*GG) - 2 * Lept0E * Lept1E * (1 + 4 * u2 * Xsi2 * GG);
00170   G4double D_Gamma   = 8 * Lept0E * Lept1E * u2 * Xsi2 * GG / I_Gamma;
00171   G4double L_Gamma   = GammaE * ((Lept0E + Lept1E) * (3 + 2 * GG) 
00172                                  - 2 * Lept1E * (1 + 4 * u2 * Xsi2 * GG))/I_Gamma;   
00173   G4double T_Gamma   = 4 * GammaE * Lept1E * Xsi * u * (2 * Xsi - 1) * GG / I_Gamma ;
00174     
00175   G4double Stokes_P1 = D_Gamma ;
00176   G4double Stokes_P2 = 0 ;
00177   G4double Stokes_P3 = (Stokes_S3*L_Gamma + Stokes_S1*T_Gamma) ;
00178 
00179   theFinalGammaPolarization.SetPhoton();
00180 
00181   theFinalGammaPolarization.setX(Stokes_P1);
00182   theFinalGammaPolarization.setY(Stokes_P2);
00183   theFinalGammaPolarization.setZ(Stokes_P3);
00184     
00185   if(theFinalGammaPolarization.mag2()>1) {
00186     G4cout<<" WARNING in pol-brem theFinalGammaPolarization \n";
00187     G4cout
00188       <<"\t"<<theFinalGammaPolarization
00189       <<"\t GG\t"<<GG
00190       <<"\t delta\t"<<delta
00191       <<G4endl;
00192   }
00193 }

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

Implements G4VPolarizedCrossSection.

Definition at line 195 of file G4PolarizedBremsstrahlungCrossSection.cc.

References G4cout.

00197 {
00198   G4cout<<"ERROR dummy routine G4PolarizedBremsstrahlungCrossSection::XSection called \n";
00199   return 0;
00200 }


The documentation for this class was generated from the following files:
Generated on Mon May 27 17:52:56 2013 for Geant4 by  doxygen 1.4.7