// @(#)root/physics:$Name:  $:$Id: TLorentzVector.cxx,v 1.7 2002/05/18 08:22:00 brun Exp $
// Author: Pasha Murat , Peter Malzacher  12/02/99
//    Oct  8 1999: changed Warning to Error and
//                 return fX in Double_t & operator()
//    Oct 20 1999: dito in Double_t operator()
//    Jan 25 2000: implemented as (fP,fE) instead of (fX,fY,fZ,fE)
//______________________________________________________________________________
//*-*-*-*-*-*-*-*-*-*-*-*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), 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

*/ // #include "TClass.h" #include "TError.h" #include "TLorentzVector.h" #include "TLorentzRotation.h" ClassImp(TLorentzVector) TLorentzVector::TLorentzVector(Double_t x, Double_t y, Double_t z, Double_t t) : fP(x,y,z), fE(t) {} TLorentzVector::TLorentzVector(const Double_t * x0) : fP(x0), fE(x0[3]) {} TLorentzVector::TLorentzVector(const Float_t * x0) : fP(x0), fE(x0[3]) {} TLorentzVector::TLorentzVector(const TVector3 & p, Double_t e) : fP(p), fE(e) {} TLorentzVector::TLorentzVector(const TLorentzVector & p) : TObject(p) , fP(p.Vect()), fE(p.T()) {} TLorentzVector::~TLorentzVector() {} Double_t TLorentzVector::operator () (int i) const { switch(i) { case kX: case kY: case kZ: return fP(i); case kT: return fE; default: Error("operator()()", "bad index (%d) returning 0",i); } return 0.; } Double_t & TLorentzVector::operator () (int i) { switch(i) { case kX: case kY: case kZ: return fP(i); case kT: return fE; default: Error("operator()()", "bad index (%d) returning &fE",i); } return fE; } void TLorentzVector::Boost(Double_t bx, Double_t by, Double_t bz) { Double_t b2 = bx*bx + by*by + bz*bz; register Double_t gamma = 1.0 / TMath::Sqrt(1.0 - b2); register Double_t bp = bx*X() + by*Y() + bz*Z(); register Double_t gamma2 = b2 > 0 ? (gamma - 1.0)/b2 : 0.0; SetX(X() + gamma2*bp*bx + gamma*bx*T()); SetY(Y() + gamma2*bp*by + gamma*by*T()); SetZ(Z() + gamma2*bp*bz + gamma*bz*T()); SetT(gamma*(T() + bp)); } Double_t TLorentzVector::Rapidity() const { return 0.5*log( (E()+Pz()) / (E()-Pz()) ); } TLorentzVector &TLorentzVector::operator *= (const TLorentzRotation & m) { return *this = m.VectorMultiplication(*this); } TLorentzVector &TLorentzVector::Transform(const TLorentzRotation & m) { return *this = m.VectorMultiplication(*this); } void TLorentzVector::Streamer(TBuffer &R__b) { // Stream an object of class TLorentzVector. Double_t x, y, z; UInt_t R__s, R__c; if (R__b.IsReading()) { Version_t R__v = R__b.ReadVersion(&R__s, &R__c); if (R__v > 3) { TLorentzVector::Class()->ReadBuffer(R__b, this, R__v, R__s, R__c); return; } //====process old versions before automatic schema evolution if (R__v != 2) TObject::Streamer(R__b); R__b >> x; R__b >> y; R__b >> z; fP.SetXYZ(x,y,z); R__b >> fE; R__b.CheckByteCount(R__s, R__c, TLorentzVector::IsA()); } else { TLorentzVector::Class()->WriteBuffer(R__b,this); } }



ROOT page - Class index - 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.