G4PolarizationHelper.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: G4PolarizationHelper.cc 69847 2013-05-16 09:36:18Z gcosmo $
00027 //
00028 // GEANT4 Class file
00029 //
00030 //
00031 // File name:     G4PolarizationHelper
00032 //
00033 // Author:        Andreas Schaelicke
00034 //
00035 // Creation date: 12.08.2006
00036 //
00037 // Modifications:
00038 //
00039 // Class Description:
00040 //
00041 // Provides some basic polarization transformation routines.
00042 //
00043 #include "G4PolarizationHelper.hh"
00044 #include "G4PhysicalConstants.hh"
00045 #include "G4StokesVector.hh"
00046 
00047 
00048 G4ThreeVector G4PolarizationHelper::GetFrame(const G4ThreeVector & mom1, const G4ThreeVector & mom2)
00049 {
00050   G4ThreeVector normal = (mom1.cross(mom2)).unit();
00051   return normal;
00052   //  return 1./normal.mag()*normal;
00053 }
00054 
00055 G4ThreeVector G4PolarizationHelper::GetParticleFrameY(const G4ThreeVector &uZ)
00056 {
00057   // compare also G4ThreeVector::rotateUz()
00058 
00059   if (uZ.x()==0. && uZ.y()==0.) {
00060     return G4ThreeVector(0.,1.,0.);
00061   }
00062 
00063   G4double invPerp = 1./std::sqrt(sqr(uZ.x())+sqr(uZ.y()));
00064   return G4ThreeVector(-uZ.y()*invPerp,uZ.x()*invPerp,0);
00065 }
00066 
00067 G4ThreeVector G4PolarizationHelper::GetParticleFrameX(const G4ThreeVector &uZ)
00068 {
00069   // compare also G4ThreeVector::rotateUz()
00070 
00071   if (uZ.x()==0. && uZ.y()==0.) {
00072     if (uZ.z()>=0.) return G4ThreeVector(1.,0.,0.);
00073     return G4ThreeVector(-1.,0.,0.);
00074   }
00075 
00076   G4double perp    = std::sqrt(sqr(uZ.x())+sqr(uZ.y()));
00077   G4double invPerp = uZ.z()/perp;
00078   return G4ThreeVector(uZ.x()*invPerp,uZ.y()*invPerp,-perp);
00079 }
00080 
00081 G4ThreeVector G4PolarizationHelper::GetRandomFrame(const G4ThreeVector & mom1)
00082 {
00083   G4double phi     =2.*pi*G4UniformRand();
00084   G4ThreeVector normal = std::cos(phi)*GetParticleFrameX(mom1) 
00085     + std::sin(phi)*G4PolarizationHelper::GetParticleFrameY(mom1);
00086   return normal;
00087 }
00088 
00089 
00090 G4ThreeVector G4PolarizationHelper::GetSpinInPRF(const G4ThreeVector &uZ, const G4ThreeVector & spin)
00091 {
00092   // compare also G4ThreeVector::rotateUz()
00093 
00094   if (uZ.x()==0. && uZ.y()==0.) {
00095     if (uZ.z()>=0.) return spin;
00096     return G4ThreeVector(-spin.x(),spin.y(),-spin.z());
00097   }
00098 
00099   G4double perp    = std::sqrt(sqr(uZ.x())+sqr(uZ.y()));
00100   G4double invPerp = 1./perp;
00101 
00102   G4ThreeVector uX(uZ.x()*uZ.z()*invPerp,uZ.y()*uZ.z()*invPerp,-perp);
00103   G4ThreeVector uY(-uZ.y()*invPerp,uZ.x()*invPerp,0); 
00104   
00105   return G4ThreeVector(spin*uX,spin*uY,spin*uZ);
00106 }
00107 
00108 void G4PolarizationHelper::TestPolarizationTransformations()
00109 {
00110   G4double theta=0.;
00111   G4cout<<"========================================\n\n";
00112   for (G4int i=0; i<=10; ++i) {
00113     theta=pi*i/10.;
00114     G4ThreeVector zAxis = G4ThreeVector(std::sin(theta),0.,std::cos(theta));
00115     if (i==5) zAxis = G4ThreeVector(1.,0.,0.);
00116     if (i==10) zAxis = G4ThreeVector(0.,0.,-1.);
00117     G4ThreeVector yAxis = GetParticleFrameY(zAxis);
00118 
00119     G4cout<<zAxis<<" "<<zAxis.mag()<<"\n";
00120     G4cout<<yAxis<<" "<<yAxis.mag()<<"\n";
00121     G4ThreeVector xAxis = yAxis.cross(zAxis);
00122     G4cout<<xAxis<<" "<<xAxis.mag()<<"\n\n";
00123   }
00124 
00125   G4cout<<"========================================\n\n";
00126 
00127   for (G4int i=0; i<=10; ++i) {
00128     theta=pi*i/10.;
00129     G4ThreeVector zAxis = G4ThreeVector(0.,std::sin(theta),std::cos(theta));
00130     if (i==5) zAxis = G4ThreeVector(0.,1.,0.);
00131     if (i==10) zAxis = G4ThreeVector(0.,0.,-1.);
00132     G4ThreeVector yAxis = GetParticleFrameY(zAxis);
00133 
00134     G4cout<<zAxis<<" "<<zAxis.mag()<<"\n";
00135     G4cout<<yAxis<<" "<<yAxis.mag()<<"\n";
00136     G4ThreeVector xAxis = yAxis.cross(zAxis);
00137     G4cout<<xAxis<<" "<<xAxis.mag()<<"\n\n";
00138 
00139     G4cout<<"spat : "<<xAxis*yAxis.cross(zAxis)<<"\n\n";
00140   }
00141   G4cout<<"========================================\n\n";
00142 }
00143 
00144 void G4PolarizationHelper::TestInteractionFrame()
00145 {
00146   // check transformation procedure for polarisation transfer 
00147   // calculation in scattering processes
00148   //  a) transfer target polarisation in beam particle reference frame (PRF)
00149   //  b) calc correct asymmetry w.r.t. scattering plane
00150   //  c) determine incomming polarisation in interaction frame (IF)
00151   //  d) transfer outgoing polarisation from IF to PRF
00152   G4cout<<"========================================\n\n";
00153 
00154   G4double theta=0.;
00155 
00156   G4ThreeVector dir0=G4ThreeVector(0.,0.,1.);
00157   G4ThreeVector dir2=G4ThreeVector(std::sin(theta),0.,std::cos(theta));
00158   
00159   G4StokesVector pol0=G4ThreeVector(0.,0.,1.);
00160   G4StokesVector pol1=G4ThreeVector(0.,0.,1.); 
00161 
00162   pol1.rotateUz(dir0);
00163 
00164   G4cout<<"========================================\n\n";
00165 
00166 
00167 }

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