G4PolarizedBhabhaCrossSection.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: G4PolarizedBhabhaCrossSection.cc 69847 2013-05-16 09:36:18Z gcosmo $
00027 // -------------------------------------------------------------------
00028 //
00029 // GEANT4 Class file
00030 //
00031 //
00032 // File name:     G4PolarizedBhabhaCrossSection
00033 //
00034 // Author:        Andreas Schaelicke
00035 //
00036 // Creation date: 12.01.2006
00037 //
00038 // Modifications:
00039 //   16-01-06 included cross section as calculated by P.Starovoitov
00040 //   24-08-06 bugfix in total cross section (A. Schaelicke)
00041 //   07-11-06 modify reference system for polarisation vectors 
00042 //            (A. Schaelicke & P.Starovoitov)
00043 //
00044 // Class Description:
00045 //   * calculates the differential cross section
00046 //     incomming positron Kpl(along positive z direction) scatters at 
00047 //     an electron Kmn at rest
00048 //   * phi denotes the angle between the scattering plane (defined by the
00049 //     outgoing electron) and X-axis
00050 //   * all stokes vectors refer to spins in the Global System (X,Y,Z)
00051 //
00052 
00053 #include "G4PolarizedBhabhaCrossSection.hh"
00054 #include "G4PhysicalConstants.hh"
00055 
00056 G4PolarizedBhabhaCrossSection::G4PolarizedBhabhaCrossSection() : phi0(1.)
00057 {
00058 }
00059 G4PolarizedBhabhaCrossSection::~G4PolarizedBhabhaCrossSection()
00060 {
00061 }
00062 void G4PolarizedBhabhaCrossSection::Initialize(
00063                           G4double e,
00064                           G4double gamma,
00065                           G4double /*phi*/,
00066                     const G4StokesVector & pol0,
00067                     const G4StokesVector & pol1,
00068                           G4int flag)
00069 {
00070   SetXmax(1.);
00071 
00072   G4double re2 = classic_electr_radius * classic_electr_radius;
00073   G4double gamma2 = gamma*gamma;
00074   G4double gamma3 = gamma2*gamma;
00075   G4double gmo = (gamma - 1.);
00076   G4double gmo2 = (gamma - 1.)*(gamma - 1.);
00077   G4double gmo3 = gmo2*(gamma - 1.);
00078   G4double gpo = (gamma + 1.);
00079   G4double gpo2 = (gamma + 1.)*(gamma + 1.);
00080   G4double gpo3 = gpo2*(gamma + 1.);
00081   G4double gpo12 = std::sqrt(gpo);
00082   G4double gpo32 = gpo*gpo12;
00083   G4double gpo52 = gpo2*gpo12;
00084 
00085   G4double pref = re2/(gamma - 1.0);
00086   G4double sqrttwo=std::sqrt(2.);
00087   G4double d = std::sqrt(1./e - 1.);
00088   G4double e2 = e*e;
00089   G4double e3 = e2*e;
00090 
00091   // *** new ***
00092   G4double gmo12 = std::sqrt(gmo);
00093   G4double gmo32 = gmo*gmo12;
00094   G4double egmp32 = std::pow(e*(2 + e*gmo)*gpo,(3./2.));
00095   G4double e32 = e*std::sqrt(e);
00096 
00097   G4bool polarized=(!pol0.IsZero())||(!pol1.IsZero());
00098 
00099   if (flag==0) polarized=false;
00100   // Unpolarised part of XS
00101   // *AS* UnpME . OK 
00102   phi0 = 0.;
00103   phi0+= e2*gmo3/gpo3;
00104   phi0+= -2.*e*gamma*gmo2/gpo3;
00105   phi0+= (3.*gamma2 + 6.*gamma + 4.)*gmo/gpo3;
00106   phi0+= -(2.*gamma2 + 4.*gamma + 1.)/(e*gpo2);
00107   phi0+= gamma2/(e2*(gamma2 - 1.));
00108   phi0*=0.25;
00109   // Initial state polarisarion dependence
00110   if (polarized) {
00111     //    G4cout<<"Polarized differential Bhabha cross section"<<G4endl;
00112     //    G4cout<<"Initial state polarisation contributions"<<G4endl;
00113     //    G4cout<<"Diagonal Matrix Elements"<<G4endl;
00114     // *** new ***
00115     G4double xx = -((e*gmo - gamma)*(-1 - gamma + e*(e*gmo - gamma)*(3 + gamma)))/(4*e*gpo3);
00116     G4double yy = (e3*gmo3 - 2*e2*gmo2*gamma - gpo*(1 + 2*gamma) + e*(-2 + gamma2 + gamma3))/(4*e*gpo3);
00117     G4double zz = ((e*gmo - gamma)*(e2*gmo*(3 + gamma) - e*gamma*(3 + gamma) + gpo*(1 + 2*gamma)))/(4*e*gpo3);
00118     // ***
00119 
00120     phi0 += xx*pol0.x()*pol1.x() + yy*pol0.y()*pol1.y() + zz*pol0.z()*pol1.z();
00121 
00122     {
00123       G4double xy = 0;
00124       G4double xz = (d*(e*gmo - gamma)*(-1 + 2*e*gmo - gamma))/(2*sqrttwo*gpo52);
00125       G4double yx = 0;
00126       G4double yz = 0;
00127       G4double zx = xz;
00128       G4double zy = 0;
00129     //    G4cout<<"Non-diagonal Matrix Elements"<<G4endl;
00130       phi0+=yx*pol0.y()*pol1.x() + xy*pol0.x()*pol1.y();
00131       phi0+=zx*pol0.z()*pol1.x() + xz*pol0.x()*pol1.z();
00132       phi0+=zy*pol0.z()*pol1.y() + yz*pol0.y()*pol1.z();
00133     }
00134   }
00135   // Final state polarisarion dependence
00136   phi2=G4ThreeVector();
00137   phi3=G4ThreeVector();
00138 
00139   if (flag>=1) {
00140     //
00141     // Final Positron Ppl
00142     //
00143         // initial positron Kpl
00144     if (!pol0.IsZero()) {
00145 
00146       G4double xxPplKpl = -((-1 + e)*(e*gmo - gamma)*(-(gamma*gpo) +  e*(-2 + gamma + gamma2)))/
00147         (4*e2*gpo*std::sqrt(gmo*gpo*(-1 + e + gamma - e*gamma)*  (1 + e + gamma - e*gamma)));
00148       G4double xyPplKpl = 0;
00149       G4double xzPplKpl = ((e*gmo - gamma)*(-1 - gamma + e*gmo*(1 + 2*gamma)))/
00150         (2*sqrttwo*e32*gmo*gpo2*std::sqrt(1 + e + gamma - e*gamma));
00151       G4double yxPplKpl = 0;
00152       G4double yyPplKpl = (gamma2*gpo + e2*gmo2*(3 + gamma) - 
00153                   e*gmo*(1 + 2*gamma*(2 + gamma)))/(4*e2*gmo*gpo2);
00154       G4double yzPplKpl = 0;
00155       G4double zxPplKpl = ((e*gmo - gamma)*(1 + e*(-1 + 2*e*gmo - 2*gamma)*gmo + gamma))/
00156         (2*sqrttwo*e*gmo*gpo2*std::sqrt(e*(1 + e + gamma - e*gamma)));
00157       G4double zyPplKpl = 0;
00158       G4double zzPplKpl = -((e*gmo - gamma)*std::sqrt((1 - e)/(e - e*gamma2 + gpo2))*
00159                    (2*e2*gmo2 + gamma + gamma2 - e*(-2 + gamma + gamma2)))/
00160         (4*e2*(-1 + gamma2));
00161 
00162       phi2[0] += xxPplKpl*pol0.x() + xyPplKpl*pol0.y() + xzPplKpl*pol0.z();
00163       phi2[1] += yxPplKpl*pol0.x() + yyPplKpl*pol0.y() + yzPplKpl*pol0.z();
00164       phi2[2] += zxPplKpl*pol0.x() + zyPplKpl*pol0.y() + zzPplKpl*pol0.z();
00165     }
00166         // initial electron Kmn
00167     if (!pol1.IsZero()) {
00168       G4double xxPplKmn = ((-1 + e)*(e*(-2 + gamma)*gmo + gamma))/(4*e*gpo32*std::sqrt(1 + e2*gmo + gamma - 2*e*gamma));
00169       G4double xyPplKmn = 0;
00170       G4double xzPplKmn = (-1 + e*gmo + gmo*gamma)/(2*sqrttwo*gpo2* std::sqrt(e*(1 + e + gamma - e*gamma)));
00171       G4double yxPplKmn = 0;
00172       G4double yyPplKmn = (-1 - 2*gamma + e*gmo*(3 + gamma))/(4*e*gpo2);
00173       G4double yzPplKmn = 0;
00174       G4double zxPplKmn = (1 + 2*e2*gmo2 + gamma + gamma2 + e*(1 + (3 - 4*gamma)*gamma))/
00175         (2*sqrttwo*gpo2*std::sqrt(e*(1 + e + gamma - e*gamma)));
00176       G4double zyPplKmn = 0;
00177       G4double zzPplKmn = -(std::sqrt((1 - e)/(e - e*gamma2 + gpo2))*
00178                    (2*e2*gmo2 + gamma + 2*gamma2 + e*(2 + gamma - 3*gamma2)))/(4*e*gpo);
00179 
00180       phi2[0] += xxPplKmn*pol1.x() + xyPplKmn*pol1.y() + xzPplKmn*pol1.z();
00181       phi2[1] += yxPplKmn*pol1.x() + yyPplKmn*pol1.y() + yzPplKmn*pol1.z();
00182       phi2[2] += zxPplKmn*pol1.x() + zyPplKmn*pol1.y() + zzPplKmn*pol1.z();
00183     }
00184 //
00185 // Final Electron Pmn
00186 //
00187         // initial positron Kpl
00188     if (!pol0.IsZero()) {
00189       G4double xxPmnKpl = ((-1 + e*gmo)*(2 + gamma))/(4*gpo* std::sqrt(e*(2 + e*gmo)*gpo));
00190       G4double xyPmnKpl = 0;
00191       G4double xzPmnKpl = (std::sqrt((-1 + e)/(-2 + e - e*gamma))*
00192             (e + gamma + e*gamma -   2*(-1 + e)*gamma2))/(2*sqrttwo*e*gpo2);
00193       G4double yxPmnKpl = 0;
00194       G4double yyPmnKpl = (-1 - 2*gamma + e*gmo*(3 + gamma))/(4*e*gpo2);
00195       G4double yzPmnKpl = 0;
00196       G4double zxPmnKpl = -((-1 + e)*(1 + 2*e*gmo)*(e*gmo - gamma))/
00197         (2*sqrttwo*e*std::sqrt(-((-1 + e)*(2 + e*gmo)))*gpo2);
00198       G4double zyPmnKpl = 0;
00199       G4double zzPmnKpl = (-2 + 2*e2*gmo2 + gamma*(-1 + 2*gamma) + 
00200                   e*(-2 + (5 - 3*gamma)*gamma))/(4*std::sqrt(e*(2 + e*gmo))*  gpo32);
00201 
00202       phi3[0] += xxPmnKpl*pol0.x() + xyPmnKpl*pol0.y() + xzPmnKpl*pol0.z();
00203       phi3[1] += yxPmnKpl*pol0.x() + yyPmnKpl*pol0.y() + yzPmnKpl*pol0.z();
00204       phi3[2] += zxPmnKpl*pol0.x() + zyPmnKpl*pol0.y() + zzPmnKpl*pol0.z();
00205     }
00206         // initial electron Kmn
00207     if (!pol1.IsZero()) {
00208       G4double xxPmnKmn = -((2 + e*gmo)*(-1 + e*gmo - gamma)*(e*gmo - gamma)*
00209                    (-2 + gamma))/(4*gmo*egmp32);
00210       G4double xyPmnKmn = 0;
00211       G4double xzPmnKmn = ((e*gmo - gamma)*
00212                            std::sqrt((-1 + e + gamma - e*gamma)/(2 + e*gmo))*
00213                            (e + gamma - e*gamma + gamma2))/
00214         (2*sqrttwo*e2*gmo32*gpo2);
00215       G4double yxPmnKmn = 0;
00216       G4double yyPmnKmn = (gamma2*gpo + e2*gmo2*(3 + gamma) - 
00217                   e*gmo*(1 + 2*gamma*(2 + gamma)))/(4*e2*gmo*gpo2);
00218       G4double yzPmnKmn = 0;
00219       G4double zxPmnKmn = -((-1 + e)*(e*gmo - gamma)*(e*gmo + 2*e2*gmo2 - gamma*gpo))/
00220         (2*sqrttwo*e2*std::sqrt(-((-1 + e)*(2 + e*gmo)))* gmo*gpo2);
00221       G4double zyPmnKmn = 0;
00222       G4double zzPmnKmn = ((e*gmo - gamma)*std::sqrt(e/((2 + e*gmo)*gpo))*
00223                   (-(e*(-2 + gamma)*gmo) + 2*e2*gmo2 + (-2 + gamma)*gpo))/(4*e2*(-1 + gamma2));
00224 
00225       phi3[0] += xxPmnKmn*pol1.x() + xyPmnKmn*pol1.y() + xzPmnKmn*pol1.z();
00226       phi3[1] += yxPmnKmn*pol1.x() + yyPmnKmn*pol1.y() + yzPmnKmn*pol1.z();
00227       phi3[2] += zxPmnKmn*pol1.x() + zyPmnKmn*pol1.y() + zzPmnKmn*pol1.z();
00228     }
00229   }
00230   phi0 *= pref;
00231   phi2 *= pref;
00232   phi3 *= pref;
00233 
00234 }
00235 
00236 G4double G4PolarizedBhabhaCrossSection::XSection(const G4StokesVector & pol2,
00237                                                  const G4StokesVector & pol3)
00238 {
00239   G4double xs=0.;
00240   xs+=phi0;
00241 
00242   G4bool polarized=(!pol2.IsZero())||(!pol3.IsZero());
00243   if (polarized) {
00244     xs+=phi2*pol2 + phi3*pol3;
00245   }
00246   return xs;
00247 }
00248 
00249 G4double G4PolarizedBhabhaCrossSection::TotalXSection(
00250   G4double xmin, G4double xmax, G4double gamma,
00251   const G4StokesVector & pol0,const G4StokesVector & pol1)
00252 {
00253   G4double xs=0.;
00254 
00255   G4double x=xmin;
00256 
00257   if (xmax != 1.) G4cout<<" warning xmax expected to be 1 but is "<<xmax<< G4endl;
00258 
00259   //  re -> electron radius^2;
00260   G4double re2 = classic_electr_radius * classic_electr_radius;
00261   G4double gamma2=gamma*gamma;
00262   G4double gmo2 = (gamma - 1.)*(gamma - 1.);
00263   G4double gpo2 = (gamma + 1.)*(gamma + 1.);
00264   G4double gpo3 = gpo2*(gamma + 1.);
00265   G4double logMEM = std::log(x);
00266   G4double pref = twopi*re2/(gamma - 1.0);
00267   // unpolarise XS
00268   G4double sigma0 = 0.;
00269         sigma0 += -gmo2*(gamma - 1.)*x*x*x/3. + gmo2*gamma*x*x;
00270         sigma0 += -(gamma - 1.)*(3.*gamma*(gamma + 2.) +4.)*x;
00271         sigma0 += (gamma*(gamma*(gamma*(4.*gamma - 1.) - 21.) - 7.)+13.)/(3.*(gamma - 1.));
00272         sigma0 /= gpo3;
00273         sigma0 += logMEM*(2. - 1./gpo2);
00274         sigma0 += gamma2/((gamma2 - 1.)*x);
00275   //    longitudinal part
00276   G4double sigma2=0.;
00277         sigma2 += logMEM*gamma*(gamma + 1.)*(2.*gamma + 1.);
00278         sigma2 += gamma*(7.*gamma*(gamma + 1.) - 2.)/3.;
00279         sigma2 += -(3.*gamma + 1.)*(gamma2 + gamma - 1.)*x;
00280         sigma2 += (gamma - 1.)*gamma*(gamma + 3.)*x*x;
00281         sigma2 += -gmo2*(gamma + 3.)*x*x*x/3.;
00282         sigma2 /= gpo3;
00283   //    transverse part
00284   G4double sigma3=0.;
00285         sigma3 += 0.5*(gamma + 1.)*(3.*gamma + 1.)*logMEM;
00286         sigma3 += (gamma*(5.*gamma - 4.) - 13.)/6.;
00287         sigma3 += 0.5*(gamma2 + 3.)*x;
00288         sigma3 += - 2.*(gamma - 1.)*gamma*x*x;  // *AS* changed sign
00289         sigma3 += 2.*gmo2*x*x*x/3.;
00290         sigma3 /= gpo3;
00291   // total cross section
00292   xs+=pref*(sigma0 + sigma2*pol0.z()*pol1.z() + sigma3*(pol0.x()*pol1.x()+pol0.y()*pol1.y()));
00293 
00294   return xs;
00295 }
00296 
00297 
00298 G4StokesVector G4PolarizedBhabhaCrossSection::GetPol2()
00299 {
00300   // Note, mean polarization can not contain correlation
00301   // effects.
00302   return  1./phi0 * phi2;
00303 }
00304 G4StokesVector G4PolarizedBhabhaCrossSection::GetPol3()
00305 {
00306   // Note, mean polarization can not contain correlation
00307   // effects.
00308   return  1./phi0 * phi3;
00309 }

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