#ifndef ROOT_Math_GenVector_PtEtaPhiM4D
#define ROOT_Math_GenVector_PtEtaPhiM4D 1
#ifndef ROOT_Math_Math
#include "Math/Math.h"
#endif
#ifndef ROOT_Math_GenVector_etaMax
#include "Math/GenVector/etaMax.h"
#endif
#ifndef ROOT_Math_GenVector_GenVector_exception
#include "Math/GenVector/GenVector_exception.h"
#endif
#ifdef TRACE_CE
#include <iostream>
#endif
namespace ROOT {
namespace Math {
template <class ScalarType>
class PtEtaPhiM4D {
public :
typedef ScalarType Scalar;
PtEtaPhiM4D() : fPt(0), fEta(0), fPhi(0), fM(0) { }
PtEtaPhiM4D(Scalar pt, Scalar eta, Scalar phi, Scalar mass) :
fPt(pt), fEta(eta), fPhi(phi), fM(mass) {
RestrictPhi();
if (fM < 0) RestrictNegMass();
}
template <class CoordSystem >
explicit PtEtaPhiM4D(const CoordSystem & c) :
fPt(c.Pt()), fEta(c.Eta()), fPhi(c.Phi()), fM(c.M()) { RestrictPhi(); }
PtEtaPhiM4D(const PtEtaPhiM4D & v) :
fPt(v.fPt), fEta(v.fEta), fPhi(v.fPhi), fM(v.fM) { }
PtEtaPhiM4D & operator = (const PtEtaPhiM4D & v) {
fPt = v.fPt;
fEta = v.fEta;
fPhi = v.fPhi;
fM = v.fM;
return *this;
}
void SetCoordinates( const Scalar src[] ) {
fPt=src[0]; fEta=src[1]; fPhi=src[2]; fM=src[3];
RestrictPhi();
if (fM <0) RestrictNegMass();
}
void GetCoordinates( Scalar dest[] ) const
{ dest[0] = fPt; dest[1] = fEta; dest[2] = fPhi; dest[3] = fM; }
void SetCoordinates(Scalar pt, Scalar eta, Scalar phi, Scalar mass) {
fPt=pt; fEta = eta; fPhi = phi; fM = mass;
RestrictPhi();
if (fM <0) RestrictNegMass();
}
void
GetCoordinates(Scalar& pt, Scalar & eta, Scalar & phi, Scalar& mass) const
{ pt=fPt; eta=fEta; phi = fPhi; mass = fM; }
Scalar Pt() const { return fPt; }
Scalar Eta() const { return fEta; }
Scalar Phi() const { return fPhi; }
Scalar M() const { return fM; }
Scalar Mag() const { return M(); }
Scalar Perp()const { return Pt(); }
Scalar Rho() const { return Pt(); }
Scalar Px() const { return fPt*cos(fPhi);}
Scalar X () const { return Px(); }
Scalar Py() const { return fPt*sin(fPhi);}
Scalar Y () const { return Py(); }
Scalar Pz() const {
return fPt > 0 ? fPt*std::sinh(fEta) :
fEta == 0 ? 0 :
fEta > 0 ? fEta - etaMax<Scalar>() :
fEta + etaMax<Scalar>() ;
}
Scalar Z () const { return Pz(); }
Scalar P() const {
return fPt > 0 ? fPt*std::cosh(fEta) :
fEta > etaMax<Scalar>() ? fEta - etaMax<Scalar>() :
fEta < -etaMax<Scalar>() ? -fEta - etaMax<Scalar>() :
0 ;
}
Scalar R() const { return P(); }
Scalar P2() const { Scalar p = P(); return p*p; }
Scalar E2() const {
Scalar e2 = P2() + M2();
return e2 > 0 ? e2 : 0;
}
Scalar E() const { return std::sqrt(E2() ); }
Scalar T() const { return E(); }
Scalar M2() const {
return ( fM >= 0 ) ? fM*fM : -fM*fM;
}
Scalar Mag2() const { return M2(); }
Scalar Pt2() const { return fPt*fPt;}
Scalar Perp2() const { return Pt2(); }
Scalar Mt2() const { return M2() + fPt*fPt; }
Scalar Mt() const {
Scalar mm = Mt2();
if (mm >= 0) {
return std::sqrt(mm);
} else {
GenVector::Throw ("PtEtaPhiM4D::Mt() - Tachyonic:\n"
" Pz^2 > E^2 so the transverse mass would be imaginary");
return -std::sqrt(-mm);
}
}
Scalar Et2() const {
return 2. * E2()/ ( std::cosh(2 * fEta) + 1 );
}
Scalar Et() const {
return E() / std::cosh(fEta);
}
private:
inline static Scalar pi() { return M_PI; }
inline void RestrictPhi() {
if ( fPhi <= -pi() || fPhi > pi() )
fPhi = fPhi - std::floor( fPhi/(2*pi()) +.5 ) * 2*pi();
return;
}
inline void RestrictNegMass() {
if ( fM >=0 ) return;
if ( P2() - fM*fM < 0 ) {
GenVector::Throw ("PtEtaPhiM4D::unphysical value of mass, set to closest physical value");
fM = - P();
}
return;
}
public:
Scalar Theta() const {
if (fPt > 0) return 2* std::atan( exp( - fEta ) );
if (fEta >= 0) return 0;
return pi();
}
void SetPt( Scalar pt) {
fPt = pt;
}
void SetEta( Scalar eta) {
fEta = eta;
}
void SetPhi( Scalar phi) {
fPhi = phi;
RestrictPhi();
}
void SetM( Scalar mass) {
fM = mass;
if (fM <0) RestrictNegMass();
}
void SetPxPyPzE(Scalar px, Scalar py, Scalar pz, Scalar e);
void Negate( ) {
fPhi = ( (fPhi > 0) ? fPhi - pi() : fPhi + pi() );
fEta = - fEta;
GenVector::Throw ("PtEtaPhiM4D::Negate - cannot negate the energy - can negate only the spatial components");
}
void Scale( Scalar a) {
if (a < 0) {
Negate(); a = -a;
}
fPt *= a;
fM *= a;
}
template <class CoordSystem >
PtEtaPhiM4D & operator = (const CoordSystem & c) {
fPt = c.Pt();
fEta = c.Eta();
fPhi = c.Phi();
fM = c.M();
return *this;
}
bool operator == (const PtEtaPhiM4D & rhs) const {
return fPt == rhs.fPt && fEta == rhs.fEta
&& fPhi == rhs.fPhi && fM == rhs.fM;
}
bool operator != (const PtEtaPhiM4D & rhs) const {return !(operator==(rhs));}
Scalar x() const { return X(); }
Scalar y() const { return Y(); }
Scalar z() const { return Z(); }
Scalar t() const { return E(); }
#if defined(__MAKECINT__) || defined(G__DICTIONARY)
void SetPx(Scalar px);
void SetPy(Scalar py);
void SetPz(Scalar pz);
void SetE(Scalar t);
#endif
private:
ScalarType fPt;
ScalarType fEta;
ScalarType fPhi;
ScalarType fM;
};
}
}
#ifndef ROOT_Math_GenVector_PxPyPzE4D
#include "Math/GenVector/PxPyPzE4D.h"
#endif
namespace ROOT {
namespace Math {
template <class ScalarType>
inline void PtEtaPhiM4D<ScalarType>::SetPxPyPzE(Scalar px, Scalar py, Scalar pz, Scalar e) {
*this = PxPyPzE4D<Scalar> (px, py, pz, e);
}
#if defined(__MAKECINT__) || defined(G__DICTIONARY)
template <class ScalarType>
void PtEtaPhiM4D<ScalarType>::SetPx(Scalar px) {
GenVector_exception e("PtEtaPhiM4D::SetPx() is not supposed to be called");
throw e;
PxPyPzE4D<Scalar> v(*this); v.SetPx(px); *this = PtEtaPhiM4D<Scalar>(v);
}
template <class ScalarType>
void PtEtaPhiM4D<ScalarType>::SetPy(Scalar py) {
GenVector_exception e("PtEtaPhiM4D::SetPx() is not supposed to be called");
throw e;
PxPyPzE4D<Scalar> v(*this); v.SetPy(py); *this = PtEtaPhiM4D<Scalar>(v);
}
template <class ScalarType>
void PtEtaPhiM4D<ScalarType>::SetPz(Scalar pz) {
GenVector_exception e("PtEtaPhiM4D::SetPx() is not supposed to be called");
throw e;
PxPyPzE4D<Scalar> v(*this); v.SetPz(pz); *this = PtEtaPhiM4D<Scalar>(v);
}
template <class ScalarType>
void PtEtaPhiM4D<ScalarType>::SetE(Scalar energy) {
GenVector_exception e("PtEtaPhiM4D::SetE() is not supposed to be called");
throw e;
PxPyPzE4D<Scalar> v(*this); v.SetE(energy); *this = PtEtaPhiM4D<Scalar>(v);
}
#endif // endif __MAKE__CINT || G__DICTIONARY
}
}
#endif // ROOT_Math_GenVector_PtEtaPhiM4D