G4PolarizedComptonCrossSection.cc

Go to the documentation of this file.
00001 //
00002 // ********************************************************************
00003 // * License and Disclaimer                                           *
00004 // *                                                                  *
00005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
00006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
00007 // * conditions of the Geant4 Software License,  included in the file *
00008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
00009 // * include a list of copyright holders.                             *
00010 // *                                                                  *
00011 // * Neither the authors of this software system, nor their employing *
00012 // * institutes,nor the agencies providing financial support for this *
00013 // * work  make  any representation or  warranty, express or implied, *
00014 // * regarding  this  software system or assume any liability for its *
00015 // * use.  Please see the license in the file  LICENSE  and URL above *
00016 // * for the full disclaimer and the limitation of liability.         *
00017 // *                                                                  *
00018 // * This  code  implementation is the result of  the  scientific and *
00019 // * technical work of the GEANT4 collaboration.                      *
00020 // * By using,  copying,  modifying or  distributing the software (or *
00021 // * any work based  on the software)  you  agree  to acknowledge its *
00022 // * use  in  resulting  scientific  publications,  and indicate your *
00023 // * acceptance of all terms of the Geant4 Software license.          *
00024 // ********************************************************************
00025 //
00026 // $Id: G4PolarizedComptonCrossSection.cc 69847 2013-05-16 09:36:18Z gcosmo $
00027 //
00028 // GEANT4 Class file
00029 //
00030 //
00031 // File name:     G4PolarizedComptonCrossSection
00032 //
00033 // Author:        Andreas Schaelicke
00034 //
00035 // Creation date: 15.05.2005
00036 //
00037 // Modifications:
00038 //
00039 // Class Description:
00040 //   determine the  polarization of the final state 
00041 //   in a Compton scattering process employing the differential 
00042 //   cross section by F.W.Lipps & H.A.Tolhoek
00043 //   ( Physica 20 (1954) 395 )
00044 //   recalculated by P.Starovoitov
00045 //
00046 #include "G4PolarizedComptonCrossSection.hh"
00047 #include "G4PhysicalConstants.hh"
00048 #include "Randomize.hh"
00049 
00050 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00051 G4PolarizedComptonCrossSection::G4PolarizedComptonCrossSection() 
00052   : gammaPol2(false), electronPol3(false)
00053 {
00054   SetYmin(0.);
00055 
00056   //  G4cout<<"G4PolarizedComptonCrossSection() init\n";
00057 
00058   re2 = classic_electr_radius * classic_electr_radius * sqr(4*pi/hbarc);
00059   //  G4double unit_conversion = hbarc_squared ;
00060   //  G4cout<<" (keV)^2* m^2 ="<<unit_conversion<<"\n";
00061   phi0 = 0.; polXS = 0.; unpXS = 0.;
00062   phi2 = G4ThreeVector(0., 0., 0.);
00063   phi3 = G4ThreeVector(0., 0., 0.);
00064   polxx = polyy = polzz = polxz = polzx = polyz = polzy = polxy = polyx = 0.;
00065   diffXSFactor = 1.;
00066   totalXSFactor = 1.;
00067 }
00068 
00069 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00070 G4PolarizedComptonCrossSection::~G4PolarizedComptonCrossSection()
00071 {}
00072 
00073 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00074 void G4PolarizedComptonCrossSection::Initialize(G4double eps, G4double X, G4double , // phi
00075                                                 const G4StokesVector & pol0,
00076                                                 const G4StokesVector & pol1,
00077                                                 G4int flag)
00078 {
00079   G4double cosT = 1. - (1./eps - 1.)/X;
00080   if(cosT > 1.+1.e-8)  cosT = 1.;
00081   if(cosT < -1.-1.e-8) cosT = -1.;
00082   G4double cosT2 = cosT*cosT;
00083   G4double cosT3 = cosT2*cosT;
00084   G4double sinT2 = 1. - cosT2;
00085   if(sinT2 > 1. + 1.e-8) sinT2 = 1.;
00086   if(sinT2 < 0.)     sinT2 = 0.;
00087   G4double sinT = std::sqrt(sinT2);
00088   G4double cos2T = 2.*cosT2 - 1.;
00089   G4double sin2T = 2.*sinT*cosT;
00090   G4double eps2 = sqr(eps);
00091   DefineCoefficients(pol0,pol1);
00092   diffXSFactor = re2/(4.*X);
00093 
00094   // unpolarized Cross Section
00095   unpXS = (eps2 + 1. - eps*sinT2)/(2.*eps);
00096   // initial polarization dependence
00097   polXS = -sinT2*pol0.x() + (1. - eps)*sinT*polzx + ((eps2 - 1.)/eps)*cosT*polzz;
00098   polXS *= 0.5;
00099 
00100   phi0 = unpXS + polXS;
00101 
00102   if (flag == 2 ){
00103   // polarization of outgoing photon
00104   G4double PHI21 = -sinT2 + 0.5*(cos2T + 3.)*pol0.x() - ((1. - eps)/eps)*sinT*polzx;
00105   PHI21 *= 0.5;
00106   G4double PHI22 = cosT*pol0.y() + ((1. - eps)/(2.*eps))*sinT*polzy;
00107   G4double PHI23 = ((eps2 + 1.)/eps)*cosT*pol0.z() - ((1. - eps)/eps)*(eps*cosT2 + 1.)*pol1.z();
00108   PHI23 += 0.5*(1. - eps)*sin2T*pol1.x();
00109   PHI23 += (eps - 1.)*(-sinT2*polxz + sinT*polyy - 0.5*sin2T*polxx);
00110   PHI23 *= 0.5;
00111   phi2 = G4ThreeVector(PHI21, PHI22, PHI23);
00112 
00113   // polarization of outgoing electron
00114   G4double PHI32 = -sinT2*polxy + ((1. - eps)/eps)*sinT*polyz + 0.5*(cos2T + 3.)*pol1.y();
00115   PHI32 *= 0.5;
00116 
00117   G4double PHI31 = 0., PHI31add = 0., PHI33 = 0., PHI33add = 0.;
00118 
00119   if ((1. - eps) > 1.e-12){
00120     G4double helpVar = std::sqrt(eps2 - 2.*cosT*eps + 1.);
00121 
00122     PHI31 = (1. - eps)*(1. + cosT)*sinT*pol0.z();
00123     PHI31 += (-eps*cosT3 + eps*cosT2 + (eps - 2.)*cosT + eps)*pol1.x();
00124     PHI31 += -(eps*cosT2 - eps*cosT + cosT + 1.)*sinT*pol1.z();
00125     PHI31 /= 2.*helpVar;
00126 
00127     PHI31add = -eps*sqr(1. - cosT)*(1. + cosT)*polxx;
00128     PHI31add += (1. - eps)*sinT2*polyy;
00129     PHI31add += -(-eps2 + cosT*(cosT*eps - eps + 1.)*eps + eps - 1.)*sinT*polxz/eps;
00130     PHI31add /= 2.*helpVar;
00131     
00132     PHI33 = ((1. - eps)/eps)*(-eps*cosT2 + eps*(eps + 1.)*cosT - 1.)*pol0.z();
00133     PHI33 += -(eps*cosT2 + (1. - eps)*eps*cosT + 1.)*sinT*pol1.x();
00134     PHI33 += -(-eps2*cosT3 + eps*(eps2 - eps + 1.)*cosT2 - cosT + eps2)*pol1.z()/eps;
00135     PHI33 /= -2.*helpVar;
00136 
00137     PHI33add = (eps*(eps - cosT - 1.)*cosT + 1.)*sinT*polxx;
00138     PHI33add += -(-eps2 + cosT*eps + eps - 1.)*sinT2*polxz;
00139     PHI33add += (eps - 1.)*(cosT - eps)*sinT*polyy;
00140     PHI33add /= -2.*helpVar;
00141   }else{
00142      PHI31 = -pol1.z() - (X - 1.)*std::sqrt(1. - eps)*pol1.x()/std::sqrt(2.*X);
00143      PHI31add = -(-X*X*pol1.z() - 2.*X*(2.*pol0.z() - pol1.z()) - (4.*pol0.x() + 5.)*pol1.z())*(1. - eps)/(4.*X);
00144      
00145      PHI33 = pol1.x() - (X - 1.)*std::sqrt(1. - eps)*pol1.z()/std::sqrt(2.*X);
00146      PHI33add = -(X*X - 2.*X + 4.*pol0.x() + 5.)*(1. - eps)*pol1.x()/(4.*X);
00147   }
00148   phi3 = G4ThreeVector(PHI31 + PHI31add, PHI32, PHI33 + PHI33add);
00149     
00150   }
00151   unpXS *= diffXSFactor;
00152   polXS *= diffXSFactor;
00153   phi0 *= diffXSFactor;
00154   phi2 *= diffXSFactor;
00155   phi3 *= diffXSFactor;
00156   
00157 }
00158 
00159 G4double G4PolarizedComptonCrossSection::XSection(const G4StokesVector & pol2,const G4StokesVector & pol3)
00160 {
00161   gammaPol2    = !(pol2==G4StokesVector::ZERO);
00162   electronPol3 = !(pol3==G4StokesVector::ZERO);
00163 
00164   G4double phi = 0.;
00165   // polarization independent part
00166   phi += phi0;
00167 
00168 
00169   if (gammaPol2) { 
00170     // part depending on the polarization of the final photon  
00171     phi += phi2*pol2;
00172   }
00173 
00174   if (electronPol3) {
00175     // part depending on the polarization of the final electron  
00176     phi += phi3*pol3;
00177   }
00178 
00179   // return cross section.
00180   return phi;
00181 } 
00182 
00183 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00184 
00185 
00186 G4double G4PolarizedComptonCrossSection::TotalXSection(G4double /*xmin*/, G4double /*xmax*/, G4double k0,
00187                                                    const G4StokesVector & pol0,
00188                                                    const G4StokesVector & pol1)
00189 {
00190   
00191   //  G4double k0 = gammaEnergy / electron_mass_c2 ;
00192   G4double k1 = 1 + 2*k0 ;
00193 
00194 //   // pi*re^2
00195 //   G4double re=2.81794e-15; //m
00196 //   G4double barn=1.e-28; //m^2
00197   G4double Z=theZ;
00198 
00199   G4double unit = Z*pi*classic_electr_radius  
00200     * classic_electr_radius ; // *1./barn;
00201 
00202   G4double pre = unit/(sqr(k0)*sqr(1.+2.*k0));
00203 
00204   G4double xs_0 = ((k0 - 2.)*k0  -2.)*sqr(k1)*std::log(k1) + 2.*k0*(k0*(k0 + 1.)*(k0 + 8.) + 2.);               
00205   G4double xs_pol = (k0 + 1.)*sqr(k1)*std::log(k1) - 2.*k0*(5.*sqr(k0) + 4.*k0 + 1.);
00206 
00207   return pre*(xs_0/k0 + pol0.p3()*pol1.z()*xs_pol);
00208 }
00209 
00210 
00211 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00212 
00213 
00214 G4StokesVector G4PolarizedComptonCrossSection::GetPol2()
00215 {
00216   // Note, mean polarization can not contain correlation 
00217   // effects.
00218   return 1./phi0 * phi2;
00219 }
00220 
00221 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00222 
00223 G4StokesVector G4PolarizedComptonCrossSection::GetPol3()
00224 {
00225   // Note, mean polarization can not contain correlation 
00226   // effects.
00227   return 1./phi0 * phi3;
00228 }
00229 
00230 
00231 
00232 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00233 void G4PolarizedComptonCrossSection::DefineCoefficients(const G4StokesVector & pol0,
00234                                                         const G4StokesVector & pol1)
00235 {
00236   polxx=pol0.x()*pol1.x();
00237   polyy=pol0.y()*pol1.y();
00238   polzz=pol0.z()*pol1.z();
00239 
00240   polxz=pol0.x()*pol1.z();
00241   polzx=pol0.z()*pol1.x();
00242 
00243   polyz=pol0.y()*pol1.z();
00244   polzy=pol0.z()*pol1.y();
00245 
00246   polxy=pol0.x()*pol1.y();
00247   polyx=pol0.y()*pol1.x();
00248 }
00249 
00250 
00251 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......

Generated on Mon May 27 17:49:23 2013 for Geant4 by  doxygen 1.4.7