Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4StokesVector Class Reference

#include <G4StokesVector.hh>

+ Inheritance diagram for G4StokesVector:

Public Member Functions

 G4StokesVector ()
 
 G4StokesVector (const G4ThreeVector &v)
 
 ~G4StokesVector ()=default
 
G4bool IsZero () const
 
G4double p1 () const
 
G4double p2 () const
 
G4double p3 () const
 
G4double Transverse () const
 
G4ThreeVector PolSqr () const
 
G4ThreeVector PolSqrt () const
 
G4ThreeVector PolError (const G4StokesVector &sum2, long n)
 
G4ThreeVector PolDiv (const G4StokesVector &)
 
void SetPhoton ()
 
void RotateAz (G4ThreeVector nInteractionFrame, G4ThreeVector particleDirection)
 
void InvRotateAz (G4ThreeVector nInteractionFrame, G4ThreeVector particleDirection)
 
void RotateAz (G4double cosphi, G4double sinphi)
 
G4double GetBeta ()
 
void DiceUniform ()
 
void DiceP1 ()
 
void DiceP2 ()
 
void DiceP3 ()
 
void FlipP3 ()
 
- Public Member Functions inherited from CLHEP::Hep3Vector
 Hep3Vector ()
 
 Hep3Vector (double x)
 
 Hep3Vector (double x, double y)
 
 Hep3Vector (double x, double y, double z)
 
 Hep3Vector (const Hep3Vector &)
 
 Hep3Vector (Hep3Vector &&)=default
 
 ~Hep3Vector ()
 
double operator() (int) const
 
double operator[] (int) const
 
double & operator() (int)
 
double & operator[] (int)
 
double x () const
 
double y () const
 
double z () const
 
void setX (double)
 
void setY (double)
 
void setZ (double)
 
void set (double x, double y, double z)
 
double phi () const
 
double theta () const
 
double cosTheta () const
 
double cos2Theta () const
 
double mag2 () const
 
double mag () const
 
void setPhi (double)
 
void setTheta (double)
 
void setMag (double)
 
double perp2 () const
 
double perp () const
 
void setPerp (double)
 
void setCylTheta (double)
 
double perp2 (const Hep3Vector &) const
 
double perp (const Hep3Vector &) const
 
Hep3Vectoroperator= (const Hep3Vector &)
 
Hep3Vectoroperator= (Hep3Vector &&)=default
 
bool operator== (const Hep3Vector &) const
 
bool operator!= (const Hep3Vector &) const
 
bool isNear (const Hep3Vector &, double epsilon=tolerance) const
 
double howNear (const Hep3Vector &v) const
 
double deltaR (const Hep3Vector &v) const
 
Hep3Vectoroperator+= (const Hep3Vector &)
 
Hep3Vectoroperator-= (const Hep3Vector &)
 
Hep3Vector operator- () const
 
Hep3Vectoroperator*= (double)
 
Hep3Vectoroperator/= (double)
 
Hep3Vector unit () const
 
Hep3Vector orthogonal () const
 
double dot (const Hep3Vector &) const
 
Hep3Vector cross (const Hep3Vector &) const
 
double angle (const Hep3Vector &) const
 
double pseudoRapidity () const
 
void setEta (double p)
 
void setCylEta (double p)
 
Hep3VectorrotateX (double)
 
Hep3VectorrotateY (double)
 
Hep3VectorrotateZ (double)
 
Hep3VectorrotateUz (const Hep3Vector &)
 
Hep3Vectorrotate (double, const Hep3Vector &)
 
Hep3Vectoroperator*= (const HepRotation &)
 
Hep3Vectortransform (const HepRotation &)
 
void setRThetaPhi (double r, double theta, double phi)
 
void setREtaPhi (double r, double eta, double phi)
 
void setRhoPhiZ (double rho, double phi, double z)
 
void setRhoPhiTheta (double rho, double phi, double theta)
 
void setRhoPhiEta (double rho, double phi, double eta)
 
double getX () const
 
double getY () const
 
double getZ () const
 
double getR () const
 
double getTheta () const
 
double getPhi () const
 
double r () const
 
double rho () const
 
double getRho () const
 
double eta () const
 
double getEta () const
 
void setR (double s)
 
void setRho (double s)
 
int compare (const Hep3Vector &v) const
 
bool operator> (const Hep3Vector &v) const
 
bool operator< (const Hep3Vector &v) const
 
bool operator>= (const Hep3Vector &v) const
 
bool operator<= (const Hep3Vector &v) const
 
double diff2 (const Hep3Vector &v) const
 
bool isParallel (const Hep3Vector &v, double epsilon=tolerance) const
 
bool isOrthogonal (const Hep3Vector &v, double epsilon=tolerance) const
 
double howParallel (const Hep3Vector &v) const
 
double howOrthogonal (const Hep3Vector &v) const
 
double beta () const
 
double gamma () const
 
double coLinearRapidity () const
 
double angle () const
 
double theta (const Hep3Vector &v2) const
 
double cosTheta (const Hep3Vector &v2) const
 
double cos2Theta (const Hep3Vector &v2) const
 
Hep3Vector project () const
 
Hep3Vector project (const Hep3Vector &v2) const
 
Hep3Vector perpPart () const
 
Hep3Vector perpPart (const Hep3Vector &v2) const
 
double rapidity () const
 
double rapidity (const Hep3Vector &v2) const
 
double eta (const Hep3Vector &v2) const
 
double polarAngle (const Hep3Vector &v2) const
 
double deltaPhi (const Hep3Vector &v2) const
 
double azimAngle (const Hep3Vector &v2) const
 
double polarAngle (const Hep3Vector &v2, const Hep3Vector &ref) const
 
double azimAngle (const Hep3Vector &v2, const Hep3Vector &ref) const
 
Hep3Vectorrotate (const Hep3Vector &axis, double delta)
 
Hep3Vectorrotate (const HepAxisAngle &ax)
 
Hep3Vectorrotate (const HepEulerAngles &e)
 
Hep3Vectorrotate (double phi, double theta, double psi)
 

Static Public Attributes

static const G4StokesVector ZERO
 
static const G4StokesVector P1
 
static const G4StokesVector P2
 
static const G4StokesVector P3
 
static const G4StokesVector M1
 
static const G4StokesVector M2
 
static const G4StokesVector M3
 
- Static Public Attributes inherited from CLHEP::Hep3Vector
static const int ToleranceTicks = 100
 

Additional Inherited Members

- Public Types inherited from CLHEP::Hep3Vector
enum  {
  X =0 , Y =1 , Z =2 , NUM_COORDINATES =3 ,
  SIZE =NUM_COORDINATES
}
 
- Static Public Member Functions inherited from CLHEP::Hep3Vector
static double setTolerance (double tol)
 
static double getTolerance ()
 
- Protected Member Functions inherited from CLHEP::Hep3Vector
void setSpherical (double r, double theta, double phi)
 
void setCylindrical (double r, double phi, double z)
 
double negativeInfinity () const
 
- Protected Attributes inherited from CLHEP::Hep3Vector
double data [3]
 
- Static Protected Attributes inherited from CLHEP::Hep3Vector
static DLL_API double tolerance = Hep3Vector::ToleranceTicks * 2.22045e-16
 

Detailed Description

Definition at line 45 of file G4StokesVector.hh.

Constructor & Destructor Documentation

◆ G4StokesVector() [1/2]

G4StokesVector::G4StokesVector ( )

Definition at line 55 of file G4StokesVector.cc.

57 , fIsPhoton(false)
58{}
CLHEP::Hep3Vector G4ThreeVector

Referenced by PolError().

◆ G4StokesVector() [2/2]

G4StokesVector::G4StokesVector ( const G4ThreeVector & v)
explicit

Definition at line 60 of file G4StokesVector.cc.

61 : G4ThreeVector(v)
62 , fIsPhoton(false)
63{}

◆ ~G4StokesVector()

G4StokesVector::~G4StokesVector ( )
default

Member Function Documentation

◆ DiceP1()

void G4StokesVector::DiceP1 ( )

Definition at line 171 of file G4StokesVector.cc.

172{
173 if(G4UniformRand() > 0.5)
174 setX(1.);
175 else
176 setX(-1.);
177 setY(0.);
178 setZ(0.);
179}
#define G4UniformRand()
Definition Randomize.hh:52
void setY(double)
void setZ(double)
void setX(double)

◆ DiceP2()

void G4StokesVector::DiceP2 ( )

Definition at line 181 of file G4StokesVector.cc.

182{
183 setX(0.);
184 if(G4UniformRand() > 0.5)
185 setY(1.);
186 else
187 setY(-1.);
188 setZ(0.);
189}

◆ DiceP3()

void G4StokesVector::DiceP3 ( )

Definition at line 191 of file G4StokesVector.cc.

192{
193 setX(0.);
194 setY(0.);
195 if(G4UniformRand() > 0.5)
196 setZ(1.);
197 else
198 setZ(-1.);
199}

◆ DiceUniform()

void G4StokesVector::DiceUniform ( )

Definition at line 161 of file G4StokesVector.cc.

162{
163 G4double costheta = 2. * G4UniformRand() - 1.;
164 G4double sintheta = std::sqrt(1. - costheta * costheta);
165 G4double aphi = 2. * CLHEP::pi * G4UniformRand();
166 setX(std::sin(aphi) * sintheta);
167 setY(std::cos(aphi) * sintheta);
168 setZ(costheta);
169}
double G4double
Definition G4Types.hh:83

◆ FlipP3()

void G4StokesVector::FlipP3 ( )

Definition at line 201 of file G4StokesVector.cc.

201{ setZ(-z()); }
double z() const

◆ GetBeta()

G4double G4StokesVector::GetBeta ( )

Definition at line 151 of file G4StokesVector.cc.

152{
153 G4double bet = getPhi();
154 if(fIsPhoton)
155 {
156 bet *= 0.5;
157 }
158 return bet;
159}
double getPhi() const

◆ InvRotateAz()

void G4StokesVector::InvRotateAz ( G4ThreeVector nInteractionFrame,
G4ThreeVector particleDirection )

Definition at line 100 of file G4StokesVector.cc.

102{
103 // note if incoming particle is on z-axis,
104 // we might encounter some nummerical problems, since
105 // nInteratonFrame and yParticleFrame are actually (almost) the same momentum
106 // and the normalization is only good to 10^-12 !
107
108 G4ThreeVector yParticleFrame =
110 G4double cosphi = yParticleFrame * nInteractionFrame;
111
112 if(cosphi > 1. + 1.e-8 || cosphi < -1. - 1.e-8)
113 {
115 ed << " warning G4StokesVector::RotateAz cosphi>1 or cosphi<-1\n";
116 G4Exception("G4StokesVector::InvRotateAz", "pol030", JustWarning, ed);
117 }
118 if(cosphi > 1.)
119 cosphi = 1.;
120 else if(cosphi < -1.)
121 cosphi = -1.;
122
123 // check sign once more!
124 G4double hel =
125 (yParticleFrame.cross(nInteractionFrame) * particleDirection) > 0. ? 1.
126 : -1.;
127 G4double sinphi = hel * std::sqrt(std::fabs(1. - cosphi * cosphi));
128 RotateAz(cosphi, -sinphi);
129}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
Hep3Vector cross(const Hep3Vector &) const
static G4ThreeVector GetParticleFrameY(const G4ThreeVector &)
void RotateAz(G4ThreeVector nInteractionFrame, G4ThreeVector particleDirection)

Referenced by G4PolarizedAnnihilationModel::SampleSecondaries(), G4PolarizedBremsstrahlungModel::SampleSecondaries(), G4PolarizedComptonModel::SampleSecondaries(), G4PolarizedGammaConversionModel::SampleSecondaries(), G4PolarizedIonisationModel::SampleSecondaries(), and G4PolarizedPhotoElectricModel::SampleSecondaries().

◆ IsZero()

G4bool G4StokesVector::IsZero ( ) const

◆ p1()

◆ p2()

◆ p3()

◆ PolDiv()

G4ThreeVector G4StokesVector::PolDiv ( const G4StokesVector & b)

Definition at line 213 of file G4StokesVector.cc.

214{
215 return G4ThreeVector(b.x() != 0. ? x() / b.x() : 11111.,
216 b.y() != 0. ? y() / b.y() : 11111.,
217 b.z() != 0. ? z() / b.z() : 11111.);
218}

◆ PolError()

G4ThreeVector G4StokesVector::PolError ( const G4StokesVector & sum2,
long n )

Definition at line 203 of file G4StokesVector.cc.

204{
205 // delta x = sqrt[ ( <x^2> - <x>^2 )/(n-1) ]
206 G4ThreeVector mean = (1. / n) * G4ThreeVector(*this);
207 G4ThreeVector polsqr = G4StokesVector(mean).PolSqr();
208 G4ThreeVector result =
209 G4StokesVector((1. / (n - 1.) * ((1. / n) * sum2 - polsqr))).PolSqrt();
210 return result;
211}

◆ PolSqr()

G4ThreeVector G4StokesVector::PolSqr ( ) const
inline

Definition at line 60 of file G4StokesVector.hh.

61 {
62 return G4ThreeVector(x() * x(), y() * y(), z() * z());
63 }

◆ PolSqrt()

G4ThreeVector G4StokesVector::PolSqrt ( ) const
inline

Definition at line 64 of file G4StokesVector.hh.

65 {
66 return G4ThreeVector(std::sqrt(x()), std::sqrt(y()), std::sqrt(z()));
67 }

◆ RotateAz() [1/2]

void G4StokesVector::RotateAz ( G4double cosphi,
G4double sinphi )

Definition at line 131 of file G4StokesVector.cc.

132{
133 if(!fIsPhoton)
134 {
135 G4double xsi1 = cosphi * p1() + sinphi * p2();
136 G4double xsi2 = -sinphi * p1() + cosphi * p2();
137 setX(xsi1);
138 setY(xsi2);
139 return;
140 }
141
142 G4double sin2phi = 2. * cosphi * sinphi;
143 G4double cos2phi = cosphi * cosphi - sinphi * sinphi;
144
145 G4double xsi1 = cos2phi * p1() + sin2phi * p2();
146 G4double xsi2 = -sin2phi * p1() + cos2phi * p2();
147 setX(xsi1);
148 setY(xsi2);
149}
G4double p1() const
G4double p2() const

◆ RotateAz() [2/2]

void G4StokesVector::RotateAz ( G4ThreeVector nInteractionFrame,
G4ThreeVector particleDirection )

Definition at line 67 of file G4StokesVector.cc.

69{
70 G4ThreeVector yParticleFrame =
72
73 G4double cosphi = yParticleFrame * nInteractionFrame;
74 if(cosphi > (1. + 1.e-8) || cosphi < (-1. - 1.e-8))
75 {
77 ed << " warning G4StokesVector::RotateAz cosphi>1 or cosphi<-1\n"
78 << " cosphi=" << cosphi << "\n"
79 << " zAxis=" << particleDirection << " (" << particleDirection.mag()
80 << ")\n"
81 << " yAxis=" << yParticleFrame << " (" << yParticleFrame.mag() << ")\n"
82 << " nAxis=" << nInteractionFrame << " (" << nInteractionFrame.mag()
83 << ")\n";
84 G4Exception("G4StokesVector::RotateAz", "pol030", JustWarning, ed);
85 }
86 if(cosphi > 1.)
87 cosphi = 1.;
88 else if(cosphi < -1.)
89 cosphi = -1.;
90
91 G4double hel =
92 (yParticleFrame.cross(nInteractionFrame) * particleDirection) > 0. ? 1.
93 : -1.;
94
95 G4double sinphi = hel * std::sqrt(1. - cosphi * cosphi);
96
97 RotateAz(cosphi, sinphi);
98}
double mag() const

Referenced by InvRotateAz(), RotateAz(), G4PolarizedAnnihilationModel::SampleSecondaries(), G4PolarizedBremsstrahlungModel::SampleSecondaries(), G4PolarizedComptonModel::SampleSecondaries(), G4PolarizedGammaConversionModel::SampleSecondaries(), G4PolarizedIonisationModel::SampleSecondaries(), and G4PolarizedPhotoElectricModel::SampleSecondaries().

◆ SetPhoton()

◆ Transverse()

G4double G4StokesVector::Transverse ( ) const
inline

Definition at line 58 of file G4StokesVector.hh.

58{ return perp(); }
double perp() const

Member Data Documentation

◆ M1

const G4StokesVector G4StokesVector::M1
static
Initial value:
=

Definition at line 94 of file G4StokesVector.hh.

◆ M2

const G4StokesVector G4StokesVector::M2
static
Initial value:
=

Definition at line 95 of file G4StokesVector.hh.

◆ M3

const G4StokesVector G4StokesVector::M3
static
Initial value:
=

Definition at line 96 of file G4StokesVector.hh.

◆ P1

const G4StokesVector G4StokesVector::P1
static

◆ P2

const G4StokesVector G4StokesVector::P2
static

◆ P3

◆ ZERO


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