G4PolarizationHelper Class Reference

#include <G4PolarizationHelper.hh>


Static Public Member Functions

static G4ThreeVector GetFrame (const G4ThreeVector &, const G4ThreeVector &)
static G4ThreeVector GetParticleFrameX (const G4ThreeVector &)
static G4ThreeVector GetParticleFrameY (const G4ThreeVector &)
static G4ThreeVector GetRandomFrame (const G4ThreeVector &)
static G4ThreeVector GetSpinInPRF (const G4ThreeVector &uZ, const G4ThreeVector &spin)
static void TestPolarizationTransformations ()
static void TestInteractionFrame ()


Detailed Description

Definition at line 49 of file G4PolarizationHelper.hh.


Member Function Documentation

G4ThreeVector G4PolarizationHelper::GetFrame ( const G4ThreeVector ,
const G4ThreeVector  
) [static]

Definition at line 48 of file G4PolarizationHelper.cc.

Referenced by G4PolarizedPEEffectModel::SampleSecondaries(), G4PolarizedMollerBhabhaModel::SampleSecondaries(), G4PolarizedGammaConversionModel::SampleSecondaries(), G4PolarizedComptonModel::SampleSecondaries(), G4PolarizedAnnihilationModel::SampleSecondaries(), and G4ePolarizedBremsstrahlungModel::SampleSecondaries().

00049 {
00050   G4ThreeVector normal = (mom1.cross(mom2)).unit();
00051   return normal;
00052   //  return 1./normal.mag()*normal;
00053 }

G4ThreeVector G4PolarizationHelper::GetParticleFrameX ( const G4ThreeVector  )  [static]

Definition at line 67 of file G4PolarizationHelper.cc.

References sqr().

Referenced by G4ePolarizedIonisation::GetMeanFreePath(), G4eplusPolarizedAnnihilation::GetMeanFreePath(), GetRandomFrame(), G4ePolarizedIonisation::PostStepGetPhysicalInteractionLength(), and G4eplusPolarizedAnnihilation::PostStepGetPhysicalInteractionLength().

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 }

G4ThreeVector G4PolarizationHelper::GetParticleFrameY ( const G4ThreeVector  )  [static]

Definition at line 55 of file G4PolarizationHelper.cc.

References sqr().

Referenced by G4ePolarizedIonisation::GetMeanFreePath(), G4eplusPolarizedAnnihilation::GetMeanFreePath(), GetRandomFrame(), G4StokesVector::InvRotateAz(), G4ePolarizedIonisation::PostStepGetPhysicalInteractionLength(), G4eplusPolarizedAnnihilation::PostStepGetPhysicalInteractionLength(), G4StokesVector::RotateAz(), and TestPolarizationTransformations().

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 }

G4ThreeVector G4PolarizationHelper::GetRandomFrame ( const G4ThreeVector  )  [static]

Definition at line 81 of file G4PolarizationHelper.cc.

References G4UniformRand, GetParticleFrameX(), GetParticleFrameY(), and G4INCL::Math::pi.

Referenced by G4PolarizedPEEffectModel::SampleSecondaries().

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 }

G4ThreeVector G4PolarizationHelper::GetSpinInPRF ( const G4ThreeVector uZ,
const G4ThreeVector spin 
) [static]

Definition at line 90 of file G4PolarizationHelper.cc.

References sqr().

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 }

void G4PolarizationHelper::TestInteractionFrame (  )  [static]

Definition at line 144 of file G4PolarizationHelper.cc.

References G4cout.

Referenced by G4PolarizationMessenger::SetNewValue().

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 }

void G4PolarizationHelper::TestPolarizationTransformations (  )  [static]

Definition at line 108 of file G4PolarizationHelper.cc.

References G4cout, GetParticleFrameY(), and G4INCL::Math::pi.

Referenced by G4PolarizationMessenger::SetNewValue().

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 }


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