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

Geometry/CLHEP/Vector/LorentzVector.h
Go to the documentation of this file.
1 // -*- C++ -*-
2 // CLASSDOC OFF
3 // $Id: LorentzVector.h,v 1.2 2003/10/23 21:29:52 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 // HepLorentzVector is a Lorentz vector consisting of Hep3Vector and
10 // double components. Lorentz transformations (rotations and boosts)
11 // of these vectors are perfomed by multiplying with objects of
12 // the HepLorenzRotation class.
13 //
14 // .SS See Also
15 // ThreeVector.h, Rotation.h, LorentzRotation.h
16 //
17 // .SS Authors
18 // Leif Lonnblad and Anders Nilsson. Modified by Evgueni Tcherniaev, Mark Fischler
19 //
20 
21 #ifndef HEP_LORENTZVECTOR_H
22 #define HEP_LORENTZVECTOR_H
23 
24 #ifdef GNUPRAGMA
25 #pragma interface
26 #endif
27 
28 #include <iostream>
29 #include "CLHEP/Vector/defs.h"
30 #include "CLHEP/Vector/ThreeVector.h"
31 
32 namespace CLHEP {
33 
34 // Declarations of classes and global methods
35 class HepLorentzVector;
36 class HepLorentzRotation;
37 class HepRotation;
38 class HepAxisAngle;
39 class HepEulerAngles;
40 class Tcomponent;
41 HepLorentzVector rotationXOf( const HepLorentzVector & vec, double delta );
42 HepLorentzVector rotationYOf( const HepLorentzVector & vec, double delta );
43 HepLorentzVector rotationZOf( const HepLorentzVector & vec, double delta );
44 HepLorentzVector rotationOf
45  ( const HepLorentzVector & vec, const Hep3Vector & axis, double delta );
46 HepLorentzVector rotationOf
47  ( const HepLorentzVector & vec, const HepAxisAngle & ax );
48 HepLorentzVector rotationOf
49  ( const HepLorentzVector & vec, const HepEulerAngles & e1 );
50 HepLorentzVector rotationOf
51  ( const HepLorentzVector & vec, double phi,
52  double theta,
53  double psi );
54 inline
55 HepLorentzVector boostXOf( const HepLorentzVector & vec, double beta );
56 inline
57 HepLorentzVector boostYOf( const HepLorentzVector & vec, double beta );
58 inline
59 HepLorentzVector boostZOf( const HepLorentzVector & vec, double beta );
60 inline HepLorentzVector boostOf
61  ( const HepLorentzVector & vec, const Hep3Vector & betaVector );
62 inline HepLorentzVector boostOf
63  ( const HepLorentzVector & vec, const Hep3Vector & axis, double beta );
64 
66 
67 
73 
74 public:
75 
76  enum { X=0, Y=1, Z=2, T=3, NUM_COORDINATES=4, SIZE=NUM_COORDINATES };
77  // Safe indexing of the coordinates when using with matrices, arrays, etc.
78  // (BaBar)
79 
80  inline HepLorentzVector(double x, double y,
81  double z, double t);
82  // Constructor giving the components x, y, z, t.
83 
84  inline HepLorentzVector(double x, double y, double z);
85  // Constructor giving the components x, y, z with t-component set to 0.0.
86 
87  explicit HepLorentzVector(double t);
88  // Constructor giving the t-component with x, y and z set to 0.0.
89 
90  inline HepLorentzVector();
91  // Default constructor with x, y, z and t set to 0.0.
92 
93  inline HepLorentzVector(const Hep3Vector & p, double e);
94  inline HepLorentzVector(double e, const Hep3Vector & p);
95  // Constructor giving a 3-Vector and a time component.
96 
97  inline HepLorentzVector(const HepLorentzVector &);
98  // Copy constructor.
99 
100  inline ~HepLorentzVector();
101  // The destructor.
102 
103  inline operator const Hep3Vector & () const;
104  inline operator Hep3Vector & ();
105  // Conversion (cast) to Hep3Vector.
106 
107  inline double x() const;
108  inline double y() const;
109  inline double z() const;
110  inline double t() const;
111  // Get position and time.
112 
113  inline void setX(double);
114  inline void setY(double);
115  inline void setZ(double);
116  inline void setT(double);
117  // Set position and time.
118 
119  inline double px() const;
120  inline double py() const;
121  inline double pz() const;
122  inline double e() const;
123  // Get momentum and energy.
124 
125  inline void setPx(double);
126  inline void setPy(double);
127  inline void setPz(double);
128  inline void setE(double);
129  // Set momentum and energy.
130 
131  inline Hep3Vector vect() const;
132  // Get spatial component.
133 
134  inline void setVect(const Hep3Vector &);
135  // Set spatial component.
136 
137  inline double theta() const;
138  inline double cosTheta() const;
139  inline double phi() const;
140  inline double rho() const;
141  // Get spatial vector components in spherical coordinate system.
142 
143  inline void setTheta(double);
144  inline void setPhi(double);
145  inline void setRho(double);
146  // Set spatial vector components in spherical coordinate system.
147 
148  double operator () (int) const;
149  inline double operator [] (int) const;
150  // Get components by index.
151 
152  double & operator () (int);
153  inline double & operator [] (int);
154  // Set components by index.
155 
157  // Assignment.
158 
159  inline HepLorentzVector operator + (const HepLorentzVector &) const;
161  // Additions.
162 
163  inline HepLorentzVector operator - (const HepLorentzVector &) const;
165  // Subtractions.
166 
167  inline HepLorentzVector operator - () const;
168  // Unary minus.
169 
170  inline HepLorentzVector & operator *= (double);
171  HepLorentzVector & operator /= (double);
172  // Scaling with real numbers.
173 
174  inline bool operator == (const HepLorentzVector &) const;
175  inline bool operator != (const HepLorentzVector &) const;
176  // Comparisons.
177 
178  inline double perp2() const;
179  // Transverse component of the spatial vector squared.
180 
181  inline double perp() const;
182  // Transverse component of the spatial vector (R in cylindrical system).
183 
184  inline void setPerp(double);
185  // Set the transverse component of the spatial vector.
186 
187  inline double perp2(const Hep3Vector &) const;
188  // Transverse component of the spatial vector w.r.t. given axis squared.
189 
190  inline double perp(const Hep3Vector &) const;
191  // Transverse component of the spatial vector w.r.t. given axis.
192 
193  inline double angle(const Hep3Vector &) const;
194  // Angle wrt. another vector.
195 
196  inline double mag2() const;
197  // Dot product of 4-vector with itself.
198  // By default the metric is TimePositive, and mag2() is the same as m2().
199 
200  inline double m2() const;
201  // Invariant mass squared.
202 
203  inline double mag() const;
204  inline double m() const;
205  // Invariant mass. If m2() is negative then -sqrt(-m2()) is returned.
206 
207  inline double mt2() const;
208  // Transverse mass squared.
209 
210  inline double mt() const;
211  // Transverse mass.
212 
213  inline double et2() const;
214  // Transverse energy squared.
215 
216  inline double et() const;
217  // Transverse energy.
218 
219  inline double dot(const HepLorentzVector &) const;
220  inline double operator * (const HepLorentzVector &) const;
221  // Scalar product.
222 
223  inline double invariantMass2( const HepLorentzVector & w ) const;
224  // Invariant mass squared of pair of 4-vectors
225 
226  double invariantMass ( const HepLorentzVector & w ) const;
227  // Invariant mass of pair of 4-vectors
228 
229  inline void setVectMag(const Hep3Vector & spatial, double magnitude);
230  inline void setVectM(const Hep3Vector & spatial, double mass);
231  // Copy spatial coordinates, and set energy = sqrt(mass^2 + spatial^2)
232 
233  inline double plus() const;
234  inline double minus() const;
235  // Returns the positive/negative light-cone component t +/- z.
236 
237  Hep3Vector boostVector() const;
238  // Boost needed from rest4Vector in rest frame to form this 4-vector
239  // Returns the spatial components divided by the time component.
240 
241  HepLorentzVector & boost(double, double, double);
242  inline HepLorentzVector & boost(const Hep3Vector &);
243  // Lorentz boost.
244 
245  HepLorentzVector & boostX( double beta );
246  HepLorentzVector & boostY( double beta );
247  HepLorentzVector & boostZ( double beta );
248  // Boost along an axis, by magnitue beta (fraction of speed of light)
249 
250  double rapidity() const;
251  // Returns the rapidity, i.e. 0.5*ln((E+pz)/(E-pz))
252 
253  inline double pseudoRapidity() const;
254  // Returns the pseudo-rapidity, i.e. -ln(tan(theta/2))
255 
256  inline bool isTimelike() const;
257  // Test if the 4-vector is timelike
258 
259  inline bool isSpacelike() const;
260  // Test if the 4-vector is spacelike
261 
262  inline bool isLightlike(double epsilon=tolerance) const;
263  // Test for lightlike is within tolerance epsilon
264 
265  HepLorentzVector & rotateX(double);
266  // Rotate the spatial component around the x-axis.
267 
268  HepLorentzVector & rotateY(double);
269  // Rotate the spatial component around the y-axis.
270 
271  HepLorentzVector & rotateZ(double);
272  // Rotate the spatial component around the z-axis.
273 
275  // Rotates the reference frame from Uz to newUz (unit vector).
276 
277  HepLorentzVector & rotate(double, const Hep3Vector &);
278  // Rotate the spatial component around specified axis.
279 
280  inline HepLorentzVector & operator *= (const HepRotation &);
281  inline HepLorentzVector & transform(const HepRotation &);
282  // Transformation with HepRotation.
283 
286  // Transformation with HepLorenzRotation.
287 
288 // = = = = = = = = = = = = = = = = = = = = = = = =
289 //
290 // Esoteric properties and operations on 4-vectors:
291 //
292 // 0 - Flexible metric convention and axial unit 4-vectors
293 // 1 - Construct and set 4-vectors in various ways
294 // 2 - Synonyms for accessing coordinates and properties
295 // 2a - Setting space coordinates in different ways
296 // 3 - Comparisions (dictionary, near-ness, and geometric)
297 // 4 - Intrinsic properties
298 // 4a - Releativistic kinematic properties
299 // 4b - Methods combining two 4-vectors
300 // 5 - Properties releative to z axis and to arbitrary directions
301 // 7 - Rotations and Boosts
302 //
303 // = = = = = = = = = = = = = = = = = = = = = = = =
304 
305 // 0 - Flexible metric convention
306 
307  static ZMpvMetric_t setMetric( ZMpvMetric_t a1 );
308  static ZMpvMetric_t getMetric();
309 
310 // 1 - Construct and set 4-vectors in various ways
311 
312  inline void set (double x, double y, double z, double t);
313  inline void set (double x, double y, double z, Tcomponent t);
314  inline HepLorentzVector(double x, double y, double z, Tcomponent t);
315  // Form 4-vector by supplying cartesian coordinate components
316 
317  inline void set (Tcomponent t, double x, double y, double z);
318  inline HepLorentzVector(Tcomponent t, double x, double y, double z);
319  // Deprecated because the 4-doubles form uses x,y,z,t, not t,x,y,z.
320 
321  inline void set ( double t );
322 
323  inline void set ( Tcomponent t );
324  inline explicit HepLorentzVector( Tcomponent t );
325  // Form 4-vector with zero space components, by supplying t component
326 
327  inline void set ( const Hep3Vector & v );
328  inline explicit HepLorentzVector( const Hep3Vector & v );
329  // Form 4-vector with zero time component, by supplying space 3-vector
330 
331  inline HepLorentzVector & operator=( const Hep3Vector & v );
332  // Form 4-vector with zero time component, equal to space 3-vector
333 
334  inline void set ( const Hep3Vector & v, double t );
335  inline void set ( double t, const Hep3Vector & v );
336  // Set using specified space vector and time component
337 
338 // 2 - Synonyms for accessing coordinates and properties
339 
340  inline double getX() const;
341  inline double getY() const;
342  inline double getZ() const;
343  inline double getT() const;
344  // Get position and time.
345 
346  inline Hep3Vector v() const;
347  inline Hep3Vector getV() const;
348  // Get spatial component. Same as vect.
349 
350  inline void setV(const Hep3Vector &);
351  // Set spatial component. Same as setVect.
352 
353 // 2a - Setting space coordinates in different ways
354 
355  inline void setV( double x, double y, double z );
356 
357  inline void setRThetaPhi( double r, double theta, double phi);
358  inline void setREtaPhi( double r, double eta, double phi);
359  inline void setRhoPhiZ( double rho, double phi, double z );
360 
361 // 3 - Comparisions (dictionary, near-ness, and geometric)
362 
363  int compare( const HepLorentzVector & w ) const;
364 
365  bool operator >( const HepLorentzVector & w ) const;
366  bool operator <( const HepLorentzVector & w ) const;
367  bool operator>=( const HepLorentzVector & w ) const;
368  bool operator<=( const HepLorentzVector & w ) const;
369 
370  bool isNear ( const HepLorentzVector & w,
371  double epsilon=tolerance ) const;
372  double howNear( const HepLorentzVector & w ) const;
373  // Is near using Euclidean measure t**2 + v**2
374 
375  bool isNearCM ( const HepLorentzVector & w,
376  double epsilon=tolerance ) const;
377  double howNearCM( const HepLorentzVector & w ) const;
378  // Is near in CM frame: Applicable only for two timelike HepLorentzVectors
379 
380  // If w1 and w2 are already in their CM frame, then w1.isNearCM(w2)
381  // is exactly equivalent to w1.isNear(w2).
382  // If w1 and w2 have T components of zero, w1.isNear(w2) is exactly
383  // equivalent to w1.getV().isNear(w2.v()).
384 
385  bool isParallel( const HepLorentzVector & w,
386  double epsilon=tolerance ) const;
387  // Test for isParallel is within tolerance epsilon
388  double howParallel (const HepLorentzVector & w) const;
389 
390  static double getTolerance();
391  static double setTolerance( double tol );
392  // Set the tolerance for HepLorentzVectors to be considered near
393  // The same tolerance is used for determining isLightlike, and isParallel
394 
395  double deltaR(const HepLorentzVector & v) const;
396  // sqrt ( (delta eta)^2 + (delta phi)^2 ) of space part
397 
398 // 4 - Intrinsic properties
399 
400  double howLightlike() const;
401  // Close to zero for almost lightlike 4-vectors; up to 1.
402 
403  inline double euclideanNorm2() const;
404  // Sum of the squares of time and space components; not Lorentz invariant.
405 
406  inline double euclideanNorm() const;
407  // Length considering the metric as (+ + + +); not Lorentz invariant.
408 
409 
410 // 4a - Relativistic kinematic properties
411 
412 // All Relativistic kinematic properties are independent of the sense of metric
413 
414  inline double restMass2() const;
415  inline double invariantMass2() const;
416  // Rest mass squared -- same as m2()
417 
418  inline double restMass() const;
419  inline double invariantMass() const;
420  // Same as m(). If m2() is negative then -sqrt(-m2()) is returned.
421 
422 // The following properties are rest-frame related,
423 // and are applicable only to non-spacelike 4-vectors
424 
426  // This 4-vector, boosted into its own rest frame: (0, 0, 0, m())
427  // The following relation holds by definition:
428  // w.rest4Vector().boost(w.boostVector()) == w
429 
430  // Beta and gamma of the boost vector
431  double beta() const;
432  // Relativistic beta of the boost vector
433 
434  double gamma() const;
435  // Relativistic gamma of the boost vector
436 
437  inline double eta() const;
438  // Pseudorapidity (of the space part)
439 
440  inline double eta(const Hep3Vector & ref) const;
441  // Pseudorapidity (of the space part) w.r.t. specified direction
442 
443  double rapidity(const Hep3Vector & ref) const;
444  // Rapidity in specified direction
445 
446  double coLinearRapidity() const;
447  // Rapidity, in the relativity textbook sense: atanh (|P|/E)
448 
449  Hep3Vector findBoostToCM() const;
450  // Boost needed to get to center-of-mass frame:
451  // w.findBoostToCM() == - w.boostVector()
452  // w.boost(w.findBoostToCM()) == w.rest4Vector()
453 
454  Hep3Vector findBoostToCM( const HepLorentzVector & w ) const;
455  // Boost needed to get to combined center-of-mass frame:
456  // w1.findBoostToCM(w2) == w2.findBoostToCM(w1)
457  // w.findBoostToCM(w) == w.findBoostToCM()
458 
459  inline double et2(const Hep3Vector &) const;
460  // Transverse energy w.r.t. given axis squared.
461 
462  inline double et(const Hep3Vector &) const;
463  // Transverse energy w.r.t. given axis.
464 
465 // 4b - Methods combining two 4-vectors
466 
467  inline double diff2( const HepLorentzVector & w ) const;
468  // (this - w).dot(this-w); sign depends on metric choice
469 
470  inline double delta2Euclidean ( const HepLorentzVector & w ) const;
471  // Euclidean norm of differnce: (delta_T)^2 + (delta_V)^2
472 
473 // 5 - Properties releative to z axis and to arbitrary directions
474 
475  double plus( const Hep3Vector & ref ) const;
476  // t + projection in reference direction
477 
478  double minus( const Hep3Vector & ref ) const;
479  // t - projection in reference direction
480 
481 // 7 - Rotations and boosts
482 
483  HepLorentzVector & rotate ( const Hep3Vector & axis, double delta );
484  // Same as rotate (delta, axis)
485 
486  HepLorentzVector & rotate ( const HepAxisAngle & ax );
487  HepLorentzVector & rotate ( const HepEulerAngles & e );
488  HepLorentzVector & rotate ( double phi,
489  double theta,
490  double psi );
491  // Rotate using these HepEuler angles - see Goldstein page 107 for conventions
492 
493  HepLorentzVector & boost ( const Hep3Vector & axis, double beta );
494  // Normalizes the Hep3Vector to define a direction, and uses beta to
495  // define the magnitude of the boost.
496 
498  ( const HepLorentzVector & vec, double delta );
500  ( const HepLorentzVector & vec, double delta );
502  ( const HepLorentzVector & vec, double delta );
504  ( const HepLorentzVector & vec, const Hep3Vector & axis, double delta );
506  ( const HepLorentzVector & vec, const HepAxisAngle & ax );
508  ( const HepLorentzVector & vec, const HepEulerAngles & e );
510  ( const HepLorentzVector & vec, double phi,
511  double theta,
512  double psi );
513 
514  inline friend HepLorentzVector boostXOf
515  ( const HepLorentzVector & vec, double beta );
516  inline friend HepLorentzVector boostYOf
517  ( const HepLorentzVector & vec, double beta );
518  inline friend HepLorentzVector boostZOf
519  ( const HepLorentzVector & vec, double beta );
520  inline friend HepLorentzVector boostOf
521  ( const HepLorentzVector & vec, const Hep3Vector & betaVector );
522  inline friend HepLorentzVector boostOf
523  ( const HepLorentzVector & vec, const Hep3Vector & axis, double beta );
524 
525 private:
526 
527  Hep3Vector pp;
528  double ee;
529 
530  static double tolerance;
531  static double metric;
532 
533 }; // HepLorentzVector
534 
535 // 8 - Axial Unit 4-vectors
536 
537 static const HepLorentzVector X_HAT4 = HepLorentzVector( 1, 0, 0, 0 );
538 static const HepLorentzVector Y_HAT4 = HepLorentzVector( 0, 1, 0, 0 );
539 static const HepLorentzVector Z_HAT4 = HepLorentzVector( 0, 0, 1, 0 );
540 static const HepLorentzVector T_HAT4 = HepLorentzVector( 0, 0, 0, 1 );
541 
542 // Global methods
543 
544 std::ostream & operator << (std::ostream &, const HepLorentzVector &);
545 // Output to a stream.
546 
547 std::istream & operator >> (std::istream &, HepLorentzVector &);
548 // Input from a stream.
549 
552 
553 inline HepLorentzVector operator * (const HepLorentzVector &, double a);
554 inline HepLorentzVector operator * (double a, const HepLorentzVector &);
555 // Scaling LorentzVector with a real number
556 
558 // Dividing LorentzVector by a real number
559 
560 // Tcomponent definition:
561 
562 // Signature protection for 4-vector constructors taking 4 components
563 class Tcomponent {
564 private:
565  double t_;
566 public:
567  explicit Tcomponent(double t) : t_(t) {}
568  operator double() const { return t_; }
569 }; // Tcomponent
570 
571 } // namespace CLHEP
572 
573 #include "CLHEP/Vector/LorentzVector.icc"
574 
575 #ifdef ENABLE_BACKWARDS_COMPATIBILITY
576 // backwards compatibility will be enabled ONLY in CLHEP 1.9
577 using namespace CLHEP;
578 #endif
579 
580 #endif /* HEP_LORENTZVECTOR_H */
delta
HepRotation delta() setPhi()
CLHEP::HepLorentzVector::rapidity
double rapidity() const
Definition: LorentzVectorK.cc:124
CLHEP::HepLorentzVector::setRThetaPhi
void setRThetaPhi(double r, double theta, double phi)
X
Definition: testSharedPtrBasic.cc:28
double
#define double(obj)
Definition: excDblThrow.cc:32
CLHEP::rotationYOf
HepLorentzVector rotationYOf(const HepLorentzVector &vec, double delta)
Definition: LorentzVectorB.cc:36
CLHEP::HepAxisAngle
Definition: Geometry/CLHEP/Vector/AxisAngle.h:37
CLHEP::ZMpvMetric_t
ZMpvMetric_t
Definition: Geometry/CLHEP/Vector/LorentzVector.h:65
CLHEP::HepLorentzVector::delta2Euclidean
double delta2Euclidean(const HepLorentzVector &w) const
CLHEP::HepLorentzVector::deltaR
double deltaR(const HepLorentzVector &v) const
Definition: LorentzVectorC.cc:196
CLHEP::HepLorentzVector::phi
double phi() const
CLHEP::HepLorentzVector::rotateY
HepLorentzVector & rotateY(double)
Definition: LorentzVector.cc:74
CLHEP::HepLorentzVector::mag
double mag() const
CLHEP::HepLorentzVector::operator*=
HepLorentzVector & operator*=(double)
a
@ a
Definition: testCategories.cc:125
CLHEP::HepLorentzVector::operator()
double operator()(int) const
Definition: LorentzVector.cc:24
CLHEP::HepLorentzVector::setRhoPhiZ
void setRhoPhiZ(double rho, double phi, double z)
CLHEP::HepLorentzVector::operator<=
bool operator<=(const HepLorentzVector &w) const
Definition: LorentzVectorC.cc:49
CLHEP::HepLorentzVector::transform
HepLorentzVector & transform(const HepRotation &)
CLHEP::HepLorentzVector::y
double y() const
CLHEP::HepLorentzVector::rotateZ
HepLorentzVector & rotateZ(double)
Definition: LorentzVector.cc:78
CLHEP::HepLorentzVector::set
void set(double x, double y, double z, double t)
CLHEP::HepLorentzVector::z
double z() const
CLHEP::HepLorentzVector::isTimelike
bool isTimelike() const
CLHEP::HepLorentzVector::px
double px() const
CLHEP::boostYOf
HepLorentzVector boostYOf(const HepLorentzVector &vec, double beta)
CLHEP::HepLorentzVector::HepLorentzVector
HepLorentzVector()
CLHEP::HepLorentzVector::x
double x() const
CLHEP::HepLorentzVector::operator-=
HepLorentzVector & operator-=(const HepLorentzVector &)
CLHEP::HepLorentzVector::operator*
double operator*(const HepLorentzVector &) const
CLHEP::HepLorentzVector::setE
void setE(double)
CLHEP::HepLorentzVector::mt
double mt() const
CLHEP::HepLorentzVector::isNearCM
bool isNearCM(const HepLorentzVector &w, double epsilon=tolerance) const
Definition: LorentzVectorC.cc:86
CLHEP::HepLorentzVectorF
HepLorentzVector HepLorentzVectorF
Definition: Geometry/CLHEP/Vector/LorentzVector.h:551
CLHEP::HepLorentzVector::dot
double dot(const HepLorentzVector &) const
CLHEP::HepLorentzVector::rotationOf
friend HepLorentzVector rotationOf(const HepLorentzVector &vec, const Hep3Vector &axis, double delta)
Definition: LorentzVectorR.cc:48
CLHEP::HepLorentzVector::~HepLorentzVector
~HepLorentzVector()
CLHEP::HepLorentzVector::vect
Hep3Vector vect() const
CLHEP::HepLorentzVector::eta
double eta() 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::HepLorentzVector::boostYOf
friend HepLorentzVector boostYOf(const HepLorentzVector &vec, double beta)
CLHEP::HepLorentzVector::angle
double angle(const Hep3Vector &) const
CLHEP::HepLorentzVector::restMass
double restMass() const
CLHEP::HepLorentzVector::howParallel
double howParallel(const HepLorentzVector &w) const
Definition: LorentzVectorC.cc:228
CLHEP::HepLorentzVector::boostOf
friend HepLorentzVector boostOf(const HepLorentzVector &vec, const Hep3Vector &betaVector)
CLHEP::HepLorentzVector::getT
double getT() const
CLHEP::HepLorentzVector::operator>
bool operator>(const HepLorentzVector &w) const
Definition: LorentzVectorC.cc:40
CLHEP::HepLorentzVector::setT
void setT(double)
CLHEP::TimeNegative
@ TimeNegative
Definition: Geometry/CLHEP/Vector/LorentzVector.h:65
CLHEP::HepLorentzVector::operator>=
bool operator>=(const HepLorentzVector &w) const
Definition: LorentzVectorC.cc:46
CLHEP::TimePositive
@ TimePositive
Definition: Geometry/CLHEP/Vector/LorentzVector.h:65
CLHEP::HepLorentzVector::rotationZOf
friend HepLorentzVector rotationZOf(const HepLorentzVector &vec, double delta)
Definition: LorentzVectorB.cc:42
CLHEP::HepLorentzVector::boostXOf
friend HepLorentzVector boostXOf(const HepLorentzVector &vec, double beta)
CLHEP::HepLorentzVector::m2
double m2() const
CLHEP::HepLorentzVector::setPy
void setPy(double)
CLHEP::boostZOf
HepLorentzVector boostZOf(const HepLorentzVector &vec, double beta)
CLHEP::HepLorentzVector::v
Hep3Vector v() const
CLHEP::HepRotation
Definition: Geometry/CLHEP/Vector/Rotation.h:48
CLHEP::HepLorentzVector::isParallel
bool isParallel(const HepLorentzVector &w, double epsilon=tolerance) const
Definition: LorentzVectorC.cc:209
CLHEP::HepLorentzVector::setPx
void setPx(double)
CLHEP::HepLorentzVector::operator=
HepLorentzVector & operator=(const HepLorentzVector &)
CLHEP::HepLorentzVector::operator+
HepLorentzVector operator+(const HepLorentzVector &) const
CLHEP::HepLorentzRotation
Definition: Geometry/CLHEP/Vector/LorentzRotation.h:54
CLHEP::HepLorentzVector::t
double t() const
CLHEP::HepLorentzVector::setX
void setX(double)
CLHEP::HepLorentzVector::getZ
double getZ() const
CLHEP::HepLorentzVector::setRho
void setRho(double)
CLHEP::HepLorentzVector::coLinearRapidity
double coLinearRapidity() const
Definition: LorentzVectorK.cc:162
CLHEP::HepLorentzVector::setY
void setY(double)
CLHEP::HepLorentzVector::operator==
bool operator==(const HepLorentzVector &) const
CLHEP::operator/
HepLorentzVector operator/(const HepLorentzVector &, double a)
Definition: LorentzVector.cc:162
CLHEP::HepLorentzVector::setVectM
void setVectM(const Hep3Vector &spatial, double mass)
CLHEP::HepLorentzVector::invariantMass2
double invariantMass2() const
CLHEP::HepLorentzVector::perp
double perp() const
CLHEP::HepLorentzVector::m
double m() const
CLHEP::HepLorentzVector::operator+=
HepLorentzVector & operator+=(const HepLorentzVector &)
Z
Definition: testSharedPtrConvertible.cc:34
psi
HepRotation psi() axis()
CLHEP::HepLorentzVector::et2
double et2() const
CLHEP::HepLorentzVector::isNear
bool isNear(const HepLorentzVector &w, double epsilon=tolerance) const
Definition: LorentzVectorC.cc:58
CLHEP::HepLorentzVector::isSpacelike
bool isSpacelike() const
CLHEP::HepLorentzVector::rest4Vector
HepLorentzVector rest4Vector() const
Definition: LorentzVectorK.cc:66
CLHEP::HepLorentzVector::e
double e() const
CLHEP::HepLorentzVector::setPz
void setPz(double)
CLHEP
Definition: ClhepVersion.h:13
CLHEP::HepLorentzVector::plus
double plus() const
CLHEP::HepLorentzVector::py
double py() const
CLHEP::HepLorentzVector::getMetric
static ZMpvMetric_t getMetric()
Definition: LorentzVectorK.cc:37
CLHEP::HepLorentzVector::boostX
HepLorentzVector & boostX(double beta)
Definition: LorentzVector.cc:192
CLHEP::HepLorentzVector::setPhi
void setPhi(double)
CLHEP::HepLorentzVector::getV
Hep3Vector getV() const
CLHEP::HepLorentzVectorD
HepLorentzVector HepLorentzVectorD
Definition: Geometry/CLHEP/Vector/LorentzVector.h:550
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::HepLorentzVector::mag2
double mag2() const
CLHEP::HepLorentzVector::et
double et() const
CLHEP::HepLorentzVector::euclideanNorm
double euclideanNorm() const
CLHEP::HepLorentzVector::T
@ T
Definition: Geometry/CLHEP/Vector/LorentzVector.h:76
CLHEP::Hep3Vector
Definition: Geometry/CLHEP/Vector/ThreeVector.h:41
CLHEP::rotationOf
HepLorentzVector rotationOf(const HepLorentzVector &vec, const Hep3Vector &axis, double delta)
Definition: LorentzVectorR.cc:48
CLHEP::HepLorentzVector::pz
double pz() const
CLHEP::HepLorentzVector::setVectMag
void setVectMag(const Hep3Vector &spatial, double magnitude)
CLHEP::HepLorentzVector::boostY
HepLorentzVector & boostY(double beta)
Definition: LorentzVector.cc:206
CLHEP::HepLorentzVector::restMass2
double restMass2() const
CLHEP::HepLorentzVector::setPerp
void setPerp(double)
CLHEP::HepLorentzVector::pseudoRapidity
double pseudoRapidity() const
CLHEP::HepLorentzVector::setTolerance
static double setTolerance(double tol)
Definition: LorentzVector.cc:234
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::HepLorentzVector::rotateX
HepLorentzVector & rotateX(double)
Definition: LorentzVector.cc:70
CLHEP::HepLorentzVector::howLightlike
double howLightlike() const
Definition: LorentzVectorC.cc:250
CLHEP::HepLorentzVector::gamma
double gamma() const
Definition: LorentzVectorK.cc:93
CLHEP::operator>>
std::istream & operator>>(std::istream &is, HepAxisAngle &aa)
Definition: AxisAngle.cc:96
CLHEP::HepLorentzVector::SIZE
@ SIZE
Definition: Geometry/CLHEP/Vector/LorentzVector.h:76
CLHEP::HepLorentzVector::boostVector
Hep3Vector boostVector() const
Definition: LorentzVector.cc:173
Y
Definition: testSharedPtrBasic.cc:34
CLHEP::boostXOf
HepLorentzVector boostXOf(const HepLorentzVector &vec, double beta)
CLHEP::HepLorentzVector::operator[]
double operator[](int) const
CLHEP::HepLorentzVector::boost
HepLorentzVector & boost(double, double, double)
Definition: LorentzVector.cc:57
CLHEP::HepLorentzVector::rotate
HepLorentzVector & rotate(double, const Hep3Vector &)
Definition: LorentzVectorR.cc:20
CLHEP::HepLorentzVector::minus
double minus() const
CLHEP::rotationXOf
HepLorentzVector rotationXOf(const HepLorentzVector &vec, double delta)
Definition: LorentzVectorB.cc:30
CLHEP::HepLorentzVector::cosTheta
double cosTheta() const
CLHEP::HepLorentzVector::setTheta
void setTheta(double)
CLHEP::boostOf
HepLorentzVector boostOf(const HepLorentzVector &vec, const Hep3Vector &betaVector)
CLHEP::HepLorentzVector::howNearCM
double howNearCM(const HepLorentzVector &w) const
Definition: LorentzVectorC.cc:134
CLHEP::HepLorentzVector::boostZ
HepLorentzVector & boostZ(double beta)
Definition: LorentzVector.cc:220
CLHEP::HepLorentzVector::beta
double beta() const
Definition: LorentzVectorK.cc:75
CLHEP::HepLorentzVector::operator-
HepLorentzVector operator-() const
CLHEP::HepLorentzVector::operator/=
HepLorentzVector & operator/=(double)
Definition: LorentzVector.cc:150
CLHEP::HepLorentzVector::operator!=
bool operator!=(const HepLorentzVector &) const
CLHEP::HepLorentzVector::rotationYOf
friend HepLorentzVector rotationYOf(const HepLorentzVector &vec, double delta)
Definition: LorentzVectorB.cc:36
CLHEP::HepLorentzVector::setZ
void setZ(double)
CLHEP::HepLorentzVector
Definition: Geometry/CLHEP/Vector/LorentzVector.h:72
CLHEP::HepLorentzVector::getY
double getY() const
CLHEP::operator<<
std::ostream & operator<<(std::ostream &os, const HepAxisAngle &aa)
Definition: AxisAngle.cc:86
CLHEP::HepLorentzVector::setMetric
static ZMpvMetric_t setMetric(ZMpvMetric_t a1)
Definition: LorentzVectorK.cc:27
CLHEP::HepEulerAngles
Definition: Geometry/CLHEP/Vector/EulerAngles.h:32
CLHEP::HepLorentzVector::euclideanNorm2
double euclideanNorm2() const
CLHEP::HepLorentzVector::invariantMass
double invariantMass() const
CLHEP::HepLorentzVector::compare
int compare(const HepLorentzVector &w) const
Definition: LorentzVectorC.cc:30
CLHEP::HepLorentzVector::howNear
double howNear(const HepLorentzVector &w) const
Definition: LorentzVectorC.cc:68
CLHEP::HepLorentzVector::perp2
double perp2() const
CLHEP::HepLorentzVector::operator<
bool operator<(const HepLorentzVector &w) const
Definition: LorentzVectorC.cc:43
CLHEP::Tcomponent
Definition: Geometry/CLHEP/Vector/LorentzVector.h:563
CLHEP::operator*
HepLorentzRotation operator*(const HepRotation &r, const HepLorentzRotation &lt)
Definition: LorentzRotation.cc:262
CLHEP::HepLorentzVector::boostZOf
friend HepLorentzVector boostZOf(const HepLorentzVector &vec, double beta)
CLHEP::HepLorentzVector::setVect
void setVect(const Hep3Vector &)
CLHEP::HepLorentzVector::theta
double theta() const
CLHEP::rotationZOf
HepLorentzVector rotationZOf(const HepLorentzVector &vec, double delta)
Definition: LorentzVectorB.cc:42
CLHEP::HepLorentzVector::NUM_COORDINATES
@ NUM_COORDINATES
Definition: Geometry/CLHEP/Vector/LorentzVector.h:76
CLHEP::HepLorentzVector::setREtaPhi
void setREtaPhi(double r, double eta, double phi)
CLHEP::Tcomponent::Tcomponent
Tcomponent(double t)
Definition: Geometry/CLHEP/Vector/LorentzVector.h:567
CLHEP::HepLorentzVector::rotateUz
HepLorentzVector & rotateUz(const Hep3Vector &)
Definition: LorentzVector.cc:83
CLHEP::HepLorentzVector::rotationXOf
friend HepLorentzVector rotationXOf(const HepLorentzVector &vec, double delta)
Definition: LorentzVectorB.cc:30
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::HepLorentzVector::isLightlike
bool isLightlike(double epsilon=tolerance) const
CLHEP::HepLorentzVector::getX
double getX() const
CLHEP::HepLorentzVector::rho
double rho() const
CLHEP::HepLorentzVector::mt2
double mt2() const
CLHEP::HepLorentzVector::findBoostToCM
Hep3Vector findBoostToCM() const
Definition: LorentzVectorK.cc:211
CLHEP::HepLorentzVector::getTolerance
static double getTolerance()
Definition: LorentzVector.cc:241
CLHEP::HepLorentzVector::setV
void setV(const Hep3Vector &)
CLHEP::HepLorentzVector::diff2
double diff2(const HepLorentzVector &w) const