Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
CLHEP::Hep3Vector Class Reference

#include <ThreeVector.h>

Inheritance diagram for CLHEP::Hep3Vector:

Public Types

enum  {
  X =0, Y =1, Z =2, NUM_COORDINATES =3,
  SIZE =NUM_COORDINATES
}
 
enum  { ToleranceTicks = 100 }
 

Public Member Functions

 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 Member Functions

static double setTolerance (double tol)
 
static double getTolerance ()
 

Protected Member Functions

void setSpherical (double r, double theta, double phi)
 
void setCylindrical (double r, double phi, double z)
 
double negativeInfinity () const
 

Protected Attributes

double dx
 
double dy
 
double dz
 

Static Protected Attributes

static DLL_API double tolerance = Hep3Vector::ToleranceTicks * 2.22045e-16
 

Detailed Description

Author

Definition at line 41 of file ThreeVector.h.

Member Enumeration Documentation

anonymous enum
Enumerator
X 
Y 
Z 
NUM_COORDINATES 
SIZE 

Definition at line 47 of file ThreeVector.h.

anonymous enum
Enumerator
ToleranceTicks 

Definition at line 298 of file ThreeVector.h.

Constructor & Destructor Documentation

CLHEP::Hep3Vector::Hep3Vector ( )
CLHEP::Hep3Vector::Hep3Vector ( double  x)
explicit
CLHEP::Hep3Vector::Hep3Vector ( double  x,
double  y 
)
CLHEP::Hep3Vector::Hep3Vector ( double  x,
double  y,
double  z 
)
CLHEP::Hep3Vector::Hep3Vector ( const Hep3Vector )
inline
CLHEP::Hep3Vector::~Hep3Vector ( )
inline

Member Function Documentation

double CLHEP::Hep3Vector::angle ( const Hep3Vector ) const

Here is the caller graph for this function:

double CLHEP::Hep3Vector::angle ( ) const
inline

Here is the caller graph for this function:

double CLHEP::Hep3Vector::azimAngle ( const Hep3Vector v2) const
double CLHEP::Hep3Vector::azimAngle ( const Hep3Vector v2,
const Hep3Vector ref 
) const

Definition at line 40 of file SpaceVectorD.cc.

41  {
42 
43  Hep3Vector vperp ( perpPart(ref) );
44  if ( vperp.mag2() == 0 ) {
45  std::cerr << "Hep3Vector::azimAngle() - "
46  << "Cannot find azimuthal angle with reference direction parallel to "
47  << "vector 1 -- will return zero" << std::endl;
48  return 0;
49  }
50 
51  Hep3Vector v2perp ( v2.perpPart(ref) );
52  if ( v2perp.mag2() == 0 ) {
53  std::cerr << "Hep3Vector::azimAngle() - "
54  << "Cannot find azimuthal angle with reference direction parallel to "
55  << "vector 2 -- will return zero" << std::endl;
56  return 0;
57  }
58 
59  double ang = vperp.angle(v2perp);
60 
61  // Now compute the sign of the answer: that of U*(VxV2) or
62  // the equivalent expression V*(V2xU).
63 
64  if ( dot(v2.cross(ref)) >= 0 ) {
65  return ang;
66  } else {
67  return -ang;
68  }
69 
70  //-| Note that if V*(V2xU) is zero, we want to return 0 or PI
71  //-| depending on whether vperp is aligned or antialigned with v2perp.
72  //-| The computed angle() expression does this properly.
73 
74 } /* azimAngle (v2, ref) */
double dot(const Hep3Vector &) const
Hep3Vector perpPart() const

Here is the call graph for this function:

double CLHEP::Hep3Vector::beta ( ) const

Definition at line 30 of file SpaceVectorP.cc.

30  {
31  double b = std::sqrt(mag2());
32 // if (b >= 1) {
33 // std::cerr << "Hep3Vector::beta() - "
34 // << "Beta taken for Hep3Vector of at least unit length" << std::endl;
35 // }
36  return b;
37 }
double mag2() const

Here is the call graph for this function:

Here is the caller graph for this function:

double CLHEP::Hep3Vector::coLinearRapidity ( ) const

Definition at line 69 of file SpaceVectorP.cc.

69  {
70  double b = beta();
71 // if (b == 1) {
72 // std::cerr << "Hep3Vector::coLinearRapidity() - "
73 // << "Co-linear Rapidity taken for Hep3Vector of unit length -- \n"
74 // << "the log should return infinity" << std::endl;
75 // }
76 // if (b > 1) {
77 // std::cerr << "Hep3Vector::coLinearRapidity() - "
78 // << "Co-linear Rapidity taken for Hep3Vector of more than unit length -- \n"
79 // << "the log would return a NAN" << std::endl;
80 // }
81  // Want inverse std::tanh(b):
82  return (.5 * std::log((1+b)/(1-b)) );
83 }
double beta() const
Definition: SpaceVectorP.cc:30

Here is the call graph for this function:

int CLHEP::Hep3Vector::compare ( const Hep3Vector v) const

Definition at line 122 of file SpaceVector.cc.

122  {
123  if ( dz > v.dz ) {
124  return 1;
125  } else if ( dz < v.dz ) {
126  return -1;
127  } else if ( dy > v.dy ) {
128  return 1;
129  } else if ( dy < v.dy ) {
130  return -1;
131  } else if ( dx > v.dx ) {
132  return 1;
133  } else if ( dx < v.dx ) {
134  return -1;
135  } else {
136  return 0;
137  }
138 } /* Compare */

Here is the caller graph for this function:

double CLHEP::Hep3Vector::cos2Theta ( ) const
inline
double CLHEP::Hep3Vector::cos2Theta ( const Hep3Vector v2) const

Definition at line 167 of file ThreeVector.cc.

167  {
168  double arg;
169  double ptot2 = mag2();
170  double qtot2 = q.mag2();
171  if ( ptot2 == 0 || qtot2 == 0 ) {
172  arg = 1.0;
173  }else{
174  double pdq = dot(q);
175  arg = (pdq/ptot2) * (pdq/qtot2);
176  // More naive methods overflow on vectors which can be squared
177  // but can't be raised to the 4th power.
178  if(arg > 1.0) arg = 1.0;
179  }
180  return arg;
181 }
double dot(const Hep3Vector &) const
double mag2() const

Here is the call graph for this function:

double CLHEP::Hep3Vector::cosTheta ( ) const
inline

Here is the caller graph for this function:

double CLHEP::Hep3Vector::cosTheta ( const Hep3Vector v2) const

Definition at line 154 of file ThreeVector.cc.

154  {
155  double arg;
156  double ptot2 = mag2()*q.mag2();
157  if(ptot2 <= 0) {
158  arg = 0.0;
159  }else{
160  arg = dot(q)/std::sqrt(ptot2);
161  if(arg > 1.0) arg = 1.0;
162  if(arg < -1.0) arg = -1.0;
163  }
164  return arg;
165 }
double dot(const Hep3Vector &) const
double mag2() const

Here is the call graph for this function:

Hep3Vector CLHEP::Hep3Vector::cross ( const Hep3Vector ) const
inline
double CLHEP::Hep3Vector::deltaPhi ( const Hep3Vector v2) const

Definition at line 138 of file ThreeVector.cc.

138  {
139  double dphi = v2.getPhi() - getPhi();
140  if ( dphi > CLHEP::pi ) {
141  dphi -= CLHEP::twopi;
142  } else if ( dphi <= -CLHEP::pi ) {
143  dphi += CLHEP::twopi;
144  }
145  return dphi;
146 } /* deltaPhi */
double getPhi() const
static constexpr double twopi
Definition: SystemOfUnits.h:55
static constexpr double pi
Definition: SystemOfUnits.h:54

Here is the call graph for this function:

Here is the caller graph for this function:

double CLHEP::Hep3Vector::deltaR ( const Hep3Vector v) const

Definition at line 148 of file ThreeVector.cc.

148  {
149  double a = eta() - v.eta();
150  double b = deltaPhi(v);
151  return std::sqrt ( a*a + b*b );
152 } /* deltaR */
double eta() const
double deltaPhi(const Hep3Vector &v2) const
Definition: ThreeVector.cc:138

Here is the call graph for this function:

double CLHEP::Hep3Vector::diff2 ( const Hep3Vector v) const
inline
double CLHEP::Hep3Vector::dot ( const Hep3Vector ) const
inline
double CLHEP::Hep3Vector::eta ( ) const

Here is the caller graph for this function:

double CLHEP::Hep3Vector::eta ( const Hep3Vector v2) const

Definition at line 117 of file SpaceVectorP.cc.

117  {
118  // Defined as -std::log ( std::tan ( .5* theta(u) ) );
119  //
120  // Quicker is to use cosTheta:
121  // std::tan (theta/2) = std::sin(theta)/(1 + std::cos(theta))
122 
123  double r1 = getR();
124  double v2r = v2.mag();
125  if ( (r1 == 0) || (v2r == 0) ) {
126  std::cerr << "Hep3Vector::eta() - "
127  << "Cannot find pseudorapidity of a zero vector relative to a vector"
128  << std::endl;
129  return 0.;
130  }
131  double c = dot(v2)/(r1*v2r);
132  if ( c >= 1 ) {
133  c = 1; //-| We don't want to return NAN because of roundoff
134  std::cerr << "Hep3Vector::eta() - "
135  << "Pseudorapidity of vector relative to parallel vector -- \n"
136  << "will give infinite result" << std::endl;
137  // We can just go on; tangent will be 0, so
138  // std::log (tangent) will be -INFINITY, so result
139  // will be +INFINITY.
140  }
141  if ( c <= -1 ) {
142  std::cerr << "Hep3Vector::eta() - "
143  << "Pseudorapidity of vector relative to anti-parallel vector -- \n"
144  << "will give negative infinite result"<< std::endl;
145  //-| We don't want to return NAN because of roundoff
146  return ( negativeInfinity() );
147  // If we just went on, the tangent would be NAN
148  // so return would be NAN. But the proper limit
149  // of tan is +Infinity, so the return should be
150  // -INFINITY.
151  }
152 
153  double tangent = std::sqrt (1-c*c) / ( 1 + c );
154  return (- std::log (tangent));
155 
156 } /* eta (u) */
double getR() const
double dot(const Hep3Vector &) const
double negativeInfinity() const
Definition: SpaceVector.cc:283

Here is the call graph for this function:

double CLHEP::Hep3Vector::gamma ( ) const

Definition at line 39 of file SpaceVectorP.cc.

39  {
40  double bbeta = std::sqrt(mag2());
41 // if (bbeta == 1) {
42 // std::cerr << "Hep3Vector::gamma() - "
43 // << "Gamma taken for Hep3Vector of unit magnitude -- infinite result"
44 // << std::endl;
45 // }
46 // if (bbeta > 1) {
47 // std::cerr << "Hep3Vector::gamma() - "
48 // << "Gamma taken for Hep3Vector of more than unit magnitude -- \n"
49 // << "the sqrt function would return NAN" << std::endl;
50 // }
51  return 1/std::sqrt(1-bbeta*bbeta);
52 }
double mag2() const

Here is the call graph for this function:

double CLHEP::Hep3Vector::getEta ( ) const
double CLHEP::Hep3Vector::getPhi ( ) const
inline

Here is the caller graph for this function:

double CLHEP::Hep3Vector::getR ( ) const
inline

Here is the caller graph for this function:

double CLHEP::Hep3Vector::getRho ( ) const
inline

Here is the caller graph for this function:

double CLHEP::Hep3Vector::getTheta ( ) const
inline

Here is the caller graph for this function:

static double CLHEP::Hep3Vector::getTolerance ( )
inlinestatic
double CLHEP::Hep3Vector::getX ( ) const
inline

Here is the caller graph for this function:

double CLHEP::Hep3Vector::getY ( ) const
inline

Here is the caller graph for this function:

double CLHEP::Hep3Vector::getZ ( ) const
inline

Here is the caller graph for this function:

double CLHEP::Hep3Vector::howNear ( const Hep3Vector v) const

Definition at line 125 of file ThreeVector.cc.

125  {
126  // | V1 - V2 | **2 / V1 dot V2, up to 1
127  double d = (*this - v).mag2();
128  double vdv = dot(v);
129  if ( (vdv > 0) && (d < vdv) ) {
130  return std::sqrt (d/vdv);
131  } else if ( (vdv == 0) && (d == 0) ) {
132  return 0;
133  } else {
134  return 1;
135  }
136 } /* howNear */
double dot(const Hep3Vector &) const
double mag2() const

Here is the call graph for this function:

double CLHEP::Hep3Vector::howOrthogonal ( const Hep3Vector v) const

Definition at line 219 of file SpaceVector.cc.

219  {
220  // | V1 dot V2 | / | V1 x V2 |
221 
222  double v1v2 = std::fabs(dot(v));
223  //-| Safe because both v1 and v2 can be squared
224  if ( v1v2 == 0 ) {
225  return 0; // Even if one or both are 0, they are considered orthogonal
226  }
227  Hep3Vector v1Xv2 ( cross(v) );
228  double abscross = v1Xv2.mag();
229  if ( v1v2 >= abscross ) {
230  return 1;
231  } else {
232  return v1v2/abscross;
233  }
234 
235 } /* howOrthogonal() */
double dot(const Hep3Vector &) const
Hep3Vector cross(const Hep3Vector &) const

Here is the call graph for this function:

Here is the caller graph for this function:

double CLHEP::Hep3Vector::howParallel ( const Hep3Vector v) const

Definition at line 168 of file SpaceVector.cc.

168  {
169  // | V1 x V2 | / | V1 dot V2 |
170  double v1v2 = std::fabs(dot(v));
171  if ( v1v2 == 0 ) {
172  // Zero is parallel to no other vector except for zero.
173  return ( (mag2() == 0) && (v.mag2() == 0) ) ? 0 : 1;
174  }
175  Hep3Vector v1Xv2 ( cross(v) );
176  double abscross = v1Xv2.mag();
177  if ( abscross >= v1v2 ) {
178  return 1;
179  } else {
180  return abscross/v1v2;
181  }
182 } /* howParallel() */
double dot(const Hep3Vector &) const
double mag2() const
Hep3Vector cross(const Hep3Vector &) const

Here is the call graph for this function:

bool CLHEP::Hep3Vector::isNear ( const Hep3Vector v,
double  epsilon = tolerance 
) const

Definition at line 120 of file ThreeVector.cc.

120  {
121  double limit = dot(v)*epsilon*epsilon;
122  return ( (*this - v).mag2() <= limit );
123 } /* isNear() */
double dot(const Hep3Vector &) const
double mag2() const
double epsilon(double density, double temperature)

Here is the call graph for this function:

bool CLHEP::Hep3Vector::isOrthogonal ( const Hep3Vector v,
double  epsilon = tolerance 
) const

Definition at line 237 of file SpaceVector.cc.

238  {
239 // | V1 x V2 | **2 <= epsilon **2 | V1 dot V2 | **2
240 // V1 is *this, V2 is v
241 
242  static const double TOOBIG = std::pow(2.0,507);
243  static const double SCALE = std::pow(2.0,-507);
244  double v1v2 = std::fabs(dot(v));
245  //-| Safe because both v1 and v2 can be squared
246  if ( v1v2 >= TOOBIG ) {
247  Hep3Vector sv1 ( *this * SCALE );
248  Hep3Vector sv2 ( v * SCALE );
249  Hep3Vector sv1Xsv2 = sv1.cross(sv2);
250  double x2 = sv1Xsv2.mag2();
251  double limit = epsilon*epsilon*x2;
252  double y2 = v1v2*SCALE*SCALE;
253  return ( y2*y2 <= limit );
254  }
255 
256  // At this point we know v1v2 can be squared.
257 
258  Hep3Vector eps_v1Xv2 ( cross(epsilon*v) );
259  if ( (std::fabs (eps_v1Xv2.dx) > TOOBIG) ||
260  (std::fabs (eps_v1Xv2.dy) > TOOBIG) ||
261  (std::fabs (eps_v1Xv2.dz) > TOOBIG) ) {
262  return true;
263  }
264 
265  // At this point we know all the math we need can be done.
266 
267  return ( v1v2*v1v2 <= eps_v1Xv2.mag2() );
268 
269 } /* isOrthogonal() */
double dot(const Hep3Vector &) const
Hep3Vector cross(const Hep3Vector &) const
double epsilon(double density, double temperature)

Here is the call graph for this function:

Here is the caller graph for this function:

bool CLHEP::Hep3Vector::isParallel ( const Hep3Vector v,
double  epsilon = tolerance 
) const

Definition at line 184 of file SpaceVector.cc.

185  {
186  // | V1 x V2 | **2 <= epsilon **2 | V1 dot V2 | **2
187  // V1 is *this, V2 is v
188 
189  static const double TOOBIG = std::pow(2.0,507);
190  static const double SCALE = std::pow(2.0,-507);
191  double v1v2 = std::fabs(dot(v));
192  if ( v1v2 == 0 ) {
193  return ( (mag2() == 0) && (v.mag2() == 0) );
194  }
195  if ( v1v2 >= TOOBIG ) {
196  Hep3Vector sv1 ( *this * SCALE );
197  Hep3Vector sv2 ( v * SCALE );
198  Hep3Vector sv1Xsv2 = sv1.cross(sv2);
199  double x2 = sv1Xsv2.mag2();
200  double limit = v1v2*SCALE*SCALE;
201  limit = epsilon*epsilon*limit*limit;
202  return ( x2 <= limit );
203  }
204 
205  // At this point we know v1v2 can be squared.
206 
207  Hep3Vector v1Xv2 ( cross(v) );
208  if ( (std::fabs (v1Xv2.dx) > TOOBIG) ||
209  (std::fabs (v1Xv2.dy) > TOOBIG) ||
210  (std::fabs (v1Xv2.dz) > TOOBIG) ) {
211  return false;
212  }
213 
214  return ( (v1Xv2.mag2()) <= ((epsilon * v1v2) * (epsilon * v1v2)) );
215 
216 } /* isParallel() */
double dot(const Hep3Vector &) const
double mag2() const
Hep3Vector cross(const Hep3Vector &) const
double epsilon(double density, double temperature)

Here is the call graph for this function:

double CLHEP::Hep3Vector::mag ( ) const
inline
double CLHEP::Hep3Vector::mag2 ( ) const
inline
double CLHEP::Hep3Vector::negativeInfinity ( ) const
protected

Definition at line 283 of file SpaceVector.cc.

283  {
284  // A byte-order-independent way to return -Infinity
285  struct Dib {
286  union {
287  double d;
288  unsigned char i[8];
289  } u;
290  };
291  Dib negOne;
292  Dib posTwo;
293  negOne.u.d = -1.0;
294  posTwo.u.d = 2.0;
295  Dib value;
296  int k;
297  for (k=0; k<8; k++) {
298  value.u.i[k] = negOne.u.i[k] | posTwo.u.i[k];
299  }
300  return value.u.d;
301 }
const XML_Char int const XML_Char * value
Definition: expat.h:331

Here is the caller graph for this function:

bool CLHEP::Hep3Vector::operator!= ( const Hep3Vector ) const
inline
double CLHEP::Hep3Vector::operator() ( int  ) const
inline
double& CLHEP::Hep3Vector::operator() ( int  )
inline
Hep3Vector& CLHEP::Hep3Vector::operator*= ( double  )
inline

Here is the caller graph for this function:

Hep3Vector & CLHEP::Hep3Vector::operator*= ( const HepRotation m1)

Definition at line 20 of file ThreeVectorR.cc.

20  {
21  return *this = m1 * (*this);
22 }
Hep3Vector& CLHEP::Hep3Vector::operator+= ( const Hep3Vector )
inline
Hep3Vector CLHEP::Hep3Vector::operator- ( ) const
inline
Hep3Vector& CLHEP::Hep3Vector::operator-= ( const Hep3Vector )
inline
Hep3Vector & CLHEP::Hep3Vector::operator/= ( double  c)

Definition at line 313 of file ThreeVector.cc.

313  {
314 // if (c == 0) {
315 // std::cerr << "Hep3Vector::operator/ () - "
316 // << "Attempt to do vector /= 0 -- "
317 // << "division by zero would produce infinite or NAN components"
318 // << std::endl;
319 // }
320  double oneOverC = 1.0/c;
321  dx *= oneOverC;
322  dy *= oneOverC;
323  dz *= oneOverC;
324  return *this;
325 }
bool CLHEP::Hep3Vector::operator< ( const Hep3Vector v) const

Definition at line 144 of file SpaceVector.cc.

144  {
145  return (compare(v) < 0);
146 }
int compare(const Hep3Vector &v) const
Definition: SpaceVector.cc:122

Here is the call graph for this function:

bool CLHEP::Hep3Vector::operator<= ( const Hep3Vector v) const

Definition at line 150 of file SpaceVector.cc.

150  {
151  return (compare(v) <= 0);
152 }
int compare(const Hep3Vector &v) const
Definition: SpaceVector.cc:122

Here is the call graph for this function:

Hep3Vector& CLHEP::Hep3Vector::operator= ( const Hep3Vector )
inline
bool CLHEP::Hep3Vector::operator== ( const Hep3Vector ) const
inline
bool CLHEP::Hep3Vector::operator> ( const Hep3Vector v) const

Definition at line 141 of file SpaceVector.cc.

141  {
142  return (compare(v) > 0);
143 }
int compare(const Hep3Vector &v) const
Definition: SpaceVector.cc:122

Here is the call graph for this function:

bool CLHEP::Hep3Vector::operator>= ( const Hep3Vector v) const

Definition at line 147 of file SpaceVector.cc.

147  {
148  return (compare(v) >= 0);
149 }
int compare(const Hep3Vector &v) const
Definition: SpaceVector.cc:122

Here is the call graph for this function:

double CLHEP::Hep3Vector::operator[] ( int  ) const
inline
double& CLHEP::Hep3Vector::operator[] ( int  )
inline
Hep3Vector CLHEP::Hep3Vector::orthogonal ( ) const
inline

Here is the caller graph for this function:

double CLHEP::Hep3Vector::perp ( ) const
inline

Here is the caller graph for this function:

double CLHEP::Hep3Vector::perp ( const Hep3Vector ) const
inline
double CLHEP::Hep3Vector::perp2 ( ) const
inline

Here is the caller graph for this function:

double CLHEP::Hep3Vector::perp2 ( const Hep3Vector ) const
inline
Hep3Vector CLHEP::Hep3Vector::perpPart ( ) const
inline

Here is the caller graph for this function:

Hep3Vector CLHEP::Hep3Vector::perpPart ( const Hep3Vector v2) const
inline
double CLHEP::Hep3Vector::phi ( ) const
inline

Here is the caller graph for this function:

double CLHEP::Hep3Vector::polarAngle ( const Hep3Vector v2) const

Definition at line 28 of file SpaceVectorD.cc.

28  {
29  return std::fabs(v2.getTheta() - getTheta());
30 } /* polarAngle */
double getTheta() const

Here is the call graph for this function:

double CLHEP::Hep3Vector::polarAngle ( const Hep3Vector v2,
const Hep3Vector ref 
) const

Definition at line 32 of file SpaceVectorD.cc.

33  {
34  return std::fabs( v2.angle(ref) - angle(ref) );
35 } /* polarAngle (v2, ref) */
double angle() const

Here is the call graph for this function:

Hep3Vector CLHEP::Hep3Vector::project ( ) const
inline

Here is the caller graph for this function:

Hep3Vector CLHEP::Hep3Vector::project ( const Hep3Vector v2) const

Definition at line 89 of file SpaceVectorP.cc.

89  {
90  double mag2v2 = v2.mag2();
91  if (mag2v2 == 0) {
92  std::cerr << "Hep3Vector::project() - "
93  << "Attempt to take projection of vector against zero reference vector"
94  << std::endl;
95  return project();
96  }
97  return ( v2 * (dot(v2)/mag2v2) );
98 }
double dot(const Hep3Vector &) const
Hep3Vector project() const

Here is the call graph for this function:

double CLHEP::Hep3Vector::pseudoRapidity ( ) const

Definition at line 58 of file ThreeVector.cc.

58  {
59  double m1 = mag();
60  if ( m1== 0 ) return 0.0;
61  if ( m1== z() ) return 1.0E72;
62  if ( m1== -z() ) return -1.0E72;
63  return 0.5*std::log( (m1+z())/(m1-z()) );
64 }
double z() const
double mag() const

Here is the call graph for this function:

double CLHEP::Hep3Vector::r ( ) const
inline
double CLHEP::Hep3Vector::rapidity ( ) const

Definition at line 54 of file SpaceVectorP.cc.

54  {
55 // if (std::fabs(dz) == 1) {
56 // std::cerr << "Hep3Vector::rapidity() - "
57 // << "Rapidity in Z direction taken for Hep3Vector with |Z| = 1 -- \n"
58 // << "the log should return infinity" <, std::endl;
59 // }
60 // if (std::fabs(dz) > 1) {
61 // std::cerr << "Hep3Vector::rapidity() - "
62 // << "Rapidity in Z direction taken for Hep3Vector with |Z| > 1 -- \n"
63 // << "the log would return a NAN" << std::endl;
64 // }
65  // Want inverse std::tanh(dz):
66  return (.5 * std::log((1+dz)/(1-dz)) );
67 }
double CLHEP::Hep3Vector::rapidity ( const Hep3Vector v2) const

Definition at line 100 of file SpaceVectorP.cc.

100  {
101  double vmag = v2.mag();
102  if ( vmag == 0 ) {
103  std::cerr << "Hep3Vector::rapidity() - "
104  << "Rapidity taken with respect to zero vector" << std::endl;
105  return 0;
106  }
107  double z1 = dot(v2)/vmag;
108 // if (std::fabs(z1) >= 1) {
109 // std::cerr << "Hep3Vector::rapidity() - "
110 // << "Rapidity taken for too large a Hep3Vector "
111 // << "-- would return infinity or NAN" << std::endl;
112 // }
113  // Want inverse std::tanh(z):
114  return (.5 * std::log((1+z1)/(1-z1)) );
115 }
double dot(const Hep3Vector &) const

Here is the call graph for this function:

double CLHEP::Hep3Vector::rho ( ) const
inline
Hep3Vector & CLHEP::Hep3Vector::rotate ( double  angle1,
const Hep3Vector aaxis 
)

Definition at line 28 of file ThreeVectorR.cc.

28  {
29  HepRotation trans;
30  trans.rotate(angle1, aaxis);
31  operator*=(trans);
32  return *this;
33 }
Hep3Vector & operator*=(double)

Here is the call graph for this function:

Here is the caller graph for this function:

Hep3Vector & CLHEP::Hep3Vector::rotate ( const Hep3Vector axis,
double  delta 
)

Definition at line 25 of file SpaceVectorR.cc.

26  {
27  double r1 = axis.mag();
28  if ( r1 == 0 ) {
29  std::cerr << "Hep3Vector::rotate() - "
30  << "Attempt to rotate around a zero vector axis! " << std::endl;
31  return *this;
32  }
33  double scale=1.0/r1;
34  double ux = scale*axis.getX();
35  double uy = scale*axis.getY();
36  double uz = scale*axis.getZ();
37  double cd = std::cos(ddelta);
38  double sd = std::sin(ddelta);
39  double ocd = 1 - cd;
40  double rx;
41  double ry;
42  double rz;
43 
44  { double ocdux = ocd * ux;
45  rx = dx * ( cd + ocdux * ux ) +
46  dy * ( ocdux * uy - sd * uz ) +
47  dz * ( ocdux * uz + sd * uy ) ;
48  }
49 
50  { double ocduy = ocd * uy;
51  ry = dy * ( cd + ocduy * uy ) +
52  dz * ( ocduy * uz - sd * ux ) +
53  dx * ( ocduy * ux + sd * uz ) ;
54  }
55 
56  { double ocduz = ocd * uz;
57  rz = dz * ( cd + ocduz * uz ) +
58  dx * ( ocduz * ux - sd * uy ) +
59  dy * ( ocduz * uy + sd * ux ) ;
60  }
61 
62  dx = rx;
63  dy = ry;
64  dz = rz;
65 
66  return *this;
67 } /* rotate */
static const G4double cd

Here is the call graph for this function:

Hep3Vector & CLHEP::Hep3Vector::rotate ( const HepAxisAngle ax)

Definition at line 112 of file SpaceVectorR.cc.

112  {
113  return rotate( ax.getAxis(), ax.delta() );
114 }
Hep3Vector & rotate(double, const Hep3Vector &)
Definition: ThreeVectorR.cc:28

Here is the call graph for this function:

Hep3Vector & CLHEP::Hep3Vector::rotate ( const HepEulerAngles e)

Definition at line 116 of file SpaceVectorR.cc.

116  {
117  return rotate( ex.phi(), ex.theta(), ex.psi() );
118 }
Hep3Vector & rotate(double, const Hep3Vector &)
Definition: ThreeVectorR.cc:28

Here is the call graph for this function:

Hep3Vector & CLHEP::Hep3Vector::rotate ( double  phi,
double  theta,
double  psi 
)

Definition at line 74 of file SpaceVectorR.cc.

76  {
77 
78  double rx;
79  double ry;
80  double rz;
81 
82  double sinPhi = std::sin( phi1 ), cosPhi = std::cos( phi1 );
83  double sinTheta = std::sin( theta1 ), cosTheta1 = std::cos( theta1 );
84  double sinPsi = std::sin( psi1 ), cosPsi = std::cos( psi1 );
85 
86  rx = (cosPsi * cosPhi - cosTheta1 * sinPsi * sinPhi) * dx +
87  (cosPsi * sinPhi + cosTheta1 * sinPsi * cosPhi) * dy +
88  (sinPsi * sinTheta) * dz ;
89 
90  ry = (- sinPsi * cosPhi - cosTheta1 * cosPsi * sinPhi) * dx +
91  (- sinPsi * sinPhi + cosTheta1 * cosPsi * cosPhi) * dy +
92  (cosPsi * sinTheta) * dz ;
93 
94  rz = (sinTheta * sinPhi) * dx +
95  (- sinTheta * cosPhi) * dy +
96  (cosTheta1) * dz ;
97 
98  dx = rx;
99  dy = ry;
100  dz = rz;
101 
102  return *this;
103 
104 } /* rotate */
Hep3Vector & CLHEP::Hep3Vector::rotateUz ( const Hep3Vector NewUzVector)

Definition at line 38 of file ThreeVector.cc.

38  {
39  // NewUzVector must be normalized !
40 
41  double u1 = NewUzVector.x();
42  double u2 = NewUzVector.y();
43  double u3 = NewUzVector.z();
44  double up = u1*u1 + u2*u2;
45 
46  if (up>0) {
47  up = std::sqrt(up);
48  double px = dx, py = dy, pz = dz;
49  dx = (u1*u3*px - u2*py)/up + u1*pz;
50  dy = (u2*u3*px + u1*py)/up + u2*pz;
51  dz = -up*px + u3*pz;
52  }
53  else if (u3 < 0.) { dx = -dx; dz = -dz; } // phi=0 teta=pi
54  else {};
55  return *this;
56 }

Here is the call graph for this function:

Hep3Vector & CLHEP::Hep3Vector::rotateX ( double  phi1)

Definition at line 90 of file ThreeVector.cc.

90  {
91  double sinphi = std::sin(phi1);
92  double cosphi = std::cos(phi1);
93  double ty;
94  ty = dy * cosphi - dz * sinphi;
95  dz = dz * cosphi + dy * sinphi;
96  dy = ty;
97  return *this;
98 } /* rotateX */

Here is the caller graph for this function:

Hep3Vector & CLHEP::Hep3Vector::rotateY ( double  phi1)

Definition at line 100 of file ThreeVector.cc.

100  {
101  double sinphi = std::sin(phi1);
102  double cosphi = std::cos(phi1);
103  double tz;
104  tz = dz * cosphi - dx * sinphi;
105  dx = dx * cosphi + dz * sinphi;
106  dz = tz;
107  return *this;
108 } /* rotateY */

Here is the caller graph for this function:

Hep3Vector & CLHEP::Hep3Vector::rotateZ ( double  phi1)

Definition at line 110 of file ThreeVector.cc.

110  {
111  double sinphi = std::sin(phi1);
112  double cosphi = std::cos(phi1);
113  double tx;
114  tx = dx * cosphi - dy * sinphi;
115  dy = dy * cosphi + dx * sinphi;
116  dx = tx;
117  return *this;
118 } /* rotateZ */

Here is the caller graph for this function:

void CLHEP::Hep3Vector::set ( double  x,
double  y,
double  z 
)
inline
void CLHEP::Hep3Vector::setCylEta ( double  p)

Definition at line 259 of file ThreeVector.cc.

259  {
260 
261  // In cylindrical coords, set eta while keeping rho and phi fixed
262 
263  double theta1 = 2 * std::atan ( std::exp (-eta1) );
264 
265  //-| The remaining code is similar to setCylTheta, The reason for
266  //-| using a copy is so as to be able to change the messages in the
267  //-| ZMthrows to say eta rather than theta. Besides, we assumedly
268  //-| need not check for theta of 0 or PI.
269 
270  if ( (dx == 0) && (dy == 0) ) {
271  if (dz == 0) {
272  std::cerr << "Hep3Vector::setCylEta() - "
273  << "Attempt to set cylEta of zero vector -- vector is unchanged"
274  << std::endl;
275  return;
276  }
277  if (theta1 == 0) {
278  dz = std::fabs(dz);
279  return;
280  }
281  if (theta1 == CLHEP::pi) {
282  dz = -std::fabs(dz);
283  return;
284  }
285  std::cerr << "Hep3Vector::setCylEta() - "
286  << "Attempt set cylindrical eta of vector along Z axis "
287  << "to a non-trivial value, while keeping rho fixed -- "
288  << "will return zero vector" << std::endl;
289  dz = 0;
290  return;
291  }
292  double phi1 (getPhi());
293  double rho1 = getRho();
294  dz = rho1 / std::tan (theta1);
295  dy = rho1 * std::sin (phi1);
296  dx = rho1 * std::cos (phi1);
297 
298 } /* setCylEta */
double getRho() const
double getPhi() const
static constexpr double pi
Definition: SystemOfUnits.h:54

Here is the call graph for this function:

void CLHEP::Hep3Vector::setCylindrical ( double  r,
double  phi,
double  z 
)
protected

Definition at line 56 of file SpaceVector.cc.

59  {
60 // if ( rho1 < 0 ) {
61 // std::cerr << "Hep3Vector::setCylindrical() - "
62 // << "Cylindrical coordinates supplied with negative Rho" << std::endl;
63 // // No special return needed if warning is ignored.
64 // }
65  dz = z1;
66  dy = rho1 * std::sin (phi1);
67  dx = rho1 * std::cos (phi1);
68  return;
69 } /* setCylindrical (r, phi, z) */
void CLHEP::Hep3Vector::setCylTheta ( double  theta1)

Definition at line 211 of file ThreeVector.cc.

211  {
212 
213  // In cylindrical coords, set theta while keeping rho and phi fixed
214 
215  if ( (dx == 0) && (dy == 0) ) {
216  if (dz == 0) {
217  std::cerr << "Hep3Vector::setCylTheta() - "
218  << "Attempt to set cylTheta of zero vector -- vector is unchanged"
219  << std::endl;
220  return;
221  }
222  if (theta1 == 0) {
223  dz = std::fabs(dz);
224  return;
225  }
226  if (theta1 == CLHEP::pi) {
227  dz = -std::fabs(dz);
228  return;
229  }
230  std::cerr << "Hep3Vector::setCylTheta() - "
231  << "Attempt set cylindrical theta of vector along Z axis "
232  << "to a non-trivial value, while keeping rho fixed -- "
233  << "will return zero vector" << std::endl;
234  dz = 0;
235  return;
236  }
237  if ( (theta1 < 0) || (theta1 > CLHEP::pi) ) {
238  std::cerr << "Hep3Vector::setCylTheta() - "
239  << "Setting Cyl theta of a vector based on a value not in [0, PI]"
240  << std::endl;
241  // No special return needed if warning is ignored.
242  }
243  double phi1 (getPhi());
244  double rho1 = getRho();
245  if ( (theta1 == 0) || (theta1 == CLHEP::pi) ) {
246  std::cerr << "Hep3Vector::setCylTheta() - "
247  << "Attempt to set cylindrical theta to 0 or PI "
248  << "while keeping rho fixed -- infinite Z will be computed"
249  << std::endl;
250  dz = (theta1==0) ? 1.0E72 : -1.0E72;
251  return;
252  }
253  dz = rho1 / std::tan (theta1);
254  dy = rho1 * std::sin (phi1);
255  dx = rho1 * std::cos (phi1);
256 
257 } /* setCylTheta */
double getRho() const
double getPhi() const
static constexpr double pi
Definition: SystemOfUnits.h:54

Here is the call graph for this function:

void CLHEP::Hep3Vector::setEta ( double  p)

Definition at line 183 of file ThreeVector.cc.

183  {
184  double phi1 = 0;
185  double r1;
186  if ( (dx == 0) && (dy == 0) ) {
187  if (dz == 0) {
188  std::cerr << "Hep3Vector::setEta() - "
189  << "Attempt to set eta of zero vector -- vector is unchanged"
190  << std::endl;
191  return;
192  }
193  std::cerr << "Hep3Vector::setEta() - "
194  << "Attempt to set eta of vector along Z axis -- will use phi = 0"
195  << std::endl;
196  r1 = std::fabs(dz);
197  } else {
198  r1 = getR();
199  phi1 = getPhi();
200  }
201  double tanHalfTheta = std::exp ( -eta1 );
202  double cosTheta1 =
203  (1 - tanHalfTheta*tanHalfTheta) / (1 + tanHalfTheta*tanHalfTheta);
204  dz = r1 * cosTheta1;
205  double rho1 = r1*std::sqrt(1 - cosTheta1*cosTheta1);
206  dy = rho1 * std::sin (phi1);
207  dx = rho1 * std::cos (phi1);
208  return;
209 }
double getR() const
double getPhi() const

Here is the call graph for this function:

void CLHEP::Hep3Vector::setMag ( double  ma)

Definition at line 25 of file ThreeVector.cc.

25  {
26  double factor = mag();
27  if (factor == 0) {
28  std::cerr << "Hep3Vector::setMag() - "
29  << "zero vector can't be stretched" << std::endl;
30  }else{
31  factor = ma/factor;
32  setX(x()*factor);
33  setY(y()*factor);
34  setZ(z()*factor);
35  }
36 }
double x() const
void setY(double)
double z() const
void setZ(double)
void setX(double)
double y() const
double mag() const

Here is the call graph for this function:

Here is the caller graph for this function:

void CLHEP::Hep3Vector::setPerp ( double  )
inline
void CLHEP::Hep3Vector::setPhi ( double  )
inline

Here is the caller graph for this function:

void CLHEP::Hep3Vector::setR ( double  s)
inline

Here is the caller graph for this function:

void CLHEP::Hep3Vector::setREtaPhi ( double  r,
double  eta,
double  phi 
)
inline
void CLHEP::Hep3Vector::setRho ( double  s)
inline
void CLHEP::Hep3Vector::setRhoPhiEta ( double  rho,
double  phi,
double  eta 
)

Definition at line 98 of file SpaceVector.cc.

101  {
102  if (rho1 == 0) {
103  std::cerr << "Hep3Vector::setRhoPhiEta() - "
104  << "Attempt set vector components rho, phi, eta with zero rho -- "
105  << "zero vector is returned, ignoring eta and phi" << std::endl;
106  dx = 0; dy = 0; dz = 0;
107  return;
108  }
109  double theta1 (2 * std::atan ( std::exp (-eta1) ));
110  dz = rho1 / std::tan (theta1);
111  dy = rho1 * std::sin (phi1);
112  dx = rho1 * std::cos (phi1);
113  return;
114 } /* setCyl (rho, phi, eta) */
void CLHEP::Hep3Vector::setRhoPhiTheta ( double  rho,
double  phi,
double  theta 
)

Definition at line 71 of file SpaceVector.cc.

74  {
75  if (rho1 == 0) {
76  std::cerr << "Hep3Vector::setRhoPhiTheta() - "
77  << "Attempt set vector components rho, phi, theta with zero rho -- "
78  << "zero vector is returned, ignoring theta and phi" << std::endl;
79  dx = 0; dy = 0; dz = 0;
80  return;
81  }
82 // if ( (theta1 == 0) || (theta1 == CLHEP::pi) ) {
83 // std::cerr << "Hep3Vector::setRhoPhiTheta() - "
84 // << "Attempt set cylindrical vector vector with finite rho and "
85 // << "theta along the Z axis: infinite Z would be computed" << std::endl;
86 // }
87 // if ( (theta1 < 0) || (theta1 > CLHEP::pi) ) {
88 // std::cerr << "Hep3Vector::setRhoPhiTheta() - "
89 // << "Rho, phi, theta set with theta not in [0, PI]" << std::endl;
90 // // No special return needed if warning is ignored.
91 // }
92  dz = rho1 / std::tan (theta1);
93  dy = rho1 * std::sin (phi1);
94  dx = rho1 * std::cos (phi1);
95  return;
96 } /* setCyl (rho, phi, theta) */
void CLHEP::Hep3Vector::setRhoPhiZ ( double  rho,
double  phi,
double  z 
)
inline
void CLHEP::Hep3Vector::setRThetaPhi ( double  r,
double  theta,
double  phi 
)
inline

Here is the caller graph for this function:

void CLHEP::Hep3Vector::setSpherical ( double  r,
double  theta,
double  phi 
)
protected

Definition at line 35 of file SpaceVector.cc.

38  {
39 // if ( r1 < 0 ) {
40 // std::cerr << "Hep3Vector::setSpherical() - "
41 // << "Spherical coordinates set with negative R" << std::endl;
42 // // No special return needed if warning is ignored.
43 // }
44 // if ( (theta1 < 0) || (theta1 > CLHEP::pi) ) {
45 // std::cerr << "Hep3Vector::setSpherical() - "
46 // << "Spherical coordinates set with theta not in [0, PI]" << std::endl;
47 // // No special return needed if warning is ignored.
48 // }
49  dz = r1 * std::cos(theta1);
50  double rho1 ( r1*std::sin(theta1));
51  dy = rho1 * std::sin (phi1);
52  dx = rho1 * std::cos (phi1);
53  return;
54 } /* setSpherical (r, theta1, phi1) */
void CLHEP::Hep3Vector::setTheta ( double  )
inline

Here is the caller graph for this function:

double CLHEP::Hep3Vector::setTolerance ( double  tol)
static

Definition at line 271 of file SpaceVector.cc.

271  {
272 // Set the tolerance for Hep3Vectors to be considered near one another
273  double oldTolerance (tolerance);
274  tolerance = tol;
275  return oldTolerance;
276 }
static DLL_API double tolerance
Definition: ThreeVector.h:399
void CLHEP::Hep3Vector::setX ( double  )
inline
void CLHEP::Hep3Vector::setY ( double  )
inline
void CLHEP::Hep3Vector::setZ ( double  )
inline
double CLHEP::Hep3Vector::theta ( ) const
inline

Here is the caller graph for this function:

double CLHEP::Hep3Vector::theta ( const Hep3Vector v2) const
inline
Hep3Vector & CLHEP::Hep3Vector::transform ( const HepRotation m1)

Definition at line 24 of file ThreeVectorR.cc.

24  {
25  return *this = m1 * (*this);
26 }

Here is the caller graph for this function:

Hep3Vector CLHEP::Hep3Vector::unit ( ) const
inline
double CLHEP::Hep3Vector::x ( ) const
inline
double CLHEP::Hep3Vector::y ( ) const
inline
double CLHEP::Hep3Vector::z ( ) const
inline

Member Data Documentation

double CLHEP::Hep3Vector::dx
protected

Definition at line 394 of file ThreeVector.h.

double CLHEP::Hep3Vector::dy
protected

Definition at line 395 of file ThreeVector.h.

double CLHEP::Hep3Vector::dz
protected

Definition at line 396 of file ThreeVector.h.

double CLHEP::Hep3Vector::tolerance = Hep3Vector::ToleranceTicks * 2.22045e-16
staticprotected

Definition at line 399 of file ThreeVector.h.


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