CLHEP VERSION Reference Documentation
   
CLHEP Home Page     CLHEP Documentation     CLHEP Bug Reports

Vector/Vector/ThreeVector.h
Go to the documentation of this file.
1 // -*- C++ -*-
2 // CLASSDOC OFF
3 // $Id: ThreeVector.h,v 1.4 2010/06/16 17:15:57 garren Exp $
4 // ---------------------------------------------------------------------------
5 // CLASSDOC ON
6 //
7 // This file is a part of the CLHEP - a Class Library for High Energy Physics.
8 //
9 // Hep3Vector is a general 3-vector class defining vectors in three
10 // dimension using double components. Rotations of these vectors are
11 // performed by multiplying with an object of the HepRotation class.
12 //
13 // .SS See Also
14 // LorentzVector.h, Rotation.h, LorentzRotation.h
15 //
16 // .SS Authors
17 // Leif Lonnblad and Anders Nilsson; Modified by Evgueni Tcherniaev;
18 // ZOOM additions by Mark Fischler
19 //
20 
21 #ifndef HEP_THREEVECTOR_H
22 #define HEP_THREEVECTOR_H
23 
24 #ifdef GNUPRAGMA
25 #pragma interface
26 #endif
27 
28 #include <iostream>
29 #include "CLHEP/Vector/defs.h"
30 
31 namespace CLHEP {
32 
33 class HepRotation;
34 class HepEulerAngles;
35 class HepAxisAngle;
36 
41 class Hep3Vector {
42 
43 public:
44 
45 // Basic properties and operations on 3-vectors:
46 
47  enum { X=0, Y=1, Z=2, NUM_COORDINATES=3, SIZE=NUM_COORDINATES };
48  // Safe indexing of the coordinates when using with matrices, arrays, etc.
49  // (BaBar)
50 
51  Hep3Vector();
52  explicit Hep3Vector(double x);
53  Hep3Vector(double x, double y);
54  Hep3Vector(double x, double y, double z);
55  // The constructor.
56 
57  inline Hep3Vector(const Hep3Vector &);
58  // The copy constructor.
59 
60  inline ~Hep3Vector();
61  // The destructor. Not virtual - inheritance from this class is dangerous.
62 
63  double operator () (int) const;
64  // Get components by index -- 0-based (Geant4)
65 
66  inline double operator [] (int) const;
67  // Get components by index -- 0-based (Geant4)
68 
69  double & operator () (int);
70  // Set components by index. 0-based.
71 
72  inline double & operator [] (int);
73  // Set components by index. 0-based.
74 
75  inline double x() const;
76  inline double y() const;
77  inline double z() const;
78  // The components in cartesian coordinate system. Same as getX() etc.
79 
80  inline void setX(double);
81  inline void setY(double);
82  inline void setZ(double);
83  // Set the components in cartesian coordinate system.
84 
85  inline void set( double x, double y, double z);
86  // Set all three components in cartesian coordinate system.
87 
88  inline double phi() const;
89  // The azimuth angle.
90 
91  inline double theta() const;
92  // The polar angle.
93 
94  inline double cosTheta() const;
95  // Cosine of the polar angle.
96 
97  inline double cos2Theta() const;
98  // Cosine squared of the polar angle - faster than cosTheta(). (ZOOM)
99 
100  inline double mag2() const;
101  // The magnitude squared (r^2 in spherical coordinate system).
102 
103  inline double mag() const;
104  // The magnitude (r in spherical coordinate system).
105 
106  inline void setPhi(double);
107  // Set phi keeping mag and theta constant (BaBar).
108 
109  inline void setTheta(double);
110  // Set theta keeping mag and phi constant (BaBar).
111 
112  void setMag(double);
113  // Set magnitude keeping theta and phi constant (BaBar).
114 
115  inline double perp2() const;
116  // The transverse component squared (rho^2 in cylindrical coordinate system).
117 
118  inline double perp() const;
119  // The transverse component (rho in cylindrical coordinate system).
120 
121  inline void setPerp(double);
122  // Set the transverse component keeping phi and z constant.
123 
124  void setCylTheta(double);
125  // Set theta while keeping transvers component and phi fixed
126 
127  inline double perp2(const Hep3Vector &) const;
128  // The transverse component w.r.t. given axis squared.
129 
130  inline double perp(const Hep3Vector &) const;
131  // The transverse component w.r.t. given axis.
132 
133  inline Hep3Vector & operator = (const Hep3Vector &);
134  // Assignment.
135 
136  inline bool operator == (const Hep3Vector &) const;
137  inline bool operator != (const Hep3Vector &) const;
138  // Comparisons (Geant4).
139 
140  bool isNear (const Hep3Vector &, double epsilon=tolerance) const;
141  // Check for equality within RELATIVE tolerance (default 2.2E-14). (ZOOM)
142  // |v1 - v2|**2 <= epsilon**2 * |v1.dot(v2)|
143 
144  double howNear(const Hep3Vector & v ) const;
145  // sqrt ( |v1-v2|**2 / v1.dot(v2) ) with a maximum of 1.
146  // If v1.dot(v2) is negative, will return 1.
147 
148  double deltaR(const Hep3Vector & v) const;
149  // sqrt( pseudorapity_difference**2 + deltaPhi **2 )
150 
151  inline Hep3Vector & operator += (const Hep3Vector &);
152  // Addition.
153 
154  inline Hep3Vector & operator -= (const Hep3Vector &);
155  // Subtraction.
156 
157  inline Hep3Vector operator - () const;
158  // Unary minus.
159 
160  inline Hep3Vector & operator *= (double);
161  // Scaling with real numbers.
162 
163  Hep3Vector & operator /= (double);
164  // Division by (non-zero) real number.
165 
166  inline Hep3Vector unit() const;
167  // Vector parallel to this, but of length 1.
168 
169  inline Hep3Vector orthogonal() const;
170  // Vector orthogonal to this (Geant4).
171 
172  inline double dot(const Hep3Vector &) const;
173  // double product.
174 
175  inline Hep3Vector cross(const Hep3Vector &) const;
176  // Cross product.
177 
178  double angle(const Hep3Vector &) const;
179  // The angle w.r.t. another 3-vector.
180 
181  double pseudoRapidity() const;
182  // Returns the pseudo-rapidity, i.e. -ln(tan(theta/2))
183 
184  void setEta ( double p );
185  // Set pseudo-rapidity, keeping magnitude and phi fixed. (ZOOM)
186 
187  void setCylEta ( double p );
188  // Set pseudo-rapidity, keeping transverse component and phi fixed. (ZOOM)
189 
190  Hep3Vector & rotateX(double);
191  // Rotates the Hep3Vector around the x-axis.
192 
193  Hep3Vector & rotateY(double);
194  // Rotates the Hep3Vector around the y-axis.
195 
196  Hep3Vector & rotateZ(double);
197  // Rotates the Hep3Vector around the z-axis.
198 
199  Hep3Vector & rotateUz(const Hep3Vector&);
200  // Rotates reference frame from Uz to newUz (unit vector) (Geant4).
201 
202  Hep3Vector & rotate(double, const Hep3Vector &);
203  // Rotates around the axis specified by another Hep3Vector.
204  // (Uses methods of HepRotation, forcing linking in of Rotation.cc.)
205 
206  Hep3Vector & operator *= (const HepRotation &);
207  Hep3Vector & transform(const HepRotation &);
208  // Transformation with a Rotation matrix.
209 
210 
211 // = = = = = = = = = = = = = = = = = = = = = = = =
212 //
213 // Esoteric properties and operations on 3-vectors:
214 //
215 // 1 - Set vectors in various coordinate systems
216 // 2 - Synonyms for accessing coordinates and properties
217 // 3 - Comparisions (dictionary, near-ness, and geometric)
218 // 4 - Intrinsic properties
219 // 5 - Properties releative to z axis and arbitrary directions
220 // 6 - Polar and azimuthal angle decomposition and deltaPhi
221 // 7 - Rotations
222 //
223 // = = = = = = = = = = = = = = = = = = = = = = = =
224 
225 // 1 - Set vectors in various coordinate systems
226 
227  inline void setRThetaPhi (double r, double theta, double phi);
228  // Set in spherical coordinates: Angles are measured in RADIANS
229 
230  inline void setREtaPhi ( double r, double eta, double phi );
231  // Set in spherical coordinates, but specify peudorapidiy to determine theta.
232 
233  inline void setRhoPhiZ (double rho, double phi, double z);
234  // Set in cylindrical coordinates: Phi angle is measured in RADIANS
235 
236  void setRhoPhiTheta ( double rho, double phi, double theta);
237  // Set in cylindrical coordinates, but specify theta to determine z.
238 
239  void setRhoPhiEta ( double rho, double phi, double eta);
240  // Set in cylindrical coordinates, but specify pseudorapidity to determine z.
241 
242 // 2 - Synonyms for accessing coordinates and properties
243 
244  inline double getX() const;
245  inline double getY() const;
246  inline double getZ() const;
247  // x(), y(), and z()
248 
249  inline double getR () const;
250  inline double getTheta() const;
251  inline double getPhi () const;
252  // mag(), theta(), and phi()
253 
254  inline double r () const;
255  // mag()
256 
257  inline double rho () const;
258  inline double getRho () const;
259  // perp()
260 
261  double eta () const;
262  double getEta () const;
263  // pseudoRapidity()
264 
265  inline void setR ( double s );
266  // setMag()
267 
268  inline void setRho ( double s );
269  // setPerp()
270 
271 // 3 - Comparisions (dictionary, near-ness, and geometric)
272 
273  int compare (const Hep3Vector & v) const;
274  bool operator > (const Hep3Vector & v) const;
275  bool operator < (const Hep3Vector & v) const;
276  bool operator>= (const Hep3Vector & v) const;
277  bool operator<= (const Hep3Vector & v) const;
278  // dictionary ordering according to z, then y, then x component
279 
280  inline double diff2 (const Hep3Vector & v) const;
281  // |v1-v2|**2
282 
283  static double setTolerance (double tol);
284  static inline double getTolerance ();
285  // Set the tolerance used in isNear() for Hep3Vectors
286 
287  bool isParallel (const Hep3Vector & v, double epsilon=tolerance) const;
288  // Are the vectors parallel, within the given tolerance?
289 
290  bool isOrthogonal (const Hep3Vector & v, double epsilon=tolerance) const;
291  // Are the vectors orthogonal, within the given tolerance?
292 
293  double howParallel (const Hep3Vector & v) const;
294  // | v1.cross(v2) / v1.dot(v2) |, to a maximum of 1.
295 
296  double howOrthogonal (const Hep3Vector & v) const;
297  // | v1.dot(v2) / v1.cross(v2) |, to a maximum of 1.
298 
299  enum { ToleranceTicks = 100 };
300 
301 // 4 - Intrinsic properties
302 
303  double beta () const;
304  // relativistic beta (considering v as a velocity vector with c=1)
305  // Same as mag() but will object if >= 1
306 
307  double gamma() const;
308  // relativistic gamma (considering v as a velocity vector with c=1)
309 
310  double coLinearRapidity() const;
311  // inverse tanh (beta)
312 
313 // 5 - Properties relative to Z axis and to an arbitrary direction
314 
315  // Note that the non-esoteric CLHEP provides
316  // theta(), cosTheta(), cos2Theta, and angle(const Hep3Vector&)
317 
318  inline double angle() const;
319  // angle against the Z axis -- synonym for theta()
320 
321  inline double theta(const Hep3Vector & v2) const;
322  // synonym for angle(v2)
323 
324  double cosTheta (const Hep3Vector & v2) const;
325  double cos2Theta(const Hep3Vector & v2) const;
326  // cos and cos^2 of the angle between two vectors
327 
328  inline Hep3Vector project () const;
329  Hep3Vector project (const Hep3Vector & v2) const;
330  // projection of a vector along a direction.
331 
332  inline Hep3Vector perpPart() const;
333  inline Hep3Vector perpPart (const Hep3Vector & v2) const;
334  // vector minus its projection along a direction.
335 
336  double rapidity () const;
337  // inverse tanh(v.z())
338 
339  double rapidity (const Hep3Vector & v2) const;
340  // rapidity with respect to specified direction:
341  // inverse tanh (v.dot(u)) where u is a unit in the direction of v2
342 
343  double eta(const Hep3Vector & v2) const;
344  // - ln tan of the angle beween the vector and the ref direction.
345 
346 // 6 - Polar and azimuthal angle decomposition and deltaPhi
347 
348  // Decomposition of an angle within reference defined by a direction:
349 
350  double polarAngle (const Hep3Vector & v2) const;
351  // The reference direction is Z: the polarAngle is abs(v.theta()-v2.theta()).
352 
353  double deltaPhi (const Hep3Vector & v2) const;
354  // v.phi()-v2.phi(), brought into the range (-PI,PI]
355 
356  double azimAngle (const Hep3Vector & v2) const;
357  // The reference direction is Z: the azimAngle is the same as deltaPhi
358 
359  double polarAngle (const Hep3Vector & v2,
360  const Hep3Vector & ref) const;
361  // For arbitrary reference direction,
362  // polarAngle is abs(v.angle(ref) - v2.angle(ref)).
363 
364  double azimAngle (const Hep3Vector & v2,
365  const Hep3Vector & ref) const;
366  // To compute azimangle, project v and v2 into the plane normal to
367  // the reference direction. Then in that plane take the angle going
368  // clockwise around the direction from projection of v to that of v2.
369 
370 // 7 - Rotations
371 
372 // These mehtods **DO NOT** use anything in the HepRotation class.
373 // Thus, use of v.rotate(axis,delta) does not force linking in Rotation.cc.
374 
375  Hep3Vector & rotate (const Hep3Vector & axis, double delta);
376  // Synonym for rotate (delta, axis)
377 
378  Hep3Vector & rotate (const HepAxisAngle & ax);
379  // HepAxisAngle is a struct holding an axis direction and an angle.
380 
381  Hep3Vector & rotate (const HepEulerAngles & e);
382  Hep3Vector & rotate (double phi,
383  double theta,
384  double psi);
385  // Rotate via Euler Angles. Our Euler Angles conventions are
386  // those of Goldstein Classical Mechanics page 107.
387 
388 protected:
389  void setSpherical (double r, double theta, double phi);
390  void setCylindrical (double r, double phi, double z);
391  double negativeInfinity() const;
392 
393 protected:
394 
395  double dx;
396  double dy;
397  double dz;
398  // The components.
399 
400  static double tolerance;
401  // default tolerance criterion for isNear() to return true.
402 }; // Hep3Vector
403 
404 // Global Methods
405 
406 Hep3Vector rotationXOf (const Hep3Vector & vec, double delta);
407 Hep3Vector rotationYOf (const Hep3Vector & vec, double delta);
408 Hep3Vector rotationZOf (const Hep3Vector & vec, double delta);
409 
410 Hep3Vector rotationOf (const Hep3Vector & vec,
411  const Hep3Vector & axis, double delta);
412 Hep3Vector rotationOf (const Hep3Vector & vec, const HepAxisAngle & ax);
413 
414 Hep3Vector rotationOf (const Hep3Vector & vec,
415  double phi, double theta, double psi);
416 Hep3Vector rotationOf (const Hep3Vector & vec, const HepEulerAngles & e);
417 // Return a new vector based on a rotation of the supplied vector
418 
419 std::ostream & operator << (std::ostream &, const Hep3Vector &);
420 // Output to a stream.
421 
422 std::istream & operator >> (std::istream &, Hep3Vector &);
423 // Input from a stream.
424 
425 extern const Hep3Vector HepXHat, HepYHat, HepZHat;
426 
429 
430 Hep3Vector operator / (const Hep3Vector &, double a);
431 // Division of 3-vectors by non-zero real number
432 
433 inline Hep3Vector operator + (const Hep3Vector &, const Hep3Vector &);
434 // Addition of 3-vectors.
435 
436 inline Hep3Vector operator - (const Hep3Vector &, const Hep3Vector &);
437 // Subtraction of 3-vectors.
438 
439 inline double operator * (const Hep3Vector &, const Hep3Vector &);
440 // double product of 3-vectors.
441 
442 inline Hep3Vector operator * (const Hep3Vector &, double a);
443 inline Hep3Vector operator * (double a, const Hep3Vector &);
444 // Scaling of 3-vectors with a real number
445 
446 } // namespace CLHEP
447 
448 #include "CLHEP/Vector/ThreeVector.icc"
449 
450 #ifdef ENABLE_BACKWARDS_COMPATIBILITY
451 // backwards compatibility will be enabled ONLY in CLHEP 1.9
452 using namespace CLHEP;
453 #endif
454 
455 #endif /* HEP_THREEVECTOR_H */
delta
HepRotation delta() setPhi()
CLHEP::Hep3Vector::cos2Theta
double cos2Theta() const
CLHEP::Hep3Vector::setEta
void setEta(double p)
Definition: ThreeVector.cc:217
CLHEP::Hep3Vector::howOrthogonal
double howOrthogonal(const Hep3Vector &v) const
Definition: SpaceVector.cc:223
CLHEP::Hep3Vector::deltaR
double deltaR(const Hep3Vector &v) const
Definition: ThreeVector.cc:182
CLHEP::Hep3Vector::setX
void setX(double)
CLHEP::Hep3Vector::getTheta
double getTheta() const
X
Definition: testSharedPtrBasic.cc:28
CLHEP::rotationYOf
HepLorentzVector rotationYOf(const HepLorentzVector &vec, double delta)
Definition: LorentzVectorB.cc:36
CLHEP::Hep3Vector::setY
void setY(double)
CLHEP::Hep3Vector::getPhi
double getPhi() const
CLHEP::Hep3Vector::diff2
double diff2(const Hep3Vector &v) const
a
@ a
Definition: testCategories.cc:125
CLHEP::Hep3Vector::operator==
bool operator==(const Hep3Vector &) const
CLHEP::Hep3Vector::operator/=
Hep3Vector & operator/=(double)
Definition: ThreeVector.cc:341
CLHEP::Hep3Vector::coLinearRapidity
double coLinearRapidity() const
Definition: SpaceVectorP.cc:70
CLHEP::Hep3Vector::compare
int compare(const Hep3Vector &v) const
Definition: SpaceVector.cc:125
CLHEP::Hep3Vector::unit
Hep3Vector unit() const
CLHEP::Hep3Vector::getY
double getY() const
CLHEP::Hep3Vector::~Hep3Vector
~Hep3Vector()
CLHEP::HepYHat
const Hep3Vector HepYHat
Definition: Geometry/CLHEP/Vector/ThreeVector.h:425
CLHEP::Hep3Vector::perpPart
Hep3Vector perpPart() const
CLHEP::Hep3Vector::eta
double eta() const
CLHEP::Hep3Vector::operator<=
bool operator<=(const Hep3Vector &v) const
Definition: SpaceVector.cc:153
CLHEP::Hep3Vector::angle
double angle() const
CLHEP::HepThreeVectorF
Hep3Vector HepThreeVectorF
Definition: Geometry/CLHEP/Vector/ThreeVector.h:428
CLHEP::Hep3Vector::getX
double getX() const
CLHEP::Hep3Vector::howParallel
double howParallel(const Hep3Vector &v) const
Definition: SpaceVector.cc:172
CLHEP::Hep3Vector::tolerance
static double tolerance
Definition: Geometry/CLHEP/Vector/ThreeVector.h:400
CLHEP::Hep3Vector::dot
double dot(const Hep3Vector &) const
CLHEP::Hep3Vector::cosTheta
double cosTheta() const
CLHEP::Hep3Vector::isNear
bool isNear(const Hep3Vector &, double epsilon=tolerance) const
Definition: ThreeVector.cc:154
CLHEP::Hep3Vector::rotateY
Hep3Vector & rotateY(double)
Definition: ThreeVector.cc:134
CLHEP::Hep3Vector::setCylTheta
void setCylTheta(double)
Definition: ThreeVector.cc:243
CLHEP::Hep3Vector::operator()
double operator()(int) const
Definition: ThreeVector.cc:40
CLHEP::Hep3Vector::phi
double phi() const
CLHEP::Hep3Vector::setPerp
void setPerp(double)
CLHEP::Hep3Vector::setRThetaPhi
void setRThetaPhi(double r, double theta, double phi)
CLHEP::Hep3Vector::getRho
double getRho() const
axis
We have the boost methods returning HepLorentzVector &rather than so things can be chained we feel the boost methods along an axis
Definition: minorMergeIssues.doc:155
CLHEP::HepThreeVectorD
Hep3Vector HepThreeVectorD
Definition: Geometry/CLHEP/Vector/ThreeVector.h:427
CLHEP::Hep3Vector::SIZE
@ SIZE
Definition: Geometry/CLHEP/Vector/ThreeVector.h:47
CLHEP::Hep3Vector::operator!=
bool operator!=(const Hep3Vector &) const
CLHEP::Hep3Vector::orthogonal
Hep3Vector orthogonal() const
CLHEP::Hep3Vector::setRho
void setRho(double s)
CLHEP::Hep3Vector::setRhoPhiZ
void setRhoPhiZ(double rho, double phi, double z)
CLHEP::Hep3Vector::rapidity
double rapidity() const
Definition: SpaceVectorP.cc:55
CLHEP::Hep3Vector::rotateX
Hep3Vector & rotateX(double)
Definition: ThreeVector.cc:124
CLHEP::Hep3Vector::dy
double dy
Definition: Geometry/CLHEP/Vector/ThreeVector.h:396
CLHEP::Hep3Vector::setRhoPhiTheta
void setRhoPhiTheta(double rho, double phi, double theta)
Definition: SpaceVector.cc:73
CLHEP::Hep3Vector::x
double x() const
CLHEP::Hep3Vector::isParallel
bool isParallel(const Hep3Vector &v, double epsilon=tolerance) const
Definition: SpaceVector.cc:188
CLHEP::Hep3Vector::azimAngle
double azimAngle(const Hep3Vector &v2) const
CLHEP::Hep3Vector::setREtaPhi
void setREtaPhi(double r, double eta, double phi)
CLHEP::operator/
HepLorentzVector operator/(const HepLorentzVector &, double a)
Definition: LorentzVector.cc:162
CLHEP::Hep3Vector::beta
double beta() const
Definition: SpaceVectorP.cc:32
CLHEP::Hep3Vector::operator=
Hep3Vector & operator=(const Hep3Vector &)
CLHEP::Hep3Vector::getTolerance
static double getTolerance()
CLHEP::Hep3Vector::cross
Hep3Vector cross(const Hep3Vector &) const
CLHEP::Hep3Vector::z
double z() const
Z
Definition: testSharedPtrConvertible.cc:34
CLHEP::Hep3Vector::project
Hep3Vector project() const
psi
HepRotation psi() axis()
CLHEP::Hep3Vector::dz
double dz
Definition: Geometry/CLHEP/Vector/ThreeVector.h:397
CLHEP
Definition: ClhepVersion.h:13
CLHEP::operator-
Hep3Vector operator-(const Hep3Vector &, const Hep3Vector &)
CLHEP::Hep3Vector::rotate
Hep3Vector & rotate(double, const Hep3Vector &)
Definition: ThreeVectorR.cc:29
CLHEP::Hep3Vector::deltaPhi
double deltaPhi(const Hep3Vector &v2) const
Definition: ThreeVector.cc:172
CLHEP::Hep3Vector::mag
double mag() const
v
they are gone ZOOM Features Discontinued The following features of the ZOOM package were felt to be extreme overkill These have been after checking that no existing user code was utilizing as in SpaceVector v
Definition: keyMergeIssues.doc:324
CLHEP::Hep3Vector::ToleranceTicks
@ ToleranceTicks
Definition: Geometry/CLHEP/Vector/ThreeVector.h:299
CLHEP::Hep3Vector::operator*=
Hep3Vector & operator*=(double)
CLHEP::Hep3Vector::setTheta
void setTheta(double)
CLHEP::Hep3Vector::rotateUz
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:72
CLHEP::Hep3Vector::pseudoRapidity
double pseudoRapidity() const
Definition: ThreeVector.cc:92
CLHEP::Hep3Vector::setPhi
void setPhi(double)
Hep3Vector
Issues Concerning the PhysicsVectors CLHEP Vector Merge The merge of ZOOM PhysicsVdectors and the CLHEP Vector package is completed The purpose of this document is to list the major issues that affected the merge of these and where relevant describe the resolutions More detailed documents describe more minor issues General Approach As agreed at the June CLHEP the approach is to combine the features of each ZOOM class with the corresponding CLHEP class expanding the interface to create a single lingua franca of what a Hep3Vector(for example) means. We are not forming SpaceVector as an class derived from Hep3Vector and enhancing it in that way. Another rule imposed by the agreement is to avoid using the Exceptions package(even though that will later go into CLHEP for other uses). A desirable goal is to avoid cluttering the interface and enlarging the code linked in when ordinary CLHEP Vector functionallity is used. To this end
CLHEP::Hep3Vector::operator+=
Hep3Vector & operator+=(const Hep3Vector &)
CLHEP::Hep3Vector::operator>=
bool operator>=(const Hep3Vector &v) const
Definition: SpaceVector.cc:150
CLHEP::Hep3Vector::mag2
double mag2() const
CLHEP::Hep3Vector::howNear
double howNear(const Hep3Vector &v) const
Definition: ThreeVector.cc:159
CLHEP::Hep3Vector::setZ
void setZ(double)
CLHEP::HepXHat
const Hep3Vector HepXHat
Definition: Matrix/CLHEP/Vector/ThreeVector.h:425
CLHEP::Hep3Vector::rho
double rho() const
CLHEP::Hep3Vector::setCylEta
void setCylEta(double p)
Definition: ThreeVector.cc:288
CLHEP::Hep3Vector::operator<
bool operator<(const Hep3Vector &v) const
Definition: SpaceVector.cc:147
CLHEP::rotationOf
HepLorentzVector rotationOf(const HepLorentzVector &vec, const Hep3Vector &axis, double delta)
Definition: LorentzVectorR.cc:48
CLHEP::Hep3Vector::NUM_COORDINATES
@ NUM_COORDINATES
Definition: Geometry/CLHEP/Vector/ThreeVector.h:47
phi
we want to make it possible for the user to use the so we provide a few new for double double phi
Definition: minorMergeIssues.doc:20
CLHEP::Hep3Vector::perp2
double perp2() const
CLHEP::Hep3Vector::operator[]
double operator[](int) const
CLHEP::operator>>
std::istream & operator>>(std::istream &is, HepAxisAngle &aa)
Definition: AxisAngle.cc:96
CLHEP::Hep3Vector::perp
double perp() const
Y
Definition: testSharedPtrBasic.cc:34
CLHEP::Hep3Vector::setR
void setR(double s)
CLHEP::Hep3Vector::dx
double dx
Definition: Geometry/CLHEP/Vector/ThreeVector.h:395
CLHEP::Hep3Vector::theta
double theta() const
CLHEP::HepZHat
const Hep3Vector HepZHat
Definition: Geometry/CLHEP/Vector/ThreeVector.h:425
CLHEP::Hep3Vector::rotateZ
Hep3Vector & rotateZ(double)
Definition: ThreeVector.cc:144
CLHEP::Hep3Vector::operator-=
Hep3Vector & operator-=(const Hep3Vector &)
CLHEP::rotationXOf
HepLorentzVector rotationXOf(const HepLorentzVector &vec, double delta)
Definition: LorentzVectorB.cc:30
CLHEP::Hep3Vector::y
double y() const
CLHEP::operator+
Hep3Vector operator+(const Hep3Vector &, const Hep3Vector &)
CLHEP::Hep3Vector::setMag
void setMag(double)
Definition: ThreeVector.cc:27
CLHEP::Hep3Vector::gamma
double gamma() const
Definition: SpaceVectorP.cc:41
CLHEP::Hep3Vector::getR
double getR() const
CLHEP::Hep3Vector::setCylindrical
void setCylindrical(double r, double phi, double z)
Definition: SpaceVector.cc:58
CLHEP::operator<<
std::ostream & operator<<(std::ostream &os, const HepAxisAngle &aa)
Definition: AxisAngle.cc:86
CLHEP::Hep3Vector::operator>
bool operator>(const Hep3Vector &v) const
Definition: SpaceVector.cc:144
CLHEP::Hep3Vector::r
double r() const
CLHEP::Hep3Vector::setSpherical
void setSpherical(double r, double theta, double phi)
Definition: SpaceVector.cc:37
CLHEP::Hep3Vector::negativeInfinity
double negativeInfinity() const
Definition: SpaceVector.cc:288
CLHEP::Hep3Vector::transform
Hep3Vector & transform(const HepRotation &)
Definition: ThreeVectorR.cc:25
CLHEP::Hep3Vector::getZ
double getZ() const
CLHEP::operator*
HepLorentzRotation operator*(const HepRotation &r, const HepLorentzRotation &lt)
Definition: LorentzRotation.cc:262
CLHEP::Hep3Vector::Hep3Vector
Hep3Vector()
CLHEP::rotationZOf
HepLorentzVector rotationZOf(const HepLorentzVector &vec, double delta)
Definition: LorentzVectorB.cc:42
CLHEP::Hep3Vector::getEta
double getEta() const
CLHEP::Hep3Vector::set
void set(double x, double y, double z)
CLHEP::Hep3Vector::setRhoPhiEta
void setRhoPhiEta(double rho, double phi, double eta)
Definition: SpaceVector.cc:100
CLHEP::Hep3Vector::setTolerance
static double setTolerance(double tol)
Definition: SpaceVector.cc:275
theta
we want to make it possible for the user to use the so we provide a few new for double theta
Definition: minorMergeIssues.doc:20
CLHEP::Hep3Vector::polarAngle
double polarAngle(const Hep3Vector &v2) const
Definition: SpaceVectorD.cc:30
CLHEP::Hep3Vector::operator-
Hep3Vector operator-() const
CLHEP::Hep3Vector::isOrthogonal
bool isOrthogonal(const Hep3Vector &v, double epsilon=tolerance) const
Definition: SpaceVector.cc:241