Geant4 10.7.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 48 of file G4PolarizationHelper.hh.

Member Function Documentation

◆ GetFrame()

G4ThreeVector G4PolarizationHelper::GetFrame ( const G4ThreeVector mom1,
const G4ThreeVector mom2 
)
static

◆ GetParticleFrameX()

G4ThreeVector G4PolarizationHelper::GetParticleFrameX ( const G4ThreeVector uZ)
static

Definition at line 66 of file G4PolarizationHelper.cc.

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

Referenced by GetRandomFrame().

◆ GetParticleFrameY()

G4ThreeVector G4PolarizationHelper::GetParticleFrameY ( const G4ThreeVector uZ)
static

Definition at line 54 of file G4PolarizationHelper.cc.

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

Referenced by GetRandomFrame(), G4StokesVector::InvRotateAz(), G4StokesVector::RotateAz(), and TestPolarizationTransformations().

◆ GetRandomFrame()

G4ThreeVector G4PolarizationHelper::GetRandomFrame ( const G4ThreeVector mom1)
static

Definition at line 80 of file G4PolarizationHelper.cc.

81{
82 G4double phi =2.*pi*G4UniformRand();
83 G4ThreeVector normal = std::cos(phi)*GetParticleFrameX(mom1)
84 + std::sin(phi)*G4PolarizationHelper::GetParticleFrameY(mom1);
85 return normal;
86}
#define G4UniformRand()
Definition: Randomize.hh:52
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 89 of file G4PolarizationHelper.cc.

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

◆ TestInteractionFrame()

void G4PolarizationHelper::TestInteractionFrame ( )
static

Definition at line 143 of file G4PolarizationHelper.cc.

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

Referenced by G4PolarizationMessenger::SetNewValue().

◆ TestPolarizationTransformations()

void G4PolarizationHelper::TestPolarizationTransformations ( )
static

Definition at line 107 of file G4PolarizationHelper.cc.

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

Referenced by G4PolarizationMessenger::SetNewValue().


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