Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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

◆ GetFrame()

◆ GetParticleFrameX()

G4ThreeVector G4PolarizationHelper::GetParticleFrameX ( const G4ThreeVector uZ)
static

Definition at line 67 of file G4PolarizationHelper.cc.

68{
69 // compare also G4ThreeVector::rotateUz()
70
71 if (uZ.x()==0. && uZ.y()==0.) {
72 if (uZ.z()>=0.) return G4ThreeVector(1.,0.,0.);
73 return G4ThreeVector(-1.,0.,0.);
74 }
75
76 G4double perp = std::sqrt(sqr(uZ.x())+sqr(uZ.y()));
77 G4double invPerp = uZ.z()/perp;
78 return G4ThreeVector(uZ.x()*invPerp,uZ.y()*invPerp,-perp);
79}
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition: G4Types.hh:64
double z() const
double x() const
double y() const
T sqr(const T &x)
Definition: templates.hh:145

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

◆ GetParticleFrameY()

G4ThreeVector G4PolarizationHelper::GetParticleFrameY ( const G4ThreeVector uZ)
static

Definition at line 55 of file G4PolarizationHelper.cc.

56{
57 // compare also G4ThreeVector::rotateUz()
58
59 if (uZ.x()==0. && uZ.y()==0.) {
60 return G4ThreeVector(0.,1.,0.);
61 }
62
63 G4double invPerp = 1./std::sqrt(sqr(uZ.x())+sqr(uZ.y()));
64 return G4ThreeVector(-uZ.y()*invPerp,uZ.x()*invPerp,0);
65}

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

◆ GetRandomFrame()

G4ThreeVector G4PolarizationHelper::GetRandomFrame ( const G4ThreeVector mom1)
static

Definition at line 81 of file G4PolarizationHelper.cc.

82{
83 G4double phi =2.*pi*G4UniformRand();
84 G4ThreeVector normal = std::cos(phi)*GetParticleFrameX(mom1)
85 + std::sin(phi)*G4PolarizationHelper::GetParticleFrameY(mom1);
86 return normal;
87}
#define G4UniformRand()
Definition: Randomize.hh:53
static G4ThreeVector GetParticleFrameY(const G4ThreeVector &)
static G4ThreeVector GetParticleFrameX(const G4ThreeVector &)
const G4double pi

Referenced by G4PolarizedPEEffectModel::SampleSecondaries().

◆ GetSpinInPRF()

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

Definition at line 90 of file G4PolarizationHelper.cc.

91{
92 // compare also G4ThreeVector::rotateUz()
93
94 if (uZ.x()==0. && uZ.y()==0.) {
95 if (uZ.z()>=0.) return spin;
96 return G4ThreeVector(-spin.x(),spin.y(),-spin.z());
97 }
98
99 G4double perp = std::sqrt(sqr(uZ.x())+sqr(uZ.y()));
100 G4double invPerp = 1./perp;
101
102 G4ThreeVector uX(uZ.x()*uZ.z()*invPerp,uZ.y()*uZ.z()*invPerp,-perp);
103 G4ThreeVector uY(-uZ.y()*invPerp,uZ.x()*invPerp,0);
104
105 return G4ThreeVector(spin*uX,spin*uY,spin*uZ);
106}

◆ TestInteractionFrame()

void G4PolarizationHelper::TestInteractionFrame ( )
static

Definition at line 144 of file G4PolarizationHelper.cc.

145{
146 // check transformation procedure for polarisation transfer
147 // calculation in scattering processes
148 // a) transfer target polarisation in beam particle reference frame (PRF)
149 // b) calc correct asymmetry w.r.t. scattering plane
150 // c) determine incomming polarisation in interaction frame (IF)
151 // d) transfer outgoing polarisation from IF to PRF
152 G4cout<<"========================================\n\n";
153
154 G4double theta=0.;
155
156 G4ThreeVector dir0=G4ThreeVector(0.,0.,1.);
157 G4ThreeVector dir2=G4ThreeVector(std::sin(theta),0.,std::cos(theta));
158
159 G4StokesVector pol0=G4ThreeVector(0.,0.,1.);
160 G4StokesVector pol1=G4ThreeVector(0.,0.,1.);
161
162 pol1.rotateUz(dir0);
163
164 G4cout<<"========================================\n\n";
165
166
167}
G4DLLIMPORT std::ostream G4cout
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:72

Referenced by G4PolarizationMessenger::SetNewValue().

◆ TestPolarizationTransformations()

void G4PolarizationHelper::TestPolarizationTransformations ( )
static

Definition at line 108 of file G4PolarizationHelper.cc.

109{
110 G4double theta=0.;
111 G4cout<<"========================================\n\n";
112 for (G4int i=0; i<=10; ++i) {
113 theta=pi*i/10.;
114 G4ThreeVector zAxis = G4ThreeVector(std::sin(theta),0.,std::cos(theta));
115 if (i==5) zAxis = G4ThreeVector(1.,0.,0.);
116 if (i==10) zAxis = G4ThreeVector(0.,0.,-1.);
117 G4ThreeVector yAxis = GetParticleFrameY(zAxis);
118
119 G4cout<<zAxis<<" "<<zAxis.mag()<<"\n";
120 G4cout<<yAxis<<" "<<yAxis.mag()<<"\n";
121 G4ThreeVector xAxis = yAxis.cross(zAxis);
122 G4cout<<xAxis<<" "<<xAxis.mag()<<"\n\n";
123 }
124
125 G4cout<<"========================================\n\n";
126
127 for (G4int i=0; i<=10; ++i) {
128 theta=pi*i/10.;
129 G4ThreeVector zAxis = G4ThreeVector(0.,std::sin(theta),std::cos(theta));
130 if (i==5) zAxis = G4ThreeVector(0.,1.,0.);
131 if (i==10) zAxis = G4ThreeVector(0.,0.,-1.);
132 G4ThreeVector yAxis = GetParticleFrameY(zAxis);
133
134 G4cout<<zAxis<<" "<<zAxis.mag()<<"\n";
135 G4cout<<yAxis<<" "<<yAxis.mag()<<"\n";
136 G4ThreeVector xAxis = yAxis.cross(zAxis);
137 G4cout<<xAxis<<" "<<xAxis.mag()<<"\n\n";
138
139 G4cout<<"spat : "<<xAxis*yAxis.cross(zAxis)<<"\n\n";
140 }
141 G4cout<<"========================================\n\n";
142}
int G4int
Definition: G4Types.hh:66
double mag() const

Referenced by G4PolarizationMessenger::SetNewValue().


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