CLHEP 2.4.6.4
C++ Class Library for High Energy Physics
Loading...
Searching...
No Matches
CLHEP::HepLorentzVector Class Reference

#include <LorentzVector.h>

Public Types

enum  {
  X =0 , Y =1 , Z =2 , T =3 ,
  NUM_COORDINATES =4 , SIZE =NUM_COORDINATES
}
 

Public Member Functions

 HepLorentzVector (double x, double y, double z, double t)
 
 HepLorentzVector (double x, double y, double z)
 
 HepLorentzVector (double t)
 
 HepLorentzVector ()
 
 HepLorentzVector (const Hep3Vector &p, double e)
 
 HepLorentzVector (double e, const Hep3Vector &p)
 
 HepLorentzVector (const HepLorentzVector &)
 
 HepLorentzVector (HepLorentzVector &&)=default
 
 ~HepLorentzVector ()
 
 operator const Hep3Vector & () const
 
 operator Hep3Vector & ()
 
double x () const
 
double y () const
 
double z () const
 
double t () const
 
void setX (double)
 
void setY (double)
 
void setZ (double)
 
void setT (double)
 
double px () const
 
double py () const
 
double pz () const
 
double e () const
 
void setPx (double)
 
void setPy (double)
 
void setPz (double)
 
void setE (double)
 
Hep3Vector vect () const
 
void setVect (const Hep3Vector &)
 
double theta () const
 
double cosTheta () const
 
double phi () const
 
double rho () const
 
void setTheta (double)
 
void setPhi (double)
 
void setRho (double)
 
double operator() (int) const
 
double operator[] (int) const
 
doubleoperator() (int)
 
doubleoperator[] (int)
 
HepLorentzVectoroperator= (const HepLorentzVector &)
 
HepLorentzVectoroperator= (HepLorentzVector &&)=default
 
HepLorentzVector operator+ (const HepLorentzVector &) const
 
HepLorentzVectoroperator+= (const HepLorentzVector &)
 
HepLorentzVector operator- (const HepLorentzVector &) const
 
HepLorentzVectoroperator-= (const HepLorentzVector &)
 
HepLorentzVector operator- () const
 
HepLorentzVectoroperator*= (double)
 
HepLorentzVectoroperator/= (double)
 
bool operator== (const HepLorentzVector &) const
 
bool operator!= (const HepLorentzVector &) const
 
double perp2 () const
 
double perp () const
 
void setPerp (double)
 
double perp2 (const Hep3Vector &) const
 
double perp (const Hep3Vector &) const
 
double angle (const Hep3Vector &) const
 
double mag2 () const
 
double m2 () const
 
double mag () const
 
double m () const
 
double mt2 () const
 
double mt () const
 
double et2 () const
 
double et () const
 
double dot (const HepLorentzVector &) const
 
double operator* (const HepLorentzVector &) const
 
double invariantMass2 (const HepLorentzVector &w) const
 
double invariantMass (const HepLorentzVector &w) const
 
void setVectMag (const Hep3Vector &spatial, double magnitude)
 
void setVectM (const Hep3Vector &spatial, double mass)
 
double plus () const
 
double minus () const
 
Hep3Vector boostVector () const
 
HepLorentzVectorboost (double, double, double)
 
HepLorentzVectorboost (const Hep3Vector &)
 
HepLorentzVectorboostX (double beta)
 
HepLorentzVectorboostY (double beta)
 
HepLorentzVectorboostZ (double beta)
 
double rapidity () const
 
double pseudoRapidity () const
 
bool isTimelike () const
 
bool isSpacelike () const
 
bool isLightlike (double epsilon=tolerance) const
 
HepLorentzVectorrotateX (double)
 
HepLorentzVectorrotateY (double)
 
HepLorentzVectorrotateZ (double)
 
HepLorentzVectorrotateUz (const Hep3Vector &)
 
HepLorentzVectorrotate (double, const Hep3Vector &)
 
HepLorentzVectoroperator*= (const HepRotation &)
 
HepLorentzVectortransform (const HepRotation &)
 
HepLorentzVectoroperator*= (const HepLorentzRotation &)
 
HepLorentzVectortransform (const HepLorentzRotation &)
 
void set (double x, double y, double z, double t)
 
void set (double x, double y, double z, Tcomponent t)
 
 HepLorentzVector (double x, double y, double z, Tcomponent t)
 
void set (Tcomponent t, double x, double y, double z)
 
 HepLorentzVector (Tcomponent t, double x, double y, double z)
 
void set (double t)
 
void set (Tcomponent t)
 
 HepLorentzVector (Tcomponent t)
 
void set (const Hep3Vector &v)
 
 HepLorentzVector (const Hep3Vector &v)
 
HepLorentzVectoroperator= (const Hep3Vector &v)
 
void set (const Hep3Vector &v, double t)
 
void set (double t, const Hep3Vector &v)
 
double getX () const
 
double getY () const
 
double getZ () const
 
double getT () const
 
Hep3Vector v () const
 
Hep3Vector getV () const
 
void setV (const Hep3Vector &)
 
void setV (double x, double y, double z)
 
void setRThetaPhi (double r, double theta, double phi)
 
void setREtaPhi (double r, double eta, double phi)
 
void setRhoPhiZ (double rho, double phi, double z)
 
int compare (const HepLorentzVector &w) const
 
bool operator> (const HepLorentzVector &w) const
 
bool operator< (const HepLorentzVector &w) const
 
bool operator>= (const HepLorentzVector &w) const
 
bool operator<= (const HepLorentzVector &w) const
 
bool isNear (const HepLorentzVector &w, double epsilon=tolerance) const
 
double howNear (const HepLorentzVector &w) const
 
bool isNearCM (const HepLorentzVector &w, double epsilon=tolerance) const
 
double howNearCM (const HepLorentzVector &w) const
 
bool isParallel (const HepLorentzVector &w, double epsilon=tolerance) const
 
double howParallel (const HepLorentzVector &w) const
 
double deltaR (const HepLorentzVector &v) const
 
double howLightlike () const
 
double euclideanNorm2 () const
 
double euclideanNorm () const
 
double restMass2 () const
 
double invariantMass2 () const
 
double restMass () const
 
double invariantMass () const
 
HepLorentzVector rest4Vector () const
 
double beta () const
 
double gamma () const
 
double eta () const
 
double eta (const Hep3Vector &ref) const
 
double rapidity (const Hep3Vector &ref) const
 
double coLinearRapidity () const
 
Hep3Vector findBoostToCM () const
 
Hep3Vector findBoostToCM (const HepLorentzVector &w) const
 
double et2 (const Hep3Vector &) const
 
double et (const Hep3Vector &) const
 
double diff2 (const HepLorentzVector &w) const
 
double delta2Euclidean (const HepLorentzVector &w) const
 
double plus (const Hep3Vector &ref) const
 
double minus (const Hep3Vector &ref) const
 
HepLorentzVectorrotate (const Hep3Vector &axis, double delta)
 
HepLorentzVectorrotate (const HepAxisAngle &ax)
 
HepLorentzVectorrotate (const HepEulerAngles &e)
 
HepLorentzVectorrotate (double phi, double theta, double psi)
 
HepLorentzVectorboost (const Hep3Vector &axis, double beta)
 

Static Public Member Functions

static ZMpvMetric_t setMetric (ZMpvMetric_t a1)
 
static ZMpvMetric_t getMetric ()
 
static double getTolerance ()
 
static double setTolerance (double tol)
 

Friends

HepLorentzVector rotationXOf (const HepLorentzVector &vec, double delta)
 
HepLorentzVector rotationYOf (const HepLorentzVector &vec, double delta)
 
HepLorentzVector rotationZOf (const HepLorentzVector &vec, double delta)
 
HepLorentzVector rotationOf (const HepLorentzVector &vec, const Hep3Vector &axis, double delta)
 
HepLorentzVector rotationOf (const HepLorentzVector &vec, const HepAxisAngle &ax)
 
HepLorentzVector rotationOf (const HepLorentzVector &vec, const HepEulerAngles &e)
 
HepLorentzVector rotationOf (const HepLorentzVector &vec, double phi, double theta, double psi)
 
HepLorentzVector boostXOf (const HepLorentzVector &vec, double beta)
 
HepLorentzVector boostYOf (const HepLorentzVector &vec, double beta)
 
HepLorentzVector boostZOf (const HepLorentzVector &vec, double beta)
 
HepLorentzVector boostOf (const HepLorentzVector &vec, const Hep3Vector &betaVector)
 
HepLorentzVector boostOf (const HepLorentzVector &vec, const Hep3Vector &axis, double beta)
 

Detailed Description

Author

Definition at line 68 of file LorentzVector.h.

Member Enumeration Documentation

◆ anonymous enum

anonymous enum
Enumerator
NUM_COORDINATES 
SIZE 

Definition at line 72 of file LorentzVector.h.

Constructor & Destructor Documentation

◆ HepLorentzVector() [1/12]

CLHEP::HepLorentzVector::HepLorentzVector ( double  x,
double  y,
double  z,
double  t 
)
inline

◆ HepLorentzVector() [2/12]

CLHEP::HepLorentzVector::HepLorentzVector ( double  x,
double  y,
double  z 
)
inline

◆ HepLorentzVector() [3/12]

CLHEP::HepLorentzVector::HepLorentzVector ( double  t)
explicit

◆ HepLorentzVector() [4/12]

CLHEP::HepLorentzVector::HepLorentzVector ( )
inline

Referenced by rest4Vector().

◆ HepLorentzVector() [5/12]

CLHEP::HepLorentzVector::HepLorentzVector ( const Hep3Vector p,
double  e 
)
inline

◆ HepLorentzVector() [6/12]

CLHEP::HepLorentzVector::HepLorentzVector ( double  e,
const Hep3Vector p 
)
inline

◆ HepLorentzVector() [7/12]

CLHEP::HepLorentzVector::HepLorentzVector ( const HepLorentzVector )
inline

◆ HepLorentzVector() [8/12]

CLHEP::HepLorentzVector::HepLorentzVector ( HepLorentzVector &&  )
inlinedefault

◆ ~HepLorentzVector()

CLHEP::HepLorentzVector::~HepLorentzVector ( )
inline

◆ HepLorentzVector() [9/12]

CLHEP::HepLorentzVector::HepLorentzVector ( double  x,
double  y,
double  z,
Tcomponent  t 
)
inline

◆ HepLorentzVector() [10/12]

CLHEP::HepLorentzVector::HepLorentzVector ( Tcomponent  t,
double  x,
double  y,
double  z 
)
inline

◆ HepLorentzVector() [11/12]

CLHEP::HepLorentzVector::HepLorentzVector ( Tcomponent  t)
inlineexplicit

◆ HepLorentzVector() [12/12]

CLHEP::HepLorentzVector::HepLorentzVector ( const Hep3Vector v)
inlineexplicit

Member Function Documentation

◆ angle()

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

◆ beta()

double CLHEP::HepLorentzVector::beta ( ) const

Definition at line 72 of file LorentzVectorK.cc.

72 {
73 if (ee == 0) {
74 if (pp.mag2() == 0) {
75 return 0;
76 } else {
77 ZMthrowA (ZMxpvInfiniteVector(
78 "beta computed for HepLorentzVector with t=0 -- infinite result"));
79 return 1./ee;
80 }
81 }
82 if (restMass2() <= 0) {
83 ZMthrowC (ZMxpvTachyonic(
84 "beta computed for a non-timelike HepLorentzVector"));
85 // result will make analytic sense but is physically meaningless
86 }
87 return std::sqrt (pp.mag2() / (ee*ee)) ;
88} /* beta */
#define ZMthrowC(A)
Definition: ZMxpv.h:133
#define ZMthrowA(A)
Definition: ZMxpv.h:128
double mag2() const
double restMass2() const

◆ boost() [1/3]

HepLorentzVector & CLHEP::HepLorentzVector::boost ( const Hep3Vector )
inline

◆ boost() [2/3]

HepLorentzVector & CLHEP::HepLorentzVector::boost ( const Hep3Vector axis,
double  beta 
)

Definition at line 49 of file LorentzVectorB.cc.

50 {
51 if (bbeta==0) {
52 return *this; // do nothing for a 0 boost
53 }
54 double r2 = aaxis.mag2();
55 if ( r2 == 0 ) {
56 ZMthrowA (ZMxpvZeroVector(
57 "A zero vector used as axis defining a boost -- no boost done"));
58 return *this;
59 }
60 double b2 = bbeta*bbeta;
61 if (b2 >= 1) {
62 ZMthrowA (ZMxpvTachyonic(
63 "LorentzVector boosted with beta >= 1 (speed of light) -- \n"
64 "no boost done"));
65 } else {
66 Hep3Vector u = aaxis.unit();
67 double ggamma = std::sqrt(1./(1.-b2));
68 double betaDotV = u.dot(pp)*bbeta;
69 double tt = ee;
70
71 ee = ggamma * (tt + betaDotV);
72 pp += ( ((ggamma-1)/b2)*betaDotV*bbeta + ggamma*bbeta*tt ) * u;
73 // Note: I have verified the behavior of this even when beta is very
74 // small -- (gamma-1)/b2 becomes inaccurate by O(1), but it is then
75 // multiplied by O(beta**2) and added to an O(beta) term, so the
76 // inaccuracy does not affect the final result.
77 }
78 return *this;
79} /* boost (axis, beta) */

◆ boost() [3/3]

HepLorentzVector & CLHEP::HepLorentzVector::boost ( double  bx,
double  by,
double  bz 
)

Definition at line 53 of file LorentzVector.cc.

54 {
55 double b2 = bx*bx + by*by + bz*bz;
56 double ggamma = 1.0 / std::sqrt(1.0 - b2);
57 double bp = bx*x() + by*y() + bz*z();
58 double gamma2 = b2 > 0 ? (ggamma - 1.0)/b2 : 0.0;
59
60 setX(x() + gamma2*bp*bx + ggamma*bx*t());
61 setY(y() + gamma2*bp*by + ggamma*by*t());
62 setZ(z() + gamma2*bp*bz + ggamma*bz*t());
63 setT(ggamma*(t() + bp));
64 return *this;
65}

Referenced by main().

◆ boostVector()

Hep3Vector CLHEP::HepLorentzVector::boostVector ( ) const

Definition at line 170 of file LorentzVector.cc.

170 {
171 if (ee == 0) {
172 if (pp.mag2() == 0) {
173 return Hep3Vector(0,0,0);
174 } else {
175 ZMthrowA (ZMxpvInfiniteVector(
176 "boostVector computed for LorentzVector with t=0 -- infinite result"));
177 return pp/ee;
178 }
179 }
180 if (restMass2() <= 0) {
181 ZMthrowC (ZMxpvTachyonic(
182 "boostVector computed for a non-timelike LorentzVector "));
183 // result will make analytic sense but is physically meaningless
184 }
185 return pp * (1./ee);
186} /* boostVector */

Referenced by findBoostToCM().

◆ boostX()

HepLorentzVector & CLHEP::HepLorentzVector::boostX ( double  beta)

Definition at line 189 of file LorentzVector.cc.

189 {
190 double b2 = bbeta*bbeta;
191 if (b2 >= 1) {
192 ZMthrowA (ZMxpvTachyonic(
193 "boost along X with beta >= 1 (speed of light) -- no boost done"));
194 } else {
195 double ggamma = std::sqrt(1./(1-b2));
196 double tt = ee;
197 ee = ggamma*(ee + bbeta*pp.getX());
198 pp.setX(ggamma*(pp.getX() + bbeta*tt));
199 }
200 return *this;
201} /* boostX */
double getX() const
void setX(double)

◆ boostY()

HepLorentzVector & CLHEP::HepLorentzVector::boostY ( double  beta)

Definition at line 203 of file LorentzVector.cc.

203 {
204 double b2 = bbeta*bbeta;
205 if (b2 >= 1) {
206 ZMthrowA (ZMxpvTachyonic(
207 "boost along Y with beta >= 1 (speed of light) -- \nno boost done"));
208 } else {
209 double ggamma = std::sqrt(1./(1-b2));
210 double tt = ee;
211 ee = ggamma*(ee + bbeta*pp.getY());
212 pp.setY(ggamma*(pp.getY() + bbeta*tt));
213 }
214 return *this;
215} /* boostY */
void setY(double)
double getY() const

◆ boostZ()

HepLorentzVector & CLHEP::HepLorentzVector::boostZ ( double  beta)

Definition at line 217 of file LorentzVector.cc.

217 {
218 double b2 = bbeta*bbeta;
219 if (b2 >= 1) {
220 ZMthrowA (ZMxpvTachyonic(
221 "boost along Z with beta >= 1 (speed of light) -- \nno boost done"));
222 } else {
223 double ggamma = std::sqrt(1./(1-b2));
224 double tt = ee;
225 ee = ggamma*(ee + bbeta*pp.getZ());
226 pp.setZ(ggamma*(pp.getZ() + bbeta*tt));
227 }
228 return *this;
229} /* boostZ */
double getZ() const
void setZ(double)

◆ coLinearRapidity()

double CLHEP::HepLorentzVector::coLinearRapidity ( ) const

Definition at line 159 of file LorentzVectorK.cc.

159 {
160 double v1 = pp.mag();
161 if (std::fabs(ee) == std::fabs(v1)) {
162 ZMthrowA (ZMxpvInfinity(
163 "co-Linear rapidity for 4-vector with |E| = |P| -- infinite result"));
164 }
165 if (std::fabs(ee) < std::fabs(v1)) {
166 ZMthrowA (ZMxpvSpacelike(
167 "co-linear rapidity for spacelike 4-vector -- undefined"));
168 return 0;
169 }
170 double q = (ee + v1) / (ee - v1);
171 return .5 * std::log(q);
172} /* rapidity */
double mag() const

◆ compare()

int CLHEP::HepLorentzVector::compare ( const HepLorentzVector w) const

Definition at line 27 of file LorentzVectorC.cc.

27 {
28 if ( ee > w.ee ) {
29 return 1;
30 } else if ( ee < w.ee ) {
31 return -1;
32 } else {
33 return ( pp.compare(w.pp) );
34 }
35} /* Compare */
int compare(const Hep3Vector &v) const
Definition: SpaceVector.cc:121

Referenced by operator<(), operator<=(), operator>(), and operator>=().

◆ cosTheta()

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

◆ delta2Euclidean()

double CLHEP::HepLorentzVector::delta2Euclidean ( const HepLorentzVector w) const
inline

◆ deltaR()

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

Definition at line 193 of file LorentzVectorC.cc.

193 {
194
195 double a = eta() - w.eta();
196 double b = pp.deltaPhi(w.getV());
197
198 return std::sqrt ( a*a + b*b );
199
200} /* deltaR */
double deltaPhi(const Hep3Vector &v2) const
Definition: ThreeVector.cc:138

◆ diff2()

double CLHEP::HepLorentzVector::diff2 ( const HepLorentzVector w) const
inline

◆ dot()

double CLHEP::HepLorentzVector::dot ( const HepLorentzVector ) const
inline

◆ e()

double CLHEP::HepLorentzVector::e ( ) const
inline

Referenced by operator()().

◆ et() [1/2]

double CLHEP::HepLorentzVector::et ( ) const
inline

◆ et() [2/2]

double CLHEP::HepLorentzVector::et ( const Hep3Vector ) const
inline

◆ et2() [1/2]

double CLHEP::HepLorentzVector::et2 ( ) const
inline

◆ et2() [2/2]

double CLHEP::HepLorentzVector::et2 ( const Hep3Vector ) const
inline

◆ eta() [1/2]

double CLHEP::HepLorentzVector::eta ( ) const
inline

Referenced by deltaR().

◆ eta() [2/2]

double CLHEP::HepLorentzVector::eta ( const Hep3Vector ref) const
inline

◆ euclideanNorm()

double CLHEP::HepLorentzVector::euclideanNorm ( ) const
inline

Referenced by howParallel(), and isParallel().

◆ euclideanNorm2()

double CLHEP::HepLorentzVector::euclideanNorm2 ( ) const
inline

Referenced by isParallel().

◆ findBoostToCM() [1/2]

Hep3Vector CLHEP::HepLorentzVector::findBoostToCM ( ) const

Definition at line 208 of file LorentzVectorK.cc.

208 {
209 return -boostVector();
210} /* boostToCM() */
Hep3Vector boostVector() const

◆ findBoostToCM() [2/2]

Hep3Vector CLHEP::HepLorentzVector::findBoostToCM ( const HepLorentzVector w) const

Definition at line 212 of file LorentzVectorK.cc.

212 {
213 double t1 = ee + w.ee;
214 Hep3Vector v1 = pp + w.pp;
215 if (t1 == 0) {
216 if (v1.mag2() == 0) {
217 return Hep3Vector(0,0,0);
218 } else {
219 ZMthrowA (ZMxpvInfiniteVector(
220 "boostToCM computed for two 4-vectors with combined t=0 -- "
221 "infinite result"));
222 return Hep3Vector(v1*(1./t1)); // Yup, 1/0 -- that is how we return infinity
223 }
224 }
225 if (t1*t1 - v1.mag2() <= 0) {
226 ZMthrowC (ZMxpvTachyonic(
227 "boostToCM computed for pair of HepLorentzVectors with non-timelike sum"));
228 // result will make analytic sense but is physically meaningless
229 }
230 return Hep3Vector(v1 * (-1./t1));
231} /* boostToCM(w) */

◆ gamma()

double CLHEP::HepLorentzVector::gamma ( ) const

Definition at line 90 of file LorentzVectorK.cc.

90 {
91 double v2 = pp.mag2();
92 double t2 = ee*ee;
93 if (ee == 0) {
94 if (pp.mag2() == 0) {
95 return 1;
96 } else {
97 ZMthrowC (ZMxpvInfiniteVector(
98 "gamma computed for HepLorentzVector with t=0 -- zero result"));
99 return 0;
100 }
101 }
102 if (t2 < v2) {
103 ZMthrowA (ZMxpvSpacelike(
104 "gamma computed for a spacelike HepLorentzVector -- imaginary result"));
105 // analytic result would be imaginary.
106 return 0;
107 } else if ( t2 == v2 ) {
108 ZMthrowA (ZMxpvInfinity(
109 "gamma computed for a lightlike HepLorentzVector -- infinite result"));
110 }
111 return 1./std::sqrt(1. - v2/t2 );
112} /* gamma */

◆ getMetric()

ZMpvMetric_t CLHEP::HepLorentzVector::getMetric ( )
static

Definition at line 34 of file LorentzVectorK.cc.

34 {
35 return ( (metric > 0) ? TimePositive : TimeNegative );
36}
@ TimePositive
Definition: LorentzVector.h:61
@ TimeNegative
Definition: LorentzVector.h:61

◆ getT()

double CLHEP::HepLorentzVector::getT ( ) const
inline

◆ getTolerance()

double CLHEP::HepLorentzVector::getTolerance ( )
static

Definition at line 238 of file LorentzVector.cc.

238 {
239// Get the tolerance for two LorentzVectors to be considered near each other
240 return tolerance;
241}

◆ getV()

Hep3Vector CLHEP::HepLorentzVector::getV ( ) const
inline

Referenced by deltaR(), and CLHEP::operator/().

◆ getX()

double CLHEP::HepLorentzVector::getX ( ) const
inline

◆ getY()

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

◆ getZ()

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

◆ howLightlike()

double CLHEP::HepLorentzVector::howLightlike ( ) const

Definition at line 247 of file LorentzVectorC.cc.

247 {
248 double m1 = std::fabs(restMass2());
249 double twoT2 = 2*ee*ee;
250 if (m1 < twoT2) {
251 return m1/twoT2;
252 } else {
253 return 1;
254 }
255} /* HowLightlike */

◆ howNear()

double CLHEP::HepLorentzVector::howNear ( const HepLorentzVector w) const

Definition at line 65 of file LorentzVectorC.cc.

65 {
66 double wdw = std::fabs(pp.dot(w.pp)) + .25*((ee+w.ee)*(ee+w.ee));
67 double delta = (pp - w.pp).mag2() + (ee-w.ee)*(ee-w.ee);
68 if ( (wdw > 0) && (delta < wdw) ) {
69 return std::sqrt (delta/wdw);
70 } else if ( (wdw == 0) && (delta == 0) ) {
71 return 0;
72 } else {
73 return 1;
74 }
75} /* howNear() */
double dot(const Hep3Vector &) const

Referenced by howNearCM().

◆ howNearCM()

double CLHEP::HepLorentzVector::howNearCM ( const HepLorentzVector w) const

Definition at line 131 of file LorentzVectorC.cc.

131 {
132
133 double tTotal = (ee + w.ee);
134 Hep3Vector vTotal (pp + w.pp);
135 double vTotal2 = vTotal.mag2();
136
137 if ( vTotal2 >= tTotal*tTotal ) {
138 // Either one or both vectors are spacelike, or the dominant T components
139 // are in opposite directions. So boosting and testing makes no sense;
140 // but we do consider two exactly equal vectors to be equal in any frame,
141 // even if they are spacelike and can't be boosted to a CM frame.
142 if (*this == w) {
143 return 0;
144 } else {
145 return 1;
146 }
147 }
148
149 if ( vTotal2 == 0 ) { // no boost needed!
150 return (howNear(w));
151 }
152
153 // Find the boost to the CM frame. We know that the total vector is timelike.
154
155 double tRecip = 1./tTotal;
156 Hep3Vector bboost ( vTotal * (-tRecip) );
157
158 //-| Note that you could do pp/t and not be terribly inefficient since
159 //-| SpaceVector/t itself takes 1/t and multiplies. The code here saves
160 //-| a redundant check for t=0.
161
162 // Boost both vectors. Since we have the same boost, there is no need
163 // to repeat the beta and gamma calculation; and there is no question
164 // about beta >= 1. That is why we don't just call w.boosted().
165
166 double b2 = vTotal2*tRecip*tRecip;
167 if ( b2 >= 1 ) { // NaN-proofing
168 ZMthrowC ( ZMxpvTachyonic (
169 "boost vector in howNearCM appears to be tachyonic"));
170 }
171 double ggamma = std::sqrt(1./(1.-b2));
172 double boostDotV1 = bboost.dot(pp);
173 double gm1_b2 = (ggamma-1)/b2;
174
175 HepLorentzVector w1 ( pp + ((gm1_b2)*boostDotV1+ggamma*ee) * bboost,
176 ggamma * (ee + boostDotV1) );
177
178 double boostDotV2 = bboost.dot(w.pp);
179 HepLorentzVector w2 ( w.pp + ((gm1_b2)*boostDotV2+ggamma*w.ee) * bboost,
180 ggamma * (w.ee + boostDotV2) );
181
182 return (w1.howNear(w2));
183
184} /* howNearCM() */
double howNear(const HepLorentzVector &w) const

◆ howParallel()

double CLHEP::HepLorentzVector::howParallel ( const HepLorentzVector w) const

Definition at line 225 of file LorentzVectorC.cc.

225 {
226
227 double norm = euclideanNorm();
228 double wnorm = w.euclideanNorm();
229 if ( norm == 0 ) {
230 if ( wnorm == 0 ) {
231 return 0;
232 } else {
233 return 1;
234 }
235 }
236 if ( wnorm == 0 ) {
237 return 1;
238 }
239
240 HepLorentzVector w1 = *this / norm;
241 HepLorentzVector w2 = w / wnorm;
242 double x1 = (w1-w2).euclideanNorm();
243 return (x1 < 1) ? x1 : 1;
244
245} /* howParallel */
double euclideanNorm() const
double norm(const HepGenMatrix &m)
Definition: GenMatrix.cc:54

◆ invariantMass() [1/2]

double CLHEP::HepLorentzVector::invariantMass ( ) const
inline

◆ invariantMass() [2/2]

double CLHEP::HepLorentzVector::invariantMass ( const HepLorentzVector w) const

Definition at line 178 of file LorentzVectorK.cc.

178 {
179 double m1 = invariantMass2(w);
180 if (m1 < 0) {
181 // We should find out why:
182 if ( ee * w.ee < 0 ) {
183 ZMthrowA (ZMxpvNegativeMass(
184 "invariant mass meaningless: \n"
185 "a negative-mass input led to spacelike 4-vector sum" ));
186 return 0;
187 } else if ( (isSpacelike() && !isLightlike()) ||
188 (w.isSpacelike() && !w.isLightlike()) ) {
189 ZMthrowA (ZMxpvSpacelike(
190 "invariant mass meaningless because of spacelike input"));
191 return 0;
192 } else {
193 // Invariant mass squared for a pair of timelike or lightlike vectors
194 // mathematically cannot be negative. If the vectors are within the
195 // tolerance of being lightlike or timelike, we can assume that prior
196 // or current roundoffs have caused the negative result, and return 0
197 // without comment.
198 return 0;
199 }
200 }
201 return (ee+w.ee >=0 ) ? std::sqrt(m1) : - std::sqrt(m1);
202} /* invariantMass */
bool isLightlike(double epsilon=tolerance) const
double invariantMass2() const
bool isSpacelike() const

◆ invariantMass2() [1/2]

double CLHEP::HepLorentzVector::invariantMass2 ( ) const
inline

Referenced by invariantMass().

◆ invariantMass2() [2/2]

double CLHEP::HepLorentzVector::invariantMass2 ( const HepLorentzVector w) const
inline

◆ isLightlike()

bool CLHEP::HepLorentzVector::isLightlike ( double  epsilon = tolerance) const
inline

Referenced by invariantMass().

◆ isNear()

bool CLHEP::HepLorentzVector::isNear ( const HepLorentzVector w,
double  epsilon = tolerance 
) const

Definition at line 55 of file LorentzVectorC.cc.

56 {
57 double limit = std::fabs(pp.dot(w.pp));
58 limit += .25*((ee+w.ee)*(ee+w.ee));
59 limit *= epsilon*epsilon;
60 double delta = (pp - w.pp).mag2();
61 delta += (ee-w.ee)*(ee-w.ee);
62 return (delta <= limit );
63} /* isNear() */

Referenced by isNearCM().

◆ isNearCM()

bool CLHEP::HepLorentzVector::isNearCM ( const HepLorentzVector w,
double  epsilon = tolerance 
) const

Definition at line 82 of file LorentzVectorC.cc.

83 {
84
85 double tTotal = (ee + w.ee);
86 Hep3Vector vTotal (pp + w.pp);
87 double vTotal2 = vTotal.mag2();
88
89 if ( vTotal2 >= tTotal*tTotal ) {
90 // Either one or both vectors are spacelike, or the dominant T components
91 // are in opposite directions. So boosting and testing makes no sense;
92 // but we do consider two exactly equal vectors to be equal in any frame,
93 // even if they are spacelike and can't be boosted to a CM frame.
94 return (*this == w);
95 }
96
97 if ( vTotal2 == 0 ) { // no boost needed!
98 return (isNear(w, epsilon));
99 }
100
101 // Find the boost to the CM frame. We know that the total vector is timelike.
102
103 double tRecip = 1./tTotal;
104 Hep3Vector bboost ( vTotal * (-tRecip) );
105
106 //-| Note that you could do pp/t and not be terribly inefficient since
107 //-| SpaceVector/t itself takes 1/t and multiplies. The code here saves
108 //-| a redundant check for t=0.
109
110 // Boost both vectors. Since we have the same boost, there is no need
111 // to repeat the beta and gamma calculation; and there is no question
112 // about beta >= 1. That is why we don't just call w.boosted().
113
114 double b2 = vTotal2*tRecip*tRecip;
115
116 double ggamma = std::sqrt(1./(1.-b2));
117 double boostDotV1 = bboost.dot(pp);
118 double gm1_b2 = (ggamma-1)/b2;
119
120 HepLorentzVector w1 ( pp + ((gm1_b2)*boostDotV1+ggamma*ee) * bboost,
121 ggamma * (ee + boostDotV1) );
122
123 double boostDotV2 = bboost.dot(w.pp);
124 HepLorentzVector w2 ( w.pp + ((gm1_b2)*boostDotV2+ggamma*w.ee) * bboost,
125 ggamma * (w.ee + boostDotV2) );
126
127 return (w1.isNear(w2, epsilon));
128
129} /* isNearCM() */
bool isNear(const HepLorentzVector &w, double epsilon=tolerance) const

◆ isParallel()

bool CLHEP::HepLorentzVector::isParallel ( const HepLorentzVector w,
double  epsilon = tolerance 
) const

Definition at line 206 of file LorentzVectorC.cc.

206 {
207 double norm = euclideanNorm();
208 double wnorm = w.euclideanNorm();
209 if ( norm == 0 ) {
210 if ( wnorm == 0 ) {
211 return true;
212 } else {
213 return false;
214 }
215 }
216 if ( wnorm == 0 ) {
217 return false;
218 }
219 HepLorentzVector w1 = *this / norm;
220 HepLorentzVector w2 = w / wnorm;
221 return ( (w1-w2).euclideanNorm2() <= epsilon*epsilon );
222} /* isParallel */
double euclideanNorm2() const

◆ isSpacelike()

bool CLHEP::HepLorentzVector::isSpacelike ( ) const
inline

Referenced by invariantMass().

◆ isTimelike()

bool CLHEP::HepLorentzVector::isTimelike ( ) const
inline

◆ m()

double CLHEP::HepLorentzVector::m ( ) const
inline

Referenced by rest4Vector().

◆ m2()

double CLHEP::HepLorentzVector::m2 ( ) const
inline

◆ mag()

double CLHEP::HepLorentzVector::mag ( ) const
inline

Referenced by main().

◆ mag2()

double CLHEP::HepLorentzVector::mag2 ( ) const
inline

Referenced by howNear(), isNear(), and main().

◆ minus() [1/2]

double CLHEP::HepLorentzVector::minus ( ) const
inline

◆ minus() [2/2]

double CLHEP::HepLorentzVector::minus ( const Hep3Vector ref) const

Definition at line 53 of file LorentzVectorK.cc.

53 {
54 double r = ref.mag();
55 if (r == 0) {
56 ZMthrowA (ZMxpvZeroVector(
57 "A zero vector used as reference to LorentzVector minus-part"));
58 return ee;
59 }
60 return ee - pp.dot(ref)/r;
61} /* plus */

◆ mt()

double CLHEP::HepLorentzVector::mt ( ) const
inline

◆ mt2()

double CLHEP::HepLorentzVector::mt2 ( ) const
inline

◆ operator const Hep3Vector &()

CLHEP::HepLorentzVector::operator const Hep3Vector & ( ) const
inline

◆ operator Hep3Vector &()

CLHEP::HepLorentzVector::operator Hep3Vector & ( )
inline

◆ operator!=()

bool CLHEP::HepLorentzVector::operator!= ( const HepLorentzVector ) const
inline

◆ operator()() [1/2]

double & CLHEP::HepLorentzVector::operator() ( int  i)

Definition at line 36 of file LorentzVector.cc.

36 {
37 static double dummy;
38 switch(i) {
39 case X:
40 case Y:
41 case Z:
42 return pp(i);
43 case T:
44 return ee;
45 default:
46 std::cerr
47 << "HepLorentzVector subscripting: bad index (" << i << ")"
48 << std::endl;
49 return dummy;
50 }
51}

◆ operator()() [2/2]

double CLHEP::HepLorentzVector::operator() ( int  i) const

Definition at line 21 of file LorentzVector.cc.

21 {
22 switch(i) {
23 case X:
24 case Y:
25 case Z:
26 return pp(i);
27 case T:
28 return e();
29 default:
30 std::cerr << "HepLorentzVector subscripting: bad index (" << i << ")"
31 << std::endl;
32 }
33 return 0.;
34}

◆ operator*()

double CLHEP::HepLorentzVector::operator* ( const HepLorentzVector ) const
inline

◆ operator*=() [1/3]

HepLorentzVector & CLHEP::HepLorentzVector::operator*= ( const HepLorentzRotation m1)

Definition at line 18 of file LorentzVectorL.cc.

18 {
19 return *this = m1.vectorMultiplication(*this);
20}

◆ operator*=() [2/3]

HepLorentzVector & CLHEP::HepLorentzVector::operator*= ( const HepRotation )
inline

◆ operator*=() [3/3]

HepLorentzVector & CLHEP::HepLorentzVector::operator*= ( double  )
inline

◆ operator+()

HepLorentzVector CLHEP::HepLorentzVector::operator+ ( const HepLorentzVector ) const
inline

◆ operator+=()

HepLorentzVector & CLHEP::HepLorentzVector::operator+= ( const HepLorentzVector )
inline

◆ operator-() [1/2]

HepLorentzVector CLHEP::HepLorentzVector::operator- ( ) const
inline

◆ operator-() [2/2]

HepLorentzVector CLHEP::HepLorentzVector::operator- ( const HepLorentzVector ) const
inline

◆ operator-=()

HepLorentzVector & CLHEP::HepLorentzVector::operator-= ( const HepLorentzVector )
inline

◆ operator/=()

HepLorentzVector & CLHEP::HepLorentzVector::operator/= ( double  c)

Definition at line 147 of file LorentzVector.cc.

147 {
148 if (c == 0) {
149 ZMthrowA (ZMxpvInfiniteVector(
150 "Attempt to do LorentzVector /= 0 -- \n"
151 "division by zero would produce infinite or NAN components"));
152 }
153 double oneOverC = 1.0/c;
154 pp *= oneOverC;
155 ee *= oneOverC;
156 return *this;
157} /* w /= c */

◆ operator<()

bool CLHEP::HepLorentzVector::operator< ( const HepLorentzVector w) const

Definition at line 40 of file LorentzVectorC.cc.

40 {
41 return (compare(w) < 0);
42}
int compare(const HepLorentzVector &w) const

◆ operator<=()

bool CLHEP::HepLorentzVector::operator<= ( const HepLorentzVector w) const

Definition at line 46 of file LorentzVectorC.cc.

46 {
47 return (compare(w) <= 0);
48}

◆ operator=() [1/3]

HepLorentzVector & CLHEP::HepLorentzVector::operator= ( const Hep3Vector v)
inline

◆ operator=() [2/3]

HepLorentzVector & CLHEP::HepLorentzVector::operator= ( const HepLorentzVector )
inline

◆ operator=() [3/3]

HepLorentzVector & CLHEP::HepLorentzVector::operator= ( HepLorentzVector &&  )
inlinedefault

◆ operator==()

bool CLHEP::HepLorentzVector::operator== ( const HepLorentzVector ) const
inline

◆ operator>()

bool CLHEP::HepLorentzVector::operator> ( const HepLorentzVector w) const

Definition at line 37 of file LorentzVectorC.cc.

37 {
38 return (compare(w) > 0);
39}

◆ operator>=()

bool CLHEP::HepLorentzVector::operator>= ( const HepLorentzVector w) const

Definition at line 43 of file LorentzVectorC.cc.

43 {
44 return (compare(w) >= 0);
45}

◆ operator[]() [1/2]

double & CLHEP::HepLorentzVector::operator[] ( int  )
inline

◆ operator[]() [2/2]

double CLHEP::HepLorentzVector::operator[] ( int  ) const
inline

◆ perp() [1/2]

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

Referenced by main().

◆ perp() [2/2]

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

◆ perp2() [1/2]

double CLHEP::HepLorentzVector::perp2 ( ) const
inline

Referenced by main().

◆ perp2() [2/2]

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

◆ phi()

double CLHEP::HepLorentzVector::phi ( ) const
inline

Referenced by main().

◆ plus() [1/2]

double CLHEP::HepLorentzVector::plus ( ) const
inline

◆ plus() [2/2]

double CLHEP::HepLorentzVector::plus ( const Hep3Vector ref) const

Definition at line 43 of file LorentzVectorK.cc.

43 {
44 double r = ref.mag();
45 if (r == 0) {
46 ZMthrowA (ZMxpvZeroVector(
47 "A zero vector used as reference to LorentzVector plus-part"));
48 return ee;
49 }
50 return ee + pp.dot(ref)/r;
51} /* plus */

◆ pseudoRapidity()

double CLHEP::HepLorentzVector::pseudoRapidity ( ) const
inline

◆ px()

double CLHEP::HepLorentzVector::px ( ) const
inline

◆ py()

double CLHEP::HepLorentzVector::py ( ) const
inline

◆ pz()

double CLHEP::HepLorentzVector::pz ( ) const
inline

◆ rapidity() [1/2]

double CLHEP::HepLorentzVector::rapidity ( ) const

Definition at line 121 of file LorentzVectorK.cc.

121 {
122 double z1 = pp.getZ();
123 if (std::fabs(ee) == std::fabs(z1)) {
124 ZMthrowA (ZMxpvInfinity(
125 "rapidity for 4-vector with |E| = |Pz| -- infinite result"));
126 }
127 if (std::fabs(ee) < std::fabs(z1)) {
128 ZMthrowA (ZMxpvSpacelike(
129 "rapidity for spacelike 4-vector with |E| < |Pz| -- undefined"));
130 return 0;
131 }
132 double q = (ee + z1) / (ee - z1);
133 //-| This cannot be negative now, since both numerator
134 //-| and denominator have the same sign as ee.
135 return .5 * std::log(q);
136} /* rapidity */

◆ rapidity() [2/2]

double CLHEP::HepLorentzVector::rapidity ( const Hep3Vector ref) const

Definition at line 138 of file LorentzVectorK.cc.

138 {
139 double r = ref.mag2();
140 if (r == 0) {
141 ZMthrowA (ZMxpvZeroVector(
142 "A zero vector used as reference to LorentzVector rapidity"));
143 return 0;
144 }
145 double vdotu = pp.dot(ref)/std::sqrt(r);
146 if (std::fabs(ee) == std::fabs(vdotu)) {
147 ZMthrowA (ZMxpvInfinity(
148 "rapidity for 4-vector with |E| = |Pu| -- infinite result"));
149 }
150 if (std::fabs(ee) < std::fabs(vdotu)) {
151 ZMthrowA (ZMxpvSpacelike(
152 "rapidity for spacelike 4-vector with |E| < |P*ref| -- undefined "));
153 return 0;
154 }
155 double q = (ee + vdotu) / (ee - vdotu);
156 return .5 * std::log(q);
157} /* rapidity(ref) */

◆ rest4Vector()

HepLorentzVector CLHEP::HepLorentzVector::rest4Vector ( ) const

Definition at line 63 of file LorentzVectorK.cc.

63 {
64 return HepLorentzVector (0, 0, 0, (t() < 0.0 ? -m() : m()));
65}

◆ restMass()

double CLHEP::HepLorentzVector::restMass ( ) const
inline

◆ restMass2()

double CLHEP::HepLorentzVector::restMass2 ( ) const
inline

Referenced by beta(), boostVector(), and howLightlike().

◆ rho()

double CLHEP::HepLorentzVector::rho ( ) const
inline

◆ rotate() [1/5]

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

Definition at line 21 of file LorentzVectorR.cc.

22 {
23 pp.rotate (aaxis, ddelta);
24 return *this;
25}
Hep3Vector & rotate(double, const Hep3Vector &)
Definition: ThreeVectorR.cc:25

◆ rotate() [2/5]

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

Definition at line 27 of file LorentzVectorR.cc.

27 {
28 pp.rotate (ax);
29 return *this;
30}

◆ rotate() [3/5]

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

Definition at line 32 of file LorentzVectorR.cc.

32 {
33 pp.rotate (e1);
34 return *this;
35}

◆ rotate() [4/5]

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

Definition at line 37 of file LorentzVectorR.cc.

39 {
40 pp.rotate (phi1, theta1, psi1);
41 return *this;
42}

◆ rotate() [5/5]

HepLorentzVector & CLHEP::HepLorentzVector::rotate ( double  a,
const Hep3Vector v1 
)

Definition at line 16 of file LorentzVectorR.cc.

16 {
17 pp.rotate(a,v1);
18 return *this;
19}

◆ rotateUz()

HepLorentzVector & CLHEP::HepLorentzVector::rotateUz ( const Hep3Vector v1)

Definition at line 80 of file LorentzVector.cc.

80 {
81 pp.rotateUz(v1);
82 return *this;
83}
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:36

◆ rotateX()

HepLorentzVector & CLHEP::HepLorentzVector::rotateX ( double  a)

Definition at line 67 of file LorentzVector.cc.

67 {
68 pp.rotateX(a);
69 return *this;
70}
Hep3Vector & rotateX(double)
Definition: ThreeVector.cc:90

◆ rotateY()

HepLorentzVector & CLHEP::HepLorentzVector::rotateY ( double  a)

Definition at line 71 of file LorentzVector.cc.

71 {
72 pp.rotateY(a);
73 return *this;
74}
Hep3Vector & rotateY(double)
Definition: ThreeVector.cc:100

Referenced by main().

◆ rotateZ()

HepLorentzVector & CLHEP::HepLorentzVector::rotateZ ( double  a)

Definition at line 75 of file LorentzVector.cc.

75 {
76 pp.rotateZ(a);
77 return *this;
78}
Hep3Vector & rotateZ(double)
Definition: ThreeVector.cc:110

Referenced by main().

◆ set() [1/8]

void CLHEP::HepLorentzVector::set ( const Hep3Vector v)
inline

◆ set() [2/8]

void CLHEP::HepLorentzVector::set ( const Hep3Vector v,
double  t 
)
inline

◆ set() [3/8]

void CLHEP::HepLorentzVector::set ( double  t)
inline

◆ set() [4/8]

void CLHEP::HepLorentzVector::set ( double  t,
const Hep3Vector v 
)
inline

◆ set() [5/8]

void CLHEP::HepLorentzVector::set ( double  x,
double  y,
double  z,
double  t 
)
inline

◆ set() [6/8]

void CLHEP::HepLorentzVector::set ( double  x,
double  y,
double  z,
Tcomponent  t 
)
inline

◆ set() [7/8]

void CLHEP::HepLorentzVector::set ( Tcomponent  t)
inline

◆ set() [8/8]

void CLHEP::HepLorentzVector::set ( Tcomponent  t,
double  x,
double  y,
double  z 
)
inline

◆ setE()

void CLHEP::HepLorentzVector::setE ( double  )
inline

◆ setMetric()

ZMpvMetric_t CLHEP::HepLorentzVector::setMetric ( ZMpvMetric_t  a1)
static

Definition at line 24 of file LorentzVectorK.cc.

24 {
25 ZMpvMetric_t oldMetric = (metric > 0) ? TimePositive : TimeNegative;
26 if ( a1 == TimeNegative ) {
27 metric = -1.0;
28 } else {
29 metric = 1.0;
30 }
31 return oldMetric;
32}

Referenced by CLHEP::HepLorentzRotation::set().

◆ setPerp()

void CLHEP::HepLorentzVector::setPerp ( double  )
inline

◆ setPhi()

void CLHEP::HepLorentzVector::setPhi ( double  )
inline

◆ setPx()

void CLHEP::HepLorentzVector::setPx ( double  )
inline

◆ setPy()

void CLHEP::HepLorentzVector::setPy ( double  )
inline

◆ setPz()

void CLHEP::HepLorentzVector::setPz ( double  )
inline

◆ setREtaPhi()

void CLHEP::HepLorentzVector::setREtaPhi ( double  r,
double  eta,
double  phi 
)
inline

◆ setRho()

void CLHEP::HepLorentzVector::setRho ( double  )
inline

◆ setRhoPhiZ()

void CLHEP::HepLorentzVector::setRhoPhiZ ( double  rho,
double  phi,
double  z 
)
inline

◆ setRThetaPhi()

void CLHEP::HepLorentzVector::setRThetaPhi ( double  r,
double  theta,
double  phi 
)
inline

◆ setT()

void CLHEP::HepLorentzVector::setT ( double  )
inline

Referenced by boost(), and CLHEP::operator>>().

◆ setTheta()

void CLHEP::HepLorentzVector::setTheta ( double  )
inline

◆ setTolerance()

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

Definition at line 231 of file LorentzVector.cc.

231 {
232// Set the tolerance for two LorentzVectors to be considered near each other
233 double oldTolerance (tolerance);
234 tolerance = tol;
235 return oldTolerance;
236}

◆ setV() [1/2]

void CLHEP::HepLorentzVector::setV ( const Hep3Vector )
inline

◆ setV() [2/2]

void CLHEP::HepLorentzVector::setV ( double  x,
double  y,
double  z 
)
inline

◆ setVect()

void CLHEP::HepLorentzVector::setVect ( const Hep3Vector )
inline

◆ setVectM()

void CLHEP::HepLorentzVector::setVectM ( const Hep3Vector spatial,
double  mass 
)
inline

◆ setVectMag()

void CLHEP::HepLorentzVector::setVectMag ( const Hep3Vector spatial,
double  magnitude 
)
inline

◆ setX()

void CLHEP::HepLorentzVector::setX ( double  )
inline

Referenced by boost(), and CLHEP::operator>>().

◆ setY()

void CLHEP::HepLorentzVector::setY ( double  )
inline

Referenced by boost(), and CLHEP::operator>>().

◆ setZ()

void CLHEP::HepLorentzVector::setZ ( double  )
inline

Referenced by boost(), and CLHEP::operator>>().

◆ t()

◆ theta()

double CLHEP::HepLorentzVector::theta ( ) const
inline

Referenced by main().

◆ transform() [1/2]

HepLorentzVector & CLHEP::HepLorentzVector::transform ( const HepLorentzRotation m1)

Definition at line 23 of file LorentzVectorL.cc.

23 {
24 return *this = m1.vectorMultiplication(*this);
25}

◆ transform() [2/2]

HepLorentzVector & CLHEP::HepLorentzVector::transform ( const HepRotation )
inline

Referenced by main().

◆ v()

Hep3Vector CLHEP::HepLorentzVector::v ( ) const
inline

◆ vect()

Hep3Vector CLHEP::HepLorentzVector::vect ( ) const
inline

Referenced by main().

◆ x()

◆ y()

◆ z()

Friends And Related Function Documentation

◆ boostOf [1/2]

HepLorentzVector boostOf ( const HepLorentzVector vec,
const Hep3Vector axis,
double  beta 
)
friend

◆ boostOf [2/2]

HepLorentzVector boostOf ( const HepLorentzVector vec,
const Hep3Vector betaVector 
)
friend

◆ boostXOf

HepLorentzVector boostXOf ( const HepLorentzVector vec,
double  beta 
)
friend

◆ boostYOf

HepLorentzVector boostYOf ( const HepLorentzVector vec,
double  beta 
)
friend

◆ boostZOf

HepLorentzVector boostZOf ( const HepLorentzVector vec,
double  beta 
)
friend

◆ rotationOf [1/4]

HepLorentzVector rotationOf ( const HepLorentzVector vec,
const Hep3Vector axis,
double  delta 
)
friend

Definition at line 44 of file LorentzVectorR.cc.

46 {
47 HepLorentzVector vv (vec);
48 return vv.rotate (aaxis, ddelta);
49}

◆ rotationOf [2/4]

HepLorentzVector rotationOf ( const HepLorentzVector vec,
const HepAxisAngle ax 
)
friend

Definition at line 51 of file LorentzVectorR.cc.

52 {
53 HepLorentzVector vv (vec);
54 return vv.rotate (ax);
55}

◆ rotationOf [3/4]

HepLorentzVector rotationOf ( const HepLorentzVector vec,
const HepEulerAngles e 
)
friend

Definition at line 57 of file LorentzVectorR.cc.

58 {
59 HepLorentzVector vv (vec);
60 return vv.rotate (e1);
61}

◆ rotationOf [4/4]

HepLorentzVector rotationOf ( const HepLorentzVector vec,
double  phi,
double  theta,
double  psi 
)
friend

Definition at line 63 of file LorentzVectorR.cc.

66 {
67 HepLorentzVector vv (vec);
68 return vv.rotate (phi1, theta1, psi1);
69}

◆ rotationXOf

HepLorentzVector rotationXOf ( const HepLorentzVector vec,
double  delta 
)
friend

Definition at line 27 of file LorentzVectorB.cc.

28 {
29 HepLorentzVector vv (vec);
30 return vv.rotateX (phi);
31}

◆ rotationYOf

HepLorentzVector rotationYOf ( const HepLorentzVector vec,
double  delta 
)
friend

Definition at line 33 of file LorentzVectorB.cc.

34 {
35 HepLorentzVector vv (vec);
36 return vv.rotateY (phi);
37}

◆ rotationZOf

HepLorentzVector rotationZOf ( const HepLorentzVector vec,
double  delta 
)
friend

Definition at line 39 of file LorentzVectorB.cc.

40 {
41 HepLorentzVector vv (vec);
42 return vv.rotateZ (phi);
43}

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