library: libPhysics
#include "TLorentzVector.h"

TLorentzVector


class description - source file - inheritance tree (.pdf)

class TLorentzVector : public TObject

Inheritance Chart:
TObject
<-
TLorentzVector

    public:
TLorentzVector(Double_t x = 0.0, Double_t y = 0.0, Double_t z = 0.0, Double_t t = 0.0) TLorentzVector(const Double_t* carray) TLorentzVector(const Float_t* carray) TLorentzVector(const TVector3& vector3, Double_t t) TLorentzVector(const TLorentzVector& lorentzvector) TLorentzVector operator+(const TLorentzVector& q) const TLorentzVector operator-(const TLorentzVector& q) const TLorentzVector operator-() const TLorentzVector operator*(Double_t a) const virtual ~TLorentzVector() Double_t Angle(const TVector3& v) const Double_t Beta() const void Boost(Double_t, Double_t, Double_t) void Boost(const TVector3& b) TVector3 BoostVector() const static TClass* Class() Double_t CosTheta() const Double_t DeltaPhi(const TLorentzVector& v) const Double_t DeltaR(const TLorentzVector& v) const Double_t Dot(const TLorentzVector& q) const Double_t DrEtaPhi(const TLorentzVector& v) const Double_t E() const Double_t Energy() const Double_t Et() const Double_t Et(const TVector3& v) const Double_t Et2() const Double_t Et2(const TVector3& v) const Double_t Eta() const TVector2 EtaPhiVector() Double_t Gamma() const void GetXYZT(Double_t* carray) const void GetXYZT(Float_t* carray) const virtual TClass* IsA() const Double_t M() const Double_t M2() const Double_t Mag() const Double_t Mag2() const Double_t Minus() const Double_t Mt() const Double_t Mt2() const Bool_t operator!=(const TLorentzVector& q) const Double_t operator()(int i) const Double_t& operator()(int i) Double_t operator*(const TLorentzVector& q) const TLorentzVector& operator*=(Double_t a) TLorentzVector& operator*=(const TRotation& m) TLorentzVector& operator*=(const TLorentzRotation&) TLorentzVector& operator+=(const TLorentzVector& q) TLorentzVector& operator-=(const TLorentzVector& q) TLorentzVector& operator=(const TLorentzVector& q) Bool_t operator==(const TLorentzVector& q) const Double_t operator[](int i) const Double_t& operator[](int i) Double_t P() const Double_t Perp() const Double_t Perp(const TVector3& v) const Double_t Perp2() const Double_t Perp2(const TVector3& v) const Double_t Phi() const Double_t Plus() const Double_t PseudoRapidity() const Double_t Pt() const Double_t Pt(const TVector3& v) const Double_t Px() const Double_t Py() const Double_t Pz() const Double_t Rapidity() const Double_t Rho() const void Rotate(Double_t a, const TVector3& v) void RotateUz(TVector3& newUzVector) void RotateX(Double_t angle) void RotateY(Double_t angle) void RotateZ(Double_t angle) void SetE(Double_t a) void SetPerp(Double_t r) void SetPhi(Double_t phi) void SetPtEtaPhiE(Double_t pt, Double_t eta, Double_t phi, Double_t e) void SetPtEtaPhiM(Double_t pt, Double_t eta, Double_t phi, Double_t m) void SetPx(Double_t a) void SetPxPyPzE(Double_t px, Double_t py, Double_t pz, Double_t e) void SetPy(Double_t a) void SetPz(Double_t a) void SetRho(Double_t rho) void SetT(Double_t a) void SetTheta(Double_t theta) void SetVect(const TVector3& vect3) void SetVectM(const TVector3& spatial, Double_t mass) void SetVectMag(const TVector3& spatial, Double_t magnitude) void SetX(Double_t a) void SetXYZM(Double_t x, Double_t y, Double_t z, Double_t m) void SetXYZT(Double_t x, Double_t y, Double_t z, Double_t t) void SetY(Double_t a) void SetZ(Double_t a) virtual void ShowMembers(TMemberInspector& insp, char* parent) virtual void Streamer(TBuffer& b) void StreamerNVirtual(TBuffer& b) Double_t T() const Double_t Theta() const TLorentzVector& Transform(const TRotation& m) TLorentzVector& Transform(const TLorentzRotation&) TVector3 Vect() const Double_t X() const Double_t Y() const Double_t Z() const

Data Members

    private:
TVector3 fP 3 vector component Double_t fE time or energy of (x,y,z,t) or (px,py,pz,e) public:
static const enum TLorentzVector:: kX static const enum TLorentzVector:: kY static const enum TLorentzVector:: kZ static const enum TLorentzVector:: kT static const enum TLorentzVector:: kNUM_COORDINATES static const enum TLorentzVector:: kSIZE

Class Description

*-*-*-*-*-*-*-*-*-*-*-*The Physics Vector package *-*-*-*-*-*-*-*-*-*-*-*
*-*                    ==========================                       *
*-* The Physics Vector package consists of five classes:                *
*-*   - TVector2                                                        *
*-*   - TVector3                                                        *
*-*   - TRotation                                                       *
*-*   - TLorentzVector                                                  *
*-*   - TLorentzRotation                                                *
*-* It is a combination of CLHEPs Vector package written by             *
*-* Leif Lonnblad, Andreas Nilsson and Evgueni Tcherniaev               *
*-* and a ROOT package written by Pasha Murat.                          *
*-* for CLHEP see:  http://wwwinfo.cern.ch/asd/lhc++/clhep/             *
*-* Adaption to ROOT by Peter Malzacher                                 *
*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*

TLorentzVector

TLorentzVector is a general four-vector class, which can be used either for the description of position and time (x,y,z,t) or momentum and energy (px,py,pz,E).
 

Declaration

TLorentzVector has been implemented as a set a TVector3 and a Double_t variable. By default all components are initialized by zero.

  TLorentzVector v1;      // initialized by (0., 0., 0., 0.)
  TLorentzVector v2(1., 1., 1., 1.);
  TLorentzVector v3(v1);
  TLorentzVector v4(TVector3(1., 2., 3.),4.);

For backward compatibility there are two constructors from an Double_t and Float_t  C array.
 

Access to the components

There are two sets of access functions to the components of a LorentzVector: X(), Y(), Z(), T() and Px(), Py(), Pz() and E(). Both sets return the same values but the first set is more relevant for use where TLorentzVector describes a combination of position and time and the second set is more relevant where TLorentzVector describes momentum and energy:

  Double_t xx =v.X();
  ...
  Double_t tt = v.T();

  Double_t px = v.Px();
  ...
  Double_t ee = v.E();

The components of TLorentzVector can also accessed by index:

  xx = v(0);       or     xx = v[0];
  yy = v(1);              yy = v[1];
  zz = v(2);              zz = v[2];
  tt = v(3);              tt = v[3];

You can use the Vect() member function to get the vector component of TLorentzVector:

  TVector3 p = v.Vect();

For setting components also two sets of member functions can be used: SetX(),.., SetPx(),..:
 
  v.SetX(1.);        or    v.SetPx(1.);
  ...                               ...
  v.SetT(1.);              v.SetE(1.);

To set more the one component by one call you can use the SetVect() function for the TVector3 part or SetXYZT(), SetPxPyPzE(). For convenience there is also a SetXYZM():

  v.SetVect(TVector3(1,2,3));
  v.SetXYZT(x,y,z,t);
  v.SetPxPyPzE(px,py,pz,e);
  v.SetXYZM(x,y,z,m);   //   ->  v=(x,y,z,e=Sqrt(x*x+y*y+z*z+m*m))

Vector components in noncartesian coordinate systems

There are a couple of memberfunctions to get and set the TVector3 part of the parameters in
sherical coordinate systems:

  Double_t m, theta, cost, phi, pp, pp2, ppv2, pp2v2;
  m = v.Rho();
  t = v.Theta();
  cost = v.CosTheta();
  phi = v.Phi();

  v.SetRho(10.);
  v.SetTheta(TMath::Pi()*.3);
  v.SetPhi(TMath::Pi());

or get infoormation about the r-coordinate in cylindrical systems:

  Double_t pp, pp2, ppv2, pp2v2;
  pp = v.Perp();         // get transvers component
  pp2 = v.Perp2();       // get transverse component squared
  ppv2 = v.Perp(v1);      // get transvers component with
                         // respect to another vector
  pp2v2 = v.Perp(v1);

for convenience there are two more set functions SetPtEtaPhiE(pt,eta,phi,e); and SetPtEtaPhiM(pt,eta,phi,m);

Arithmetic and comparison operators

The TLorentzVector class provides operators to add, subtract or compare four-vectors:

  v3 = -v1;
  v1 = v2+v3;
  v1+= v3;
  v1 = v2 + v3;
  v1-= v3;

  if (v1 == v2) {...}
  if(v1 != v3) {...}

Magnitude/Invariant mass, beta, gamma, scalar product

The scalar product of two four-vectors is calculated with the (-,-,-,+) metric,
   i.e.   s = v1*v2 = t1*t2-x1*x2-y1*y2-z1*z2
The magnitude squared mag2 of a four-vector is therfore:
          mag2 = v*v = t*t-x*x-y*y-z*z
It mag2 is negative mag = -Sqrt(-mag*mag). The member functions are:

  Double_t s, s2;
  s  = v1.Dot(v2);     // scalar product
  s  = v1*v2;         // scalar product
  s2 = v.Mag2();   or    s2 = v.M2();
  s  = v.Mag();          s  = v.M();

Since in case of momentum and energy the magnitude has the meaning of invariant mass TLorentzVector provides the more meaningful aliases M2() and M();

The member functions Beta() and Gamma() returns beta and gamma = 1/Sqrt(1-beta*beta).

Lorentz boost

A boost in a general direction can be parameterized with three parameters which can be taken as the components of a three vector b = (bx,by,bz). With
  x = (x,y,z) and gamma = 1/Sqrt(1-beta*beta) (beta being the module of vector b), an arbitary active Lorentz boost transformation (from the rod frame to the original frame) can be written as:
          x = x' + (gamma-1)/(beta*beta) * (b*x') * b + gamma * t' * b
          t = gamma (t'+ b*x').

The member function Boost() performs a boost transformation from the rod frame to the original frame. BoostVector() returns a TVector3 of the spatial components divided by the time component:

  TVector3 b;
  v.Boost(bx,by,bz);
  v.Boost(b);
  b = v.BoostVector();   // b=(x/t,y/t,z/t)

Rotations

There are four sets of functions to rotate the TVector3 component of a TLorentzVector:
rotation around axes
  v.RotateX(TMath::Pi()/2.);
  v.RotateY(.5);
  v.RotateZ(.99);
rotation around an arbitary axis
  v.Rotate(TMath::Pi()/4., v1); // rotation around v1
transformation from rotated frame
  v.RotateUz(direction); //  direction must be a unit TVector3
by TRotation (see TRotation)
  TRotation r;
  v.Transform(r);    or     v *= r; // Attention v=M*v

Misc

Angle between two vectors
  Double_t a = v1.Angle(v2);  // get angle between v1 and v2
Light-cone components
Member functions Plus() and Minus() return the positive and negative light-cone components:

  Double_t pcone = v.Plus();
  Double_t mcone = v.Minus();

CAVEAT: The values returned are T{+,-}Z. It is known that some authors find it easier to define these components as (T{+,-}Z)/sqrt(2). Thus check what definition is used in the physics you're working in and adapt your code accordingly.

Transformation by TLorentzRotation
A general Lorentz transformation see class TLorentzRotation can be used by the Transform() member function, the *= or * operator of the TLorentzRotation class:

  TLorentzRotation l;
  v.Transform(l);
  v = l*v;     or     v *= l;  // Attention v = l*v


TLorentzVector(Double_t x, Double_t y, Double_t z, Double_t t) : fP(x,y,z), fE(t)

TLorentzVector(const Double_t * x0) : fP(x0), fE(x0[3])

TLorentzVector(const Float_t * x0) : fP(x0), fE(x0[3])

TLorentzVector(const TVector3 & p, Double_t e) : fP(p), fE(e)

TLorentzVector(const TLorentzVector & p) : TObject(p) , fP(p.Vect()), fE(p.T())

~TLorentzVector()

void Boost(Double_t bx, Double_t by, Double_t bz)

Double_t Rapidity() const

TLorentzVector& Transform(const TLorentzRotation & m)

void Streamer(TBuffer &R__b)
 Stream an object of class TLorentzVector.



Inline Functions


               Double_t X() const
               Double_t Y() const
               Double_t Z() const
               Double_t T() const
                   void SetX(Double_t a)
                   void SetY(Double_t a)
                   void SetZ(Double_t a)
                   void SetT(Double_t a)
               Double_t Px() const
               Double_t Py() const
               Double_t Pz() const
               Double_t P() const
               Double_t E() const
               Double_t Energy() const
                   void SetPx(Double_t a)
                   void SetPy(Double_t a)
                   void SetPz(Double_t a)
                   void SetE(Double_t a)
               TVector3 Vect() const
                   void SetVect(const TVector3& vect3)
               Double_t Theta() const
               Double_t CosTheta() const
               Double_t Phi() const
               Double_t Rho() const
                   void SetTheta(Double_t theta)
                   void SetPhi(Double_t phi)
                   void SetRho(Double_t rho)
                   void SetPxPyPzE(Double_t px, Double_t py, Double_t pz, Double_t e)
                   void SetXYZT(Double_t x, Double_t y, Double_t z, Double_t t)
                   void SetXYZM(Double_t x, Double_t y, Double_t z, Double_t m)
                   void SetPtEtaPhiM(Double_t pt, Double_t eta, Double_t phi, Double_t m)
                   void SetPtEtaPhiE(Double_t pt, Double_t eta, Double_t phi, Double_t e)
                   void GetXYZT(Double_t* carray) const
                   void GetXYZT(Float_t* carray) const
               Double_t operator()(int i) const
               Double_t operator[](int i) const
              Double_t& operator()(int i)
              Double_t& operator[](int i)
        TLorentzVector& operator=(const TLorentzVector& q)
         TLorentzVector operator+(const TLorentzVector& q) const
        TLorentzVector& operator+=(const TLorentzVector& q)
         TLorentzVector operator-(const TLorentzVector& q) const
        TLorentzVector& operator-=(const TLorentzVector& q)
         TLorentzVector operator-() const
         TLorentzVector operator*(Double_t a) const
        TLorentzVector& operator*=(Double_t a)
                 Bool_t operator==(const TLorentzVector& q) const
                 Bool_t operator!=(const TLorentzVector& q) const
               Double_t Perp2() const
               Double_t Pt() const
               Double_t Perp() const
                   void SetPerp(Double_t r)
               Double_t Perp2(const TVector3& v) const
               Double_t Pt(const TVector3& v) const
               Double_t Perp(const TVector3& v) const
               Double_t Et2() const
               Double_t Et() const
               Double_t Et2(const TVector3& v) const
               Double_t Et(const TVector3& v) const
               Double_t DeltaPhi(const TLorentzVector& v) const
               Double_t DeltaR(const TLorentzVector& v) const
               Double_t DrEtaPhi(const TLorentzVector& v) const
               TVector2 EtaPhiVector()
               Double_t Angle(const TVector3& v) const
               Double_t Mag2() const
               Double_t M2() const
               Double_t Mag() const
               Double_t M() const
               Double_t Mt2() const
               Double_t Mt() const
               Double_t Beta() const
               Double_t Gamma() const
               Double_t Dot(const TLorentzVector& q) const
               Double_t operator*(const TLorentzVector& q) const
                   void SetVectMag(const TVector3& spatial, Double_t magnitude)
                   void SetVectM(const TVector3& spatial, Double_t mass)
               Double_t Plus() const
               Double_t Minus() const
               TVector3 BoostVector() const
                   void Boost(const TVector3& b)
               Double_t Eta() const
               Double_t PseudoRapidity() const
                   void RotateX(Double_t angle)
                   void RotateY(Double_t angle)
                   void RotateZ(Double_t angle)
                   void RotateUz(TVector3& newUzVector)
                   void Rotate(Double_t a, const TVector3& v)
        TLorentzVector& operator*=(const TRotation& m)
        TLorentzVector& operator*=(const TLorentzRotation&)
        TLorentzVector& Transform(const TLorentzRotation&)
                TClass* Class()
                TClass* IsA() const
                   void ShowMembers(TMemberInspector& insp, char* parent)
                   void StreamerNVirtual(TBuffer& b)


Author: Pasha Murat, Peter Malzacher 12/02/99
Last update: root/physics:$Name: $:$Id: TLorentzVector.cxx,v 1.9 2004/04/20 09:29:56 brun Exp $


ROOT page - Class index - Class Hierarchy - Top of the page

This page has been automatically generated. If you have any comments or suggestions about the page layout send a mail to ROOT support, or contact the developers with any questions or problems regarding ROOT.