Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions
G4PolarizedPEEffectCrossSection Class Reference

#include <G4PolarizedPEEffectCrossSection.hh>

Inheritance diagram for G4PolarizedPEEffectCrossSection:
G4VPolarizedCrossSection

Public Member Functions

 G4PolarizedPEEffectCrossSection ()
 
virtual ~G4PolarizedPEEffectCrossSection ()
 
virtual void Initialize (G4double aGammaE, G4double aLept0E, G4double sintheta, const G4StokesVector &beamPol, const G4StokesVector &, G4int flag=0)
 
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)
 

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 50 of file G4PolarizedPEEffectCrossSection.hh.

Constructor & Destructor Documentation

G4PolarizedPEEffectCrossSection::G4PolarizedPEEffectCrossSection ( )

Definition at line 48 of file G4PolarizedPEEffectCrossSection.cc.

49  {
50  cout<<"G4PolarizedPEEffectCrossSection() init\n";
51 
52 }
G4PolarizedPEEffectCrossSection::~G4PolarizedPEEffectCrossSection ( )
virtual

Definition at line 55 of file G4PolarizedPEEffectCrossSection.cc.

56 {}

Member Function Documentation

G4StokesVector G4PolarizedPEEffectCrossSection::GetPol2 ( )
virtual

Reimplemented from G4VPolarizedCrossSection.

Definition at line 150 of file G4PolarizedPEEffectCrossSection.cc.

Referenced by G4PolarizedPEEffectModel::SampleSecondaries().

151 {
152  return theFinalElectronPolarization;
153 }
G4StokesVector G4PolarizedPEEffectCrossSection::GetPol3 ( )
virtual

Reimplemented from G4VPolarizedCrossSection.

Definition at line 157 of file G4PolarizedPEEffectCrossSection.cc.

158 {
159  return G4StokesVector();
160 }
void G4PolarizedPEEffectCrossSection::Initialize ( G4double  aGammaE,
G4double  aLept0E,
G4double  sintheta,
const G4StokesVector beamPol,
const G4StokesVector ,
G4int  flag = 0 
)
virtual

Reimplemented from G4VPolarizedCrossSection.

Definition at line 59 of file G4PolarizedPEEffectCrossSection.cc.

References python.hepunit::electron_mass_c2, and CLHEP::Hep3Vector::z().

Referenced by G4PolarizedPEEffectModel::SampleSecondaries().

65 {
66  // cout<<"G4PolarizedPEEffectCrossSection::Initialize()\n";
67 
68 // G4StokesVector PolarizedPhotoElectricEffect::Transfer_G4StokesVector(
69 // G4double aGammaE, // Incoming Primary Gamma Energy.
70 // G4ThreeVector aGammaDir, // Incoming Primary Gamma Direction.
71 // G4StokesVector beamPol, // Incoming Primary Gamma polarization.
72 // G4double aLept0E, // The Lepton e- of interest Total energy.
73 // G4ThreeVector aParticl_01_Dir, // The Lepton e- of interest direction.
74 // G4double cos_aTetha_Angle // The lepton of interest Scattering angle.
75 // )
76 
77 // ***********************************************************
78 // ************ added by Karim Polarization transfer to e- in PhotoelectricEffect.
79 // ************
80 // ***********************************************************
81  G4double Gfactor = aLept0E/electron_mass_c2+1;
82  G4double Gfactor_2 = Gfactor * Gfactor;
83 
84  G4double BETA = sqrt(1. - 1./(Gfactor_2));
85 
86  G4double Stokes_P3 = beamPol.z() ;
87 
88  G4double m0_c2 = electron_mass_c2;
89  G4double Lept0E = aLept0E/m0_c2+1, Lept0E2 = Lept0E * Lept0E ;
90  G4double GammaE = aGammaE/m0_c2;
91 
92 
93 // G4double cosTheta = cos_aTetha_Angle;
94 // G4double sinTheta = sqrt(1- cos_aTetha_Angle * cos_aTetha_Angle);
95  G4double cosTheta = std::sqrt(1. - sinTheta*sinTheta);
96 
97  G4double D_Lepton0 = (1./GammaE) * ((2./(GammaE*Lept0E*(1-BETA*cosTheta)))-1.);
98 
99  G4double I_Lepton0 = 1.0+D_Lepton0;
100 
101  G4double A_Lepton0 = (Lept0E/(Lept0E+1))*(2.0/(GammaE*Lept0E)
102  + BETA*cosTheta
103  +(2.0/((GammaE*Lept0E2)*(1.0-BETA*cosTheta)))) / I_Lepton0 ;
104 
105  G4double B_Lepton0 = (Lept0E/(Lept0E+1.0)) * BETA * sinTheta * (2.0/(GammaE*Lept0E*(1-BETA*cosTheta))-1.0)/I_Lepton0;
106 
107  G4double Stokes_S1 = (Stokes_P3 * B_Lepton0) ;
108  G4double Stokes_S2 = 0.;
109  G4double Stokes_S3 = (Stokes_P3 * A_Lepton0) ;
110 
111 
112  theFinalElectronPolarization.setX(Stokes_S1);
113  theFinalElectronPolarization.setY(Stokes_S2);
114  theFinalElectronPolarization.setZ(Stokes_S3);
115 
116  if((theFinalElectronPolarization.x()*theFinalElectronPolarization.x()
117  + theFinalElectronPolarization.y()* theFinalElectronPolarization.y()
118  + theFinalElectronPolarization.z()* theFinalElectronPolarization.z())>1)
119 
120  {
121  cout<<"Warning: PhotoelectricEffect Problem in pol-transfer photon to lepton:Px2 + Py2 + Pz2 > 1"<<endl;
122  cout<<"Polarization transfer forced to be total and similar as incoming Photo"<<endl;
123  // *KL* Surprising if it arrives (never seen it up to now)
124  theFinalElectronPolarization = beamPol; // suplement de securite
125 // cout<<"PhotoEffect okay :"
126 // <<"\t"<<(aLept0E-m0_c2)/aGammaE
127 // <<"\t"<<aGammaE
128 // <<"\t"<<aLept0E
129 // <<"\t"<<cos_aTetha_Angle
130 // <<"\t"<<beamPol
131 // <<"\t"<<theFinalElectronPolarization
132 // <<"\t"<<A_Lepton0
133 // <<endl;
134  }
135 }
double x() const
void setY(double)
double z() const
void setZ(double)
void setX(double)
float electron_mass_c2
Definition: hepunit.py:274
double y() const
double G4double
Definition: G4Types.hh:76
G4double G4PolarizedPEEffectCrossSection::XSection ( const G4StokesVector pol2,
const G4StokesVector pol3 
)
virtual

Implements G4VPolarizedCrossSection.

Definition at line 141 of file G4PolarizedPEEffectCrossSection.cc.

143 {
144  cout<<"ERROR dummy routine G4PolarizedPEEffectCrossSection::XSection() called\n";
145  return 0;
146 }

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