Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions | Static Public Attributes
G4StokesVector Class Reference

#include <G4StokesVector.hh>

Inheritance diagram for G4StokesVector:
CLHEP::Hep3Vector

Public Member Functions

 G4StokesVector ()
 
 G4StokesVector (const G4ThreeVector &v)
 
virtual ~G4StokesVector ()
 
G4double p1 () const
 
G4double p2 () const
 
G4double p3 () const
 
G4bool IsZero () 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 ()
 
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 &)
 
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 =G4ThreeVector(0.,0.,0.)
 
static const G4StokesVector P1 =G4ThreeVector(1.,0.,0.)
 
static const G4StokesVector P2 =G4ThreeVector(0.,1.,0.)
 
static const G4StokesVector P3 =G4ThreeVector(0.,0.,1.)
 
static const G4StokesVector M1 =G4ThreeVector(-1.,0.,0.)
 
static const G4StokesVector M2 =G4ThreeVector(0.,-1.,0.)
 
static const G4StokesVector M3 =G4ThreeVector(0.,0.,-1.)
 

Additional Inherited Members

- Public Types inherited from CLHEP::Hep3Vector
enum  {
  X =0, Y =1, Z =2, NUM_COORDINATES =3,
  SIZE =NUM_COORDINATES
}
 
enum  { ToleranceTicks = 100 }
 
- 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 dx
 
double dy
 
double dz
 
- Static Protected Attributes inherited from CLHEP::Hep3Vector
static DLL_API double tolerance = Hep3Vector::ToleranceTicks * 2.22045e-16
 

Detailed Description

Definition at line 58 of file G4StokesVector.hh.

Constructor & Destructor Documentation

G4StokesVector::G4StokesVector ( )

Definition at line 57 of file G4StokesVector.cc.

Referenced by PolError().

58  : G4ThreeVector(),isPhoton(false)
59 {
60 }
CLHEP::Hep3Vector G4ThreeVector
G4StokesVector::G4StokesVector ( const G4ThreeVector v)

Definition at line 62 of file G4StokesVector.cc.

63  : G4ThreeVector(v),isPhoton(false)
64 {
65 }
CLHEP::Hep3Vector G4ThreeVector
G4StokesVector::~G4StokesVector ( )
virtual

Definition at line 67 of file G4StokesVector.cc.

68 {
69 }

Member Function Documentation

void G4StokesVector::DiceP1 ( )

Definition at line 168 of file G4StokesVector.cc.

References G4UniformRand, CLHEP::Hep3Vector::setX(), CLHEP::Hep3Vector::setY(), and CLHEP::Hep3Vector::setZ().

169 {
170  if (G4UniformRand()>0.5) setX(1.);
171  else setX(-1.);
172  setY(0.);
173  setZ(0.);
174 }
void setY(double)
void setZ(double)
void setX(double)
#define G4UniformRand()
Definition: Randomize.hh:87
void G4StokesVector::DiceP2 ( )

Definition at line 176 of file G4StokesVector.cc.

References G4UniformRand, CLHEP::Hep3Vector::setX(), CLHEP::Hep3Vector::setY(), and CLHEP::Hep3Vector::setZ().

177 {
178  setX(0.);
179  if (G4UniformRand()>0.5) setY(1.);
180  else setY(-1.);
181  setZ(0.);
182 }
void setY(double)
void setZ(double)
void setX(double)
#define G4UniformRand()
Definition: Randomize.hh:87
void G4StokesVector::DiceP3 ( )

Definition at line 184 of file G4StokesVector.cc.

References G4UniformRand, CLHEP::Hep3Vector::setX(), CLHEP::Hep3Vector::setY(), and CLHEP::Hep3Vector::setZ().

185 {
186  setX(0.);
187  setY(0.);
188  if (G4UniformRand()>0.5) setZ(1.);
189  else setZ(-1.);
190 }
void setY(double)
void setZ(double)
void setX(double)
#define G4UniformRand()
Definition: Randomize.hh:87
void G4StokesVector::DiceUniform ( )

Definition at line 158 of file G4StokesVector.cc.

References G4UniformRand, python.hepunit::pi, CLHEP::Hep3Vector::setX(), CLHEP::Hep3Vector::setY(), and CLHEP::Hep3Vector::setZ().

159 {
160  G4double costheta=2.*G4UniformRand()-1.;
161  G4double sintheta=std::sqrt(1.-costheta*costheta);
162  G4double aphi =2.*pi*G4UniformRand();
163  setX(std::sin(aphi)*sintheta);
164  setY(std::cos(aphi)*sintheta);
165  setZ(costheta);
166 }
void setY(double)
void setZ(double)
void setX(double)
#define G4UniformRand()
Definition: Randomize.hh:87
double G4double
Definition: G4Types.hh:76
void G4StokesVector::FlipP3 ( )

Definition at line 192 of file G4StokesVector.cc.

References CLHEP::Hep3Vector::setZ(), and CLHEP::Hep3Vector::z().

193 {
194  setZ(-z());
195 }
double z() const
void setZ(double)
G4double G4StokesVector::GetBeta ( )

Definition at line 151 of file G4StokesVector.cc.

References CLHEP::Hep3Vector::getPhi().

152 {
153  G4double bet=getPhi();
154  if (isPhoton) { bet *= 0.5; }
155  return bet;
156 }
double getPhi() const
double G4double
Definition: G4Types.hh:76
void G4StokesVector::InvRotateAz ( G4ThreeVector  nInteractionFrame,
G4ThreeVector  particleDirection 
)

Definition at line 108 of file G4StokesVector.cc.

References CLHEP::Hep3Vector::cross(), G4cout, G4PolarizationHelper::GetParticleFrameY(), and RotateAz().

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

110 {
111  // note if incomming particle is on z-axis,
112  // we might encounter some nummerical problems, since
113  // nInteratonFrame and yParticleFrame are actually (almost) the same momentum
114  // and the normalization is only good to 10^-12 !
115 
116  G4ThreeVector yParticleFrame =
117  G4PolarizationHelper::GetParticleFrameY(particleDirection);
118  G4double cosphi=yParticleFrame*nInteractionFrame;
119 
120  if (cosphi>1.+1.e-8 || cosphi<-1.-1.e-8) {
121  G4cout<<" warning G4StokesVector::RotateAz cosphi>1 or cosphi<-1\n";
122  }
123  if (cosphi>1) cosphi=1.;
124  else if (cosphi<-1)cosphi=-1.;
125 
126  // check sign once more!
127  G4double hel=(yParticleFrame.cross(nInteractionFrame)*particleDirection)>0?1.:-1.;
128  G4double sinphi=hel*std::sqrt(std::fabs(1.-cosphi*cosphi));
129  RotateAz(cosphi,-sinphi);
130 }
G4GLOB_DLL std::ostream G4cout
static G4ThreeVector GetParticleFrameY(const G4ThreeVector &)
Hep3Vector cross(const Hep3Vector &) const
double G4double
Definition: G4Types.hh:76
void RotateAz(G4ThreeVector nInteractionFrame, G4ThreeVector particleDirection)
G4bool G4StokesVector::IsZero ( ) const
inline
G4double G4StokesVector::p1 ( ) const
inline
G4double G4StokesVector::p2 ( ) const
inline
G4double G4StokesVector::p3 ( ) const
inline
G4ThreeVector G4StokesVector::PolDiv ( const G4StokesVector b)

Definition at line 204 of file G4StokesVector.cc.

References CLHEP::Hep3Vector::x(), CLHEP::Hep3Vector::y(), and CLHEP::Hep3Vector::z().

205  {return G4ThreeVector(b.x()!=0. ? x()/b.x() : 11111.,
206  b.y()!=0. ? y()/b.y() : 11111.,
207  b.z()!=0. ? z()/b.z() : 11111.);}
CLHEP::Hep3Vector G4ThreeVector
double x() const
double z() const
double y() const
G4ThreeVector G4StokesVector::PolError ( const G4StokesVector sum2,
long  n 
)

Definition at line 197 of file G4StokesVector.cc.

References G4StokesVector(), n, PolSqr(), and PolSqrt().

198 {
199  // delta x = sqrt[ ( <x^2> - <x>^2 )/(n-1) ]
200  G4StokesVector mean=(1./n)*(*this);
201  return G4StokesVector((1./(n-1.)*((1./n)*sum2 - mean.PolSqr()))).PolSqrt();
202 }
G4ThreeVector PolSqr() const
const G4int n
G4ThreeVector PolSqrt() const
G4ThreeVector G4StokesVector::PolSqr ( ) const
inline

Definition at line 81 of file G4StokesVector.hh.

References CLHEP::Hep3Vector::x(), CLHEP::Hep3Vector::y(), and CLHEP::Hep3Vector::z().

Referenced by PolError().

81  {
82  return G4ThreeVector(x()*x(),y()*y(),z()*z());
83  }
CLHEP::Hep3Vector G4ThreeVector
double x() const
double z() const
double y() const
G4ThreeVector G4StokesVector::PolSqrt ( ) const
inline

Definition at line 84 of file G4StokesVector.hh.

References CLHEP::Hep3Vector::x(), CLHEP::Hep3Vector::y(), and CLHEP::Hep3Vector::z().

Referenced by PolError().

84  {
85  return G4ThreeVector(std::sqrt(x()),std::sqrt(y()),std::sqrt(z()));
86  }
CLHEP::Hep3Vector G4ThreeVector
double x() const
double z() const
double y() const
void G4StokesVector::RotateAz ( G4ThreeVector  nInteractionFrame,
G4ThreeVector  particleDirection 
)

Definition at line 71 of file G4StokesVector.cc.

References CLHEP::Hep3Vector::cross(), G4cout, G4endl, G4PolarizationHelper::GetParticleFrameY(), and CLHEP::Hep3Vector::mag().

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

73 {
74  G4ThreeVector yParticleFrame =
75  G4PolarizationHelper::GetParticleFrameY(particleDirection);
76 
77 
78  G4double cosphi=yParticleFrame*nInteractionFrame;
79  if (cosphi>(1.+1.e-8) || cosphi<(-1.-1.e-8)) {
80  G4cout<<" warning G4StokesVector::RotateAz cosphi>1 or cosphi<-1\n"
81  <<" cosphi="<<cosphi<<"\n"
82  <<" zAxis="<<particleDirection<<" ("<<particleDirection.mag()<<")\n"
83  <<" yAxis="<<yParticleFrame<<" ("<<yParticleFrame.mag()<<")\n"
84  <<" nAxis="<<nInteractionFrame<<" ("
85  <<nInteractionFrame.mag()<<")"<<G4endl;
86  }
87  if (cosphi>1.) cosphi=1.;
88  else if (cosphi<-1.) cosphi=-1.;
89 
90 // G4cout<<" cosphi="<<cosphi<<"\n"
91 // <<" zAxis="<<particleDirection<<" ("<<particleDirection.mag()<<")\n"
92 // <<" yAxis="<<yParticleFrame<<" ("<<yParticleFrame.mag()<<","<<(yParticleFrame*particleDirection)<<")\n"
93 // <<" nAxis="<<nInteractionFrame<<" ("
94 // <<nInteractionFrame.mag()<<")"<<G4endl;
95 
96  // G4double hel=sgn(cross(yParticleFrame*nInteractionFrame)*zInteractionFrame);
97  // Why not particleDirection instead of zInteractionFrame ???!!!
98  // -> is the same, since SYSIN is called with p1, and p2 as first parameter!
99  G4double hel=(yParticleFrame.cross(nInteractionFrame)*particleDirection)>0?1.:-1.;
100 
101  G4double sinphi=hel*std::sqrt(1.-cosphi*cosphi);
102  // G4cout<<" sin2 + cos2 -1 = "<<(sinphi*sinphi+cosphi*cosphi-1)<<"\n";
103 
104  RotateAz(cosphi,sinphi);
105 }
G4GLOB_DLL std::ostream G4cout
static G4ThreeVector GetParticleFrameY(const G4ThreeVector &)
#define G4endl
Definition: G4ios.hh:61
Hep3Vector cross(const Hep3Vector &) const
double G4double
Definition: G4Types.hh:76
double mag() const
void RotateAz(G4ThreeVector nInteractionFrame, G4ThreeVector particleDirection)
void G4StokesVector::RotateAz ( G4double  cosphi,
G4double  sinphi 
)

Definition at line 132 of file G4StokesVector.cc.

References p1(), p2(), CLHEP::Hep3Vector::setX(), and CLHEP::Hep3Vector::setY().

133 {
134  if (!isPhoton) {
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 p2() const
void setY(double)
void setX(double)
G4double p1() const
double G4double
Definition: G4Types.hh:76
void G4StokesVector::SetPhoton ( )
inline
G4double G4StokesVector::Transverse ( ) const
inline

Definition at line 79 of file G4StokesVector.hh.

References CLHEP::Hep3Vector::perp().

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

Field Documentation

const G4StokesVector G4StokesVector::M1 =G4ThreeVector(-1.,0.,0.)
static

Definition at line 66 of file G4StokesVector.hh.

const G4StokesVector G4StokesVector::M2 =G4ThreeVector(0.,-1.,0.)
static

Definition at line 67 of file G4StokesVector.hh.

const G4StokesVector G4StokesVector::M3 =G4ThreeVector(0.,0.,-1.)
static

Definition at line 68 of file G4StokesVector.hh.

const G4StokesVector G4StokesVector::P1 =G4ThreeVector(1.,0.,0.)
static
const G4StokesVector G4StokesVector::P2 =G4ThreeVector(0.,1.,0.)
static
const G4StokesVector G4StokesVector::P3 =G4ThreeVector(0.,0.,1.)
static
const G4StokesVector G4StokesVector::ZERO =G4ThreeVector(0.,0.,0.)
static

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