//*CMZ : 2.22/09 12/07/99 18.44.19 by Rene Brun //*CMZ : 2.22/07 02/07/99 08.20.08 by Rene Brun //*CMZ : 2.22/06 23/06/99 08.24.47 by Rene Brun //*CMZ : 2.21/06 15/02/99 09.36.40 by Rene Brun //*-- Author : Pasha Murat , Peter Malzacher 12/02/99 //______________________________________________________________________________ //*-*-*-*-*-*-*-*-*-*-*-*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 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.
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))
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);
v3 = -v1;
v1 = v2+v3;
v1+= v3;
v1 = v2 + v3;
v1-= v3;
if (v1 == v2) {...}
if(v1 != v3) {...}
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).
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)
Double_t pcone = v.Plus();
Double_t mcone = v.Minus();
TLorentzRotation l;
v.Transform(l);
v = l*v; or
v *= l; // Attention v = l*v
*/
//
//*KEEP,TError.
#include "TError.h"
//*KEEP,TLorentzVector,T=C++.
#include "TLorentzVector.h"
//*KEEP,TLorentzRotation,T=C++.
#include "TLorentzRotation.h"
//*KEND.
ClassImp(TLorentzVector)
TLorentzVector::TLorentzVector(Double_t x, Double_t y, Double_t z, Double_t t)
: fX(x), fY(y), fZ(z), fE(t) {}
TLorentzVector::TLorentzVector(Double_t * x0)
: fX(x0[0]), fY(x0[1]), fZ(x0[2]), fE(x0[3]) {}
TLorentzVector::TLorentzVector(Float_t * x0)
: fX(x0[0]), fY(x0[1]), fZ(x0[2]), fE(x0[3]) {}
TLorentzVector::TLorentzVector(const TVector3 & p, Double_t e)
: fX(p.X()), fY(p.Y()), fZ(p.Z()), fE(e) {}
TLorentzVector::TLorentzVector(const TLorentzVector & p)
: fX(p.X()), fY(p.Y()), fZ(p.Z()), fE(p.T()) {}
Double_t TLorentzVector::operator () (int i) const
{
switch(i) {
case kX: return fX;
case kY: return fY;
case kZ: return fZ;
case kT: return fE;
default:
Warning("operator()()", "bad index (%d)",i);
}
return 0.;
}
Double_t & TLorentzVector::operator () (int i)
{
switch(i) {
case kX: return fX;
case kY: return fY;
case kZ: return fZ;
case kT: return fE;
default:
Warning("operator()()", "bad index (%d)",i);
}
Double_t dummy;
Double_t & rdummy = dummy;
return rdummy;
}
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::RotateX(Double_t angle) {
Double_t s = TMath::Sin(angle);
Double_t c = TMath::Cos(angle);
Double_t yy = fY;
fY = c*yy - s*fZ;
fZ = s*yy + c*fZ;
}
void TLorentzVector::RotateY(Double_t angle) {
Double_t s = TMath::Sin(angle);
Double_t c = TMath::Cos(angle);
Double_t zz = fZ;
fZ = c*zz - s*fX;
fX = s*zz + c*fX;
}
void TLorentzVector::RotateZ(Double_t angle) {
Double_t s = TMath::Sin(angle);
Double_t c = TMath::Cos(angle);
Double_t xx = fX;
fX = c*xx - s*fY;
fY = s*xx + c*fY;
}
void TLorentzVector::RotateUz(TVector3 &newUzVector) {
// NewUzVector must be normalized !
Double_t u1 = newUzVector.X();
Double_t u2 = newUzVector.Y();
Double_t u3 = newUzVector.Z();
Double_t up = u1*u1 + u2*u2;
if (up) {
up = TMath::Sqrt(up);
Double_t px = fX, py = fY, pz = fZ;
fX = (u1*u3*px - u2*py + u1*up*pz)/up;
fY = (u2*u3*px + u1*py + u2*up*pz)/up;
fZ = (u3*u3*px - px + u3*up*pz)/up;
}
else if (u3 < 0.) { fX = -fX; fZ = -fZ; } // phi=0 teta=pi
else {};
}