Geant4 11.2.2
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 41 of file G4PolarizationHelper.hh.

Member Function Documentation

◆ GetFrame()

◆ GetParticleFrameX()

G4ThreeVector G4PolarizationHelper::GetParticleFrameX ( const G4ThreeVector & uZ)
static

Definition at line 59 of file G4PolarizationHelper.cc.

60{
61 if(uZ.x() == 0. && uZ.y() == 0.)
62 {
63 if(uZ.z() >= 0.)
64 return G4ThreeVector(1., 0., 0.);
65 return G4ThreeVector(-1., 0., 0.);
66 }
67
68 G4double perp = std::sqrt(sqr(uZ.x()) + sqr(uZ.y()));
69 G4double invPerp = uZ.z() / perp;
70 return G4ThreeVector(uZ.x() * invPerp, uZ.y() * invPerp, -perp);
71}
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 48 of file G4PolarizationHelper.cc.

49{
50 if(uZ.x() == 0. && uZ.y() == 0.)
51 {
52 return G4ThreeVector(0., 1., 0.);
53 }
54
55 G4double invPerp = 1. / std::sqrt(sqr(uZ.x()) + sqr(uZ.y()));
56 return G4ThreeVector(-uZ.y() * invPerp, uZ.x() * invPerp, 0);
57}

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

◆ GetRandomFrame()

G4ThreeVector G4PolarizationHelper::GetRandomFrame ( const G4ThreeVector & mom1)
static

Definition at line 73 of file G4PolarizationHelper.cc.

74{
75 G4double phi = 2. * pi * G4UniformRand();
76 G4ThreeVector normal =
77 std::cos(phi) * GetParticleFrameX(mom1) +
78 std::sin(phi) * G4PolarizationHelper::GetParticleFrameY(mom1);
79 return normal;
80}
#define G4UniformRand()
Definition Randomize.hh:52
static G4ThreeVector GetParticleFrameY(const G4ThreeVector &)
static G4ThreeVector GetParticleFrameX(const G4ThreeVector &)

Referenced by G4PolarizedPhotoElectricModel::SampleSecondaries().

◆ GetSpinInPRF()

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

Definition at line 82 of file G4PolarizationHelper.cc.

84{
85 if(uZ.x() == 0. && uZ.y() == 0.)
86 {
87 if(uZ.z() >= 0.)
88 return spin;
89 return G4ThreeVector(-spin.x(), spin.y(), -spin.z());
90 }
91
92 G4double perp = std::sqrt(sqr(uZ.x()) + sqr(uZ.y()));
93 G4double invPerp = 1. / perp;
94
95 G4ThreeVector uX(uZ.x() * uZ.z() * invPerp, uZ.y() * uZ.z() * invPerp, -perp);
96 G4ThreeVector uY(-uZ.y() * invPerp, uZ.x() * invPerp, 0);
97
98 return G4ThreeVector(spin * uX, spin * uY, spin * uZ);
99}

◆ 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 incoming 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
160
161 pol1.rotateUz(dir0);
162
163 G4cout << "========================================\n\n";
164}
G4GLOB_DLL std::ostream G4cout
Hep3Vector & rotateUz(const Hep3Vector &)
static const G4StokesVector P3

Referenced by G4PolarizationMessenger::SetNewValue().

◆ TestPolarizationTransformations()

void G4PolarizationHelper::TestPolarizationTransformations ( )
static

Definition at line 101 of file G4PolarizationHelper.cc.

102{
103 G4double theta = 0.;
104 G4cout << "========================================\n\n";
105 for(G4int i = 0; i <= 10; ++i)
106 {
107 theta = pi * i / 10.;
108 G4ThreeVector zAxis = G4ThreeVector(std::sin(theta), 0., std::cos(theta));
109 if(i == 5)
110 zAxis = G4ThreeVector(1., 0., 0.);
111 if(i == 10)
112 zAxis = G4ThreeVector(0., 0., -1.);
113 G4ThreeVector yAxis = GetParticleFrameY(zAxis);
114
115 G4cout << zAxis << " " << zAxis.mag() << "\n";
116 G4cout << yAxis << " " << yAxis.mag() << "\n";
117 G4ThreeVector xAxis = yAxis.cross(zAxis);
118 G4cout << xAxis << " " << xAxis.mag() << "\n\n";
119 }
120
121 G4cout << "========================================\n\n";
122
123 for(G4int i = 0; i <= 10; ++i)
124 {
125 theta = pi * i / 10.;
126 G4ThreeVector zAxis = G4ThreeVector(0., std::sin(theta), std::cos(theta));
127 if(i == 5)
128 zAxis = G4ThreeVector(0., 1., 0.);
129 if(i == 10)
130 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: