// @(#)root/mathcore:$Id$
// Authors: W. Brown, M. Fischler, L. Moneta    2005

/**********************************************************************
 *                                                                    *
 * Copyright (c) 2005 , LCG ROOT MathLib Team                         *
 *                                                                    *
 *                                                                    *
 **********************************************************************/

// Header file for class LorentzVector
//
// Created by:    moneta   at Tue May 31 17:06:09 2005
// Major mods by: fischler at Wed Jul 20   2005
//
// Last update: $Id$
//
#ifndef ROOT_Math_GenVector_LorentzVector
#define ROOT_Math_GenVector_LorentzVector  1

#ifndef ROOT_Math_GenVector_PxPyPzE4D
#include "Math/GenVector/PxPyPzE4D.h"
#endif

#ifndef ROOT_Math_GenVector_DisplacementVector3D
#include "Math/GenVector/DisplacementVector3D.h"
#endif

#ifndef ROOT_Math_GenVector_GenVectorIO
#include "Math/GenVector/GenVectorIO.h"
#endif



namespace ROOT {

  namespace Math {

//__________________________________________________________________________________________
    /**
        Class describing a generic LorentzVector in the 4D space-time,
        using the specified coordinate system for the spatial vector part.
        The metric used for the LorentzVector is (-,-,-,+).
        In the case of LorentzVector we don't distinguish the concepts
        of points and displacement vectors as in the 3D case,
        since the main use case for 4D Vectors is to describe the kinematics of
        relativistic particles. A LorentzVector behaves like a
        DisplacementVector in 4D.  The Minkowski components could be viewed as
        v and t, or for kinematic 4-vectors, as p and E.

     @ingroup GenVector
    */
    template< class CoordSystem >
    class LorentzVector {

    public:

       // ------ ctors ------

       typedef typename CoordSystem::Scalar Scalar;
       typedef CoordSystem CoordinateType;

       /**
          default constructor of an empty vector (Px = Py = Pz = E = 0 )
       */
       LorentzVector ( ) : fCoordinates() { }

       /**
          generic constructors from four scalar values.
          The association between values and coordinate depends on the
          coordinate system.  For PxPyPzE4D,
          \param a scalar value (Px)
          \param b scalar value (Py)
          \param c scalar value (Pz)
          \param d scalar value (E)
       */
       LorentzVector(const Scalar & a,
                     const Scalar & b,
                     const Scalar & c,
                     const Scalar & d) :
          fCoordinates(a , b,  c, d)  { }

       /**
          constructor from a LorentzVector expressed in different
          coordinates, or using a different Scalar type
       */
       template< class Coords >
       explicit LorentzVector(const LorentzVector<Coords> & v ) :
          fCoordinates( v.Coordinates() ) { }

       /**
          Construct from a foreign 4D vector type, for example, HepLorentzVector
          Precondition: v must implement methods x(), y(), z(), and t()
       */
       template<class ForeignLorentzVector>
       explicit LorentzVector( const ForeignLorentzVector & v) :
          fCoordinates(PxPyPzE4D<Scalar>( v.x(), v.y(), v.z(), v.t()  ) ) { }

#ifdef LATER
       /**
          construct from a generic linear algebra  vector implementing operator []
          and with a size of at least 4. This could be also a C array
          In this case v[0] is the first data member
          ( Px for a PxPyPzE4D base)
          \param v LA vector
          \param index0 index of first vector element (Px)
       */
       template< class LAVector >
       explicit LorentzVector(const LAVector & v, size_t index0 ) {
          fCoordinates = CoordSystem ( v[index0], v[index0+1], v[index0+2], v[index0+3] );
       }
#endif


       // ------ assignment ------

       /**
          Assignment operator from a lorentz vector of arbitrary type
       */
       template< class OtherCoords >
       LorentzVector & operator= ( const LorentzVector<OtherCoords> & v) {
          fCoordinates = v.Coordinates();
          return *this;
       }

       /**
          assignment from any other Lorentz vector  implementing
          x(), y(), z() and t()
       */
       template<class ForeignLorentzVector>
       LorentzVector & operator = ( const ForeignLorentzVector & v) {
          SetXYZT( v.x(), v.y(), v.z(), v.t() );
          return *this;
       }

#ifdef LATER
       /**
          assign from a generic linear algebra  vector implementing operator []
          and with a size of at least 4
          In this case v[0] is the first data member
          ( Px for a PxPyPzE4D base)
          \param v LA vector
          \param index0 index of first vector element (Px)
       */
       template< class LAVector >
       LorentzVector & AssignFrom(const LAVector & v, size_t index0=0 ) {
          fCoordinates.SetCoordinates( v[index0], v[index0+1], v[index0+2], v[index0+3] );
          return *this;
       }
#endif

       // ------ Set, Get, and access coordinate data ------

       /**
          Retrieve a const reference to  the coordinates object
       */
       const CoordSystem & Coordinates() const {
          return fCoordinates;
       }

       /**
          Set internal data based on an array of 4 Scalar numbers
       */
       LorentzVector<CoordSystem>& SetCoordinates( const Scalar src[] ) {
          fCoordinates.SetCoordinates(src);
          return *this;
       }

       /**
          Set internal data based on 4 Scalar numbers
       */
       LorentzVector<CoordSystem>& SetCoordinates( Scalar a, Scalar b, Scalar c, Scalar d ) {
          fCoordinates.SetCoordinates(a, b, c, d);
          return *this;
       }

       /**
          Set internal data based on 4 Scalars at *begin to *end
       */
//#ifdef NDEBUG
          //this does not compile in CINT
//        template< class IT >
//        LorentzVector<CoordSystem>& SetCoordinates( IT begin, IT /* end */  ) {
// #endif
       template< class IT >
#ifndef NDEBUG
       LorentzVector<CoordSystem>& SetCoordinates( IT begin, IT end  ) {
#else
       LorentzVector<CoordSystem>& SetCoordinates( IT begin, IT /* end */  ) {
#endif
          IT a = begin; IT b = ++begin; IT c = ++begin; IT d = ++begin;
          assert (++begin==end);
          SetCoordinates (*a,*b,*c,*d);
          return *this;
       }

       /**
          get internal data into 4 Scalar numbers
       */
       void GetCoordinates( Scalar& a, Scalar& b, Scalar& c, Scalar & d ) const
       { fCoordinates.GetCoordinates(a, b, c, d);  }

       /**
          get internal data into an array of 4 Scalar numbers
       */
       void GetCoordinates( Scalar dest[] ) const
       { fCoordinates.GetCoordinates(dest);  }

       /**
          get internal data into 4 Scalars at *begin to *end
       */
       template <class IT>
#ifndef NDEBUG
       void GetCoordinates( IT begin, IT end ) const
#else
       void GetCoordinates( IT begin, IT /* end */ ) const
#endif
       { IT a = begin; IT b = ++begin; IT c = ++begin; IT d = ++begin;
       assert (++begin==end);
       GetCoordinates (*a,*b,*c,*d);
       }

       /**
          get internal data into 4 Scalars at *begin
       */
       template <class IT>
       void GetCoordinates( IT begin ) const {
          Scalar a,b,c,d = 0;
          GetCoordinates (a,b,c,d);
          *begin++ = a;
          *begin++ = b;
          *begin++ = c;
          *begin   = d;
       }

       /**
          set the values of the vector from the cartesian components (x,y,z,t)
          (if the vector is held in another coordinates, like (Pt,eta,phi,m)
          then (x, y, z, t) are converted to that form)
       */
       LorentzVector<CoordSystem>& SetXYZT (Scalar xx, Scalar yy, Scalar zz, Scalar tt) {
          fCoordinates.SetPxPyPzE(xx,yy,zz,tt);
          return *this;
       }
       LorentzVector<CoordSystem>& SetPxPyPzE (Scalar xx, Scalar yy, Scalar zz, Scalar ee) {
          fCoordinates.SetPxPyPzE(xx,yy,zz,ee);
          return *this;
       }

       // ------------------- Equality -----------------

       /**
          Exact equality
       */
       bool operator==(const LorentzVector & rhs) const {
          return fCoordinates==rhs.fCoordinates;
       }
       bool operator!= (const LorentzVector & rhs) const {
          return !(operator==(rhs));
       }

       // ------ Individual element access, in various coordinate systems ------

       // individual coordinate accessors in various coordinate systems

       /**
          spatial X component
       */
       Scalar Px() const  { return fCoordinates.Px(); }
       Scalar X()  const  { return fCoordinates.Px(); }
       /**
          spatial Y component
       */
       Scalar Py() const { return fCoordinates.Py(); }
       Scalar Y()  const { return fCoordinates.Py(); }
       /**
          spatial Z component
       */
       Scalar Pz() const { return fCoordinates.Pz(); }
       Scalar Z()  const { return fCoordinates.Pz(); }
       /**
          return 4-th component (time, or energy for a 4-momentum vector)
       */
       Scalar E()  const { return fCoordinates.E(); }
       Scalar T()  const { return fCoordinates.E(); }
       /**
          return magnitude (mass) squared  M2 = T**2 - X**2 - Y**2 - Z**2
          (we use -,-,-,+ metric)
       */
       Scalar M2()   const { return fCoordinates.M2(); }
       /**
          return magnitude (mass) using the  (-,-,-,+)  metric.
          If M2 is negative (space-like vector) a GenVector_exception
          is suggested and if continuing, - sqrt( -M2) is returned
       */
       Scalar M() const    { return fCoordinates.M();}
       /**
          return the spatial (3D) magnitude ( sqrt(X**2 + Y**2 + Z**2) )
       */
       Scalar R() const { return fCoordinates.R(); }
       Scalar P() const { return fCoordinates.R(); }
       /**
          return the square of the spatial (3D) magnitude ( X**2 + Y**2 + Z**2 )
       */
       Scalar P2() const { return P() * P(); }
       /**
          return the square of the transverse spatial component ( X**2 + Y**2 )
       */
       Scalar Perp2( ) const { return fCoordinates.Perp2();}

       /**
          return the  transverse spatial component sqrt ( X**2 + Y**2 )
       */
       Scalar Pt()  const { return fCoordinates.Pt(); }
       Scalar Rho() const { return fCoordinates.Pt(); }

       /**
          return the transverse mass squared
          \f[ m_t^2 = E^2 - p{_z}^2 \f]
       */
       Scalar Mt2() const { return fCoordinates.Mt2(); }

       /**
          return the transverse mass
          \f[ \sqrt{ m_t^2 = E^2 - p{_z}^2} X sign(E^ - p{_z}^2) \f]
       */
       Scalar Mt() const { return fCoordinates.Mt(); }

       /**
          return the transverse energy squared
          \f[ e_t = \frac{E^2 p_{\perp}^2 }{ |p|^2 } \f]
       */
       Scalar Et2() const { return fCoordinates.Et2(); }

       /**
          return the transverse energy
          \f[ e_t = \sqrt{ \frac{E^2 p_{\perp}^2 }{ |p|^2 } } X sign(E) \f]
       */
       Scalar Et() const { return fCoordinates.Et(); }

       /**
          azimuthal  Angle
       */
       Scalar Phi() const  { return fCoordinates.Phi();}

       /**
          polar Angle
       */
       Scalar Theta() const { return fCoordinates.Theta(); }

       /**
          pseudorapidity
          \f[ \eta = - \ln { \tan { \frac { \theta} {2} } } \f]
       */
       Scalar Eta() const { return fCoordinates.Eta(); }

       /**
          get the spatial components of the Vector in a
          DisplacementVector based on Cartesian Coordinates
       */
       ::ROOT::Math::DisplacementVector3D<Cartesian3D<Scalar> > Vect() const {
          return ::ROOT::Math::DisplacementVector3D<Cartesian3D<Scalar> >( X(), Y(), Z() );
       }

       // ------ Operations combining two Lorentz vectors ------

       /**
          scalar (Dot) product of two LorentzVector vectors (metric is -,-,-,+)
          Enable the product using any other LorentzVector implementing
          the x(), y() , y() and t() member functions
          \param  q  any LorentzVector implementing the x(), y() , z() and t()
          member functions
          \return the result of v.q of type according to the base scalar type of v
       */

       template< class OtherLorentzVector >
       Scalar Dot(const OtherLorentzVector & q) const {
          return t()*q.t() - x()*q.x() - y()*q.y() - z()*q.z();
       }

       /**
          Self addition with another Vector ( v+= q )
          Enable the addition with any other LorentzVector
          \param  q  any LorentzVector implementing the x(), y() , z() and t()
          member functions
       */
      template< class OtherLorentzVector >
      inline LorentzVector & operator += ( const OtherLorentzVector & q)
       {
          SetXYZT( x() + q.x(), y() + q.y(), z() + q.z(), t() + q.t()  );
          return *this;
       }

       /**
          Self subtraction of another Vector from this ( v-= q )
          Enable the addition with any other LorentzVector
          \param  q  any LorentzVector implementing the x(), y() , z() and t()
          member functions
       */
       template< class OtherLorentzVector >
       LorentzVector & operator -= ( const OtherLorentzVector & q) {
          SetXYZT( x() - q.x(), y() - q.y(), z() - q.z(), t() - q.t()  );
          return *this;
       }

       /**
          addition of two LorentzVectors (v3 = v1 + v2)
          Enable the addition with any other LorentzVector
          \param  v2  any LorentzVector implementing the x(), y() , z() and t()
          member functions
          \return a new LorentzVector of the same type as v1
       */
       template<class OtherLorentzVector>
       LorentzVector  operator +  ( const OtherLorentzVector & v2) const
       {
          LorentzVector<CoordinateType> v3(*this);
          v3 += v2;
          return v3;
       }

       /**
          subtraction of two LorentzVectors (v3 = v1 - v2)
          Enable the subtraction of any other LorentzVector
          \param  v2  any LorentzVector implementing the x(), y() , z() and t()
          member functions
          \return a new LorentzVector of the same type as v1
       */
       template<class OtherLorentzVector>
       LorentzVector  operator -  ( const OtherLorentzVector & v2) const {
          LorentzVector<CoordinateType> v3(*this);
          v3 -= v2;
          return v3;
       }

       //--- scaling operations ------

       /**
          multiplication by a scalar quantity v *= a
       */
       LorentzVector & operator *= (Scalar a) {
          fCoordinates.Scale(a);
          return *this;
       }

       /**
          division by a scalar quantity v /= a
       */
       LorentzVector & operator /= (Scalar a) {
          fCoordinates.Scale(1/a);
          return *this;
       }

       /**
          product of a LorentzVector by a scalar quantity
          \param a  scalar quantity of type a
          \return a new mathcoreLorentzVector q = v * a same type as v
       */
       LorentzVector operator * ( const Scalar & a) const {
          LorentzVector tmp(*this);
          tmp *= a;
          return tmp;
       }

       /**
          Divide a LorentzVector by a scalar quantity
          \param a  scalar quantity of type a
          \return a new mathcoreLorentzVector q = v / a same type as v
       */
       LorentzVector<CoordSystem> operator / ( const Scalar & a) const {
          LorentzVector<CoordSystem> tmp(*this);
          tmp /= a;
          return tmp;
       }

       /**
          Negative of a LorentzVector (q = - v )
          \return a new LorentzVector with opposite direction and time
       */
       LorentzVector operator - () const {
          //LorentzVector<CoordinateType> v(*this);
          //v.Negate();
          return operator*( Scalar(-1) );
       }
       LorentzVector operator + () const {
          return *this;
       }

       // ---- Relativistic Properties ----

       /**
          Rapidity relative to the Z axis:  .5 log [(E+Pz)/(E-Pz)]
       */
       Scalar Rapidity() const {
          // TODO - It would be good to check that E > Pz and use the Throw()
          //        mechanism or at least load a NAN if not.
          //        We should then move the code to a .cpp file.
          Scalar ee = E();
          Scalar ppz = Pz();
          return .5f* std::log( (ee+ppz)/(ee-ppz) );
       }

       /**
          Rapidity in the direction of travel: atanh (|P|/E)=.5 log[(E+P)/(E-P)]
       */
       Scalar ColinearRapidity() const {
          // TODO - It would be good to check that E > P and use the Throw()
          //        mechanism or at least load a NAN if not.
          Scalar ee = E();
          Scalar pp = P();
          return .5f* std::log( (ee+pp)/(ee-pp) );
       }

       /**
          Determine if momentum-energy can represent a physical massive particle
       */
       bool isTimelike( ) const {
          Scalar ee = E(); Scalar pp = P(); return ee*ee > pp*pp;
       }

       /**
          Determine if momentum-energy can represent a massless particle
       */
       bool isLightlike( Scalar tolerance
                         = 100*std::numeric_limits<Scalar>::epsilon() ) const {
          Scalar ee = E(); Scalar pp = P(); Scalar delta = ee-pp;
          if ( ee==0 ) return pp==0;
          return delta*delta < tolerance * ee*ee;
       }

       /**
          Determine if momentum-energy is spacelike, and represents a tachyon
       */
       bool isSpacelike( ) const {
          Scalar ee = E(); Scalar pp = P(); return ee*ee < pp*pp;
       }

       typedef DisplacementVector3D< Cartesian3D<Scalar> > BetaVector;

       /**
          The beta vector for the boost that would bring this vector into
          its center of mass frame (zero momentum)
       */
       BetaVector BoostToCM( ) const {
          if (E() == 0) {
             if (P() == 0) {
                return BetaVector();
             } else {
                // TODO - should attempt to Throw with msg about
                // boostVector computed for LorentzVector with t=0
                return -Vect()/E();
             }
          }
          if (M2() <= 0) {
             // TODO - should attempt to Throw with msg about
             // boostVector computed for a non-timelike LorentzVector
          }
          return -Vect()/E();
       }

       /**
          The beta vector for the boost that would bring this vector into
          its center of mass frame (zero momentum)
       */
       template <class Other4Vector>
       BetaVector BoostToCM(const Other4Vector& v ) const {
          Scalar eSum = E() + v.E();
          DisplacementVector3D< Cartesian3D<Scalar> > vecSum = Vect() + v.Vect();
          if (eSum == 0) {
             if (vecSum.Mag2() == 0) {
                return BetaVector();
             } else {
                // TODO - should attempt to Throw with msg about
                // boostToCM computed for two 4-vectors with combined t=0
                return BetaVector(vecSum/eSum);
             }
             // TODO - should attempt to Throw with msg about
             // boostToCM computed for two 4-vectors with combined e=0
          }
          return BetaVector (vecSum * (-1./eSum));
       }

       //beta and gamma

       /**
           Return beta scalar value
       */
       Scalar Beta() const {
          if ( E() == 0 ) {
             if ( P2() == 0)
                // to avoid Nan
                return 0;
             else {
                GenVector::Throw ("LorentzVector::Beta() - beta computed for LorentzVector with t = 0. Return an Infinite result");
                return 1./E();
             }
          }
          if ( M2() <= 0 ) {
             GenVector::Throw ("LorentzVector::Beta() - beta computed for non-timelike LorentzVector . Result is physically meaningless" );
          }
          return P() / E();
       }
       /**
           Return Gamma scalar value
       */
       Scalar Gamma() const {
          Scalar v2 = P2();
          Scalar t2 = E()*E();
          if (E() == 0) {
             if ( P2() == 0) {
                return 1;
             } else {
                GenVector::Throw ("LorentzVector::Gamma() - gamma computed for LorentzVector with t = 0. Return a zero result");

             }
          }
          if ( t2 < v2 ) {
             GenVector::Throw ("LorentzVector::Gamma() - gamma computed for a spacelike LorentzVector. Imaginary result");
             return 0;
          }
          else if ( t2 == v2 ) {
             GenVector::Throw ("LorentzVector::Gamma() - gamma computed for a lightlike LorentzVector. Infinite result");
          }
          return 1./std::sqrt(1. - v2/t2 );
       } /* gamma */


       // Method providing limited backward name compatibility with CLHEP ----

       Scalar x()     const { return fCoordinates.Px();     }
       Scalar y()     const { return fCoordinates.Py();     }
       Scalar z()     const { return fCoordinates.Pz();     }
       Scalar t()     const { return fCoordinates.E();      }
       Scalar px()    const { return fCoordinates.Px();     }
       Scalar py()    const { return fCoordinates.Py();     }
       Scalar pz()    const { return fCoordinates.Pz();     }
       Scalar e()     const { return fCoordinates.E();      }
       Scalar r()     const { return fCoordinates.R();      }
       Scalar theta() const { return fCoordinates.Theta();  }
       Scalar phi()   const { return fCoordinates.Phi();    }
       Scalar rho()   const { return fCoordinates.Rho();    }
       Scalar eta()   const { return fCoordinates.Eta();    }
       Scalar pt()    const { return fCoordinates.Pt();     }
       Scalar perp2() const { return fCoordinates.Perp2();  }
       Scalar mag2()  const { return fCoordinates.M2();     }
       Scalar mag()   const { return fCoordinates.M();      }
       Scalar mt()    const { return fCoordinates.Mt();     }
       Scalar mt2()   const { return fCoordinates.Mt2();    }


       // Methods  requested by CMS ---
       Scalar energy() const { return fCoordinates.E();      }
       Scalar mass()   const { return fCoordinates.M();      }
       Scalar mass2()  const { return fCoordinates.M2();     }


       /**
          Methods setting a Single-component
          Work only if the component is one of which the vector is represented.
          For example SetE will work for a PxPyPzE Vector but not for a PxPyPzM Vector.
       */
       LorentzVector<CoordSystem>& SetE  ( Scalar a )  { fCoordinates.SetE  (a); return *this; }
       LorentzVector<CoordSystem>& SetEta( Scalar a )  { fCoordinates.SetEta(a); return *this; }
       LorentzVector<CoordSystem>& SetM  ( Scalar a )  { fCoordinates.SetM  (a); return *this; }
       LorentzVector<CoordSystem>& SetPhi( Scalar a )  { fCoordinates.SetPhi(a); return *this; }
       LorentzVector<CoordSystem>& SetPt ( Scalar a )  { fCoordinates.SetPt (a); return *this; }
       LorentzVector<CoordSystem>& SetPx ( Scalar a )  { fCoordinates.SetPx (a); return *this; }
       LorentzVector<CoordSystem>& SetPy ( Scalar a )  { fCoordinates.SetPy (a); return *this; }
       LorentzVector<CoordSystem>& SetPz ( Scalar a )  { fCoordinates.SetPz (a); return *this; }

    private:

       CoordSystem  fCoordinates;    // internal coordinate system


    };  // LorentzVector<>



  // global nethods

  /**
     Scale of a LorentzVector with a scalar quantity a
     \param a  scalar quantity of typpe a
     \param v  mathcore::LorentzVector based on any coordinate system
     \return a new mathcoreLorentzVector q = v * a same type as v
   */
    template< class CoordSystem >
    inline LorentzVector<CoordSystem> operator *
    ( const typename  LorentzVector<CoordSystem>::Scalar & a,
      const LorentzVector<CoordSystem>& v) {
       LorentzVector<CoordSystem> tmp(v);
       tmp *= a;
       return tmp;
    }

    // ------------- I/O to/from streams -------------

    template< class char_t, class traits_t, class Coords >
    inline
    std::basic_ostream<char_t,traits_t> &
    operator << ( std::basic_ostream<char_t,traits_t> & os
                  , LorentzVector<Coords> const & v
       )
    {
       if( !os )  return os;

       typename Coords::Scalar a, b, c, d;
       v.GetCoordinates(a, b, c, d);

       if( detail::get_manip( os, detail::bitforbit ) )  {
        detail::set_manip( os, detail::bitforbit, '\00' );
        // TODO: call MF's bitwise-accurate functions on each of a, b, c, d
       }
       else  {
          os << detail::get_manip( os, detail::open  ) << a
             << detail::get_manip( os, detail::sep   ) << b
             << detail::get_manip( os, detail::sep   ) << c
             << detail::get_manip( os, detail::sep   ) << d
             << detail::get_manip( os, detail::close );
       }

       return os;

    }  // op<< <>()


     template< class char_t, class traits_t, class Coords >
     inline
     std::basic_istream<char_t,traits_t> &
     operator >> ( std::basic_istream<char_t,traits_t> & is
                   , LorentzVector<Coords> & v
        )
     {
        if( !is )  return is;

        typename Coords::Scalar a, b, c, d;

        if( detail::get_manip( is, detail::bitforbit ) )  {
           detail::set_manip( is, detail::bitforbit, '\00' );
           // TODO: call MF's bitwise-accurate functions on each of a, b, c
        }
        else  {
           detail::require_delim( is, detail::open  );  is >> a;
           detail::require_delim( is, detail::sep   );  is >> b;
           detail::require_delim( is, detail::sep   );  is >> c;
           detail::require_delim( is, detail::sep   );  is >> d;
           detail::require_delim( is, detail::close );
        }

        if( is )
           v.SetCoordinates(a, b, c, d);
        return is;

     }  // op>> <>()



  } // end namespace Math

} // end namespace ROOT



#endif

//#include "Math/GenVector/LorentzVectorOperations.h"



 LorentzVector.h:1
 LorentzVector.h:2
 LorentzVector.h:3
 LorentzVector.h:4
 LorentzVector.h:5
 LorentzVector.h:6
 LorentzVector.h:7
 LorentzVector.h:8
 LorentzVector.h:9
 LorentzVector.h:10
 LorentzVector.h:11
 LorentzVector.h:12
 LorentzVector.h:13
 LorentzVector.h:14
 LorentzVector.h:15
 LorentzVector.h:16
 LorentzVector.h:17
 LorentzVector.h:18
 LorentzVector.h:19
 LorentzVector.h:20
 LorentzVector.h:21
 LorentzVector.h:22
 LorentzVector.h:23
 LorentzVector.h:24
 LorentzVector.h:25
 LorentzVector.h:26
 LorentzVector.h:27
 LorentzVector.h:28
 LorentzVector.h:29
 LorentzVector.h:30
 LorentzVector.h:31
 LorentzVector.h:32
 LorentzVector.h:33
 LorentzVector.h:34
 LorentzVector.h:35
 LorentzVector.h:36
 LorentzVector.h:37
 LorentzVector.h:38
 LorentzVector.h:39
 LorentzVector.h:40
 LorentzVector.h:41
 LorentzVector.h:42
 LorentzVector.h:43
 LorentzVector.h:44
 LorentzVector.h:45
 LorentzVector.h:46
 LorentzVector.h:47
 LorentzVector.h:48
 LorentzVector.h:49
 LorentzVector.h:50
 LorentzVector.h:51
 LorentzVector.h:52
 LorentzVector.h:53
 LorentzVector.h:54
 LorentzVector.h:55
 LorentzVector.h:56
 LorentzVector.h:57
 LorentzVector.h:58
 LorentzVector.h:59
 LorentzVector.h:60
 LorentzVector.h:61
 LorentzVector.h:62
 LorentzVector.h:63
 LorentzVector.h:64
 LorentzVector.h:65
 LorentzVector.h:66
 LorentzVector.h:67
 LorentzVector.h:68
 LorentzVector.h:69
 LorentzVector.h:70
 LorentzVector.h:71
 LorentzVector.h:72
 LorentzVector.h:73
 LorentzVector.h:74
 LorentzVector.h:75
 LorentzVector.h:76
 LorentzVector.h:77
 LorentzVector.h:78
 LorentzVector.h:79
 LorentzVector.h:80
 LorentzVector.h:81
 LorentzVector.h:82
 LorentzVector.h:83
 LorentzVector.h:84
 LorentzVector.h:85
 LorentzVector.h:86
 LorentzVector.h:87
 LorentzVector.h:88
 LorentzVector.h:89
 LorentzVector.h:90
 LorentzVector.h:91
 LorentzVector.h:92
 LorentzVector.h:93
 LorentzVector.h:94
 LorentzVector.h:95
 LorentzVector.h:96
 LorentzVector.h:97
 LorentzVector.h:98
 LorentzVector.h:99
 LorentzVector.h:100
 LorentzVector.h:101
 LorentzVector.h:102
 LorentzVector.h:103
 LorentzVector.h:104
 LorentzVector.h:105
 LorentzVector.h:106
 LorentzVector.h:107
 LorentzVector.h:108
 LorentzVector.h:109
 LorentzVector.h:110
 LorentzVector.h:111
 LorentzVector.h:112
 LorentzVector.h:113
 LorentzVector.h:114
 LorentzVector.h:115
 LorentzVector.h:116
 LorentzVector.h:117
 LorentzVector.h:118
 LorentzVector.h:119
 LorentzVector.h:120
 LorentzVector.h:121
 LorentzVector.h:122
 LorentzVector.h:123
 LorentzVector.h:124
 LorentzVector.h:125
 LorentzVector.h:126
 LorentzVector.h:127
 LorentzVector.h:128
 LorentzVector.h:129
 LorentzVector.h:130
 LorentzVector.h:131
 LorentzVector.h:132
 LorentzVector.h:133
 LorentzVector.h:134
 LorentzVector.h:135
 LorentzVector.h:136
 LorentzVector.h:137
 LorentzVector.h:138
 LorentzVector.h:139
 LorentzVector.h:140
 LorentzVector.h:141
 LorentzVector.h:142
 LorentzVector.h:143
 LorentzVector.h:144
 LorentzVector.h:145
 LorentzVector.h:146
 LorentzVector.h:147
 LorentzVector.h:148
 LorentzVector.h:149
 LorentzVector.h:150
 LorentzVector.h:151
 LorentzVector.h:152
 LorentzVector.h:153
 LorentzVector.h:154
 LorentzVector.h:155
 LorentzVector.h:156
 LorentzVector.h:157
 LorentzVector.h:158
 LorentzVector.h:159
 LorentzVector.h:160
 LorentzVector.h:161
 LorentzVector.h:162
 LorentzVector.h:163
 LorentzVector.h:164
 LorentzVector.h:165
 LorentzVector.h:166
 LorentzVector.h:167
 LorentzVector.h:168
 LorentzVector.h:169
 LorentzVector.h:170
 LorentzVector.h:171
 LorentzVector.h:172
 LorentzVector.h:173
 LorentzVector.h:174
 LorentzVector.h:175
 LorentzVector.h:176
 LorentzVector.h:177
 LorentzVector.h:178
 LorentzVector.h:179
 LorentzVector.h:180
 LorentzVector.h:181
 LorentzVector.h:182
 LorentzVector.h:183
 LorentzVector.h:184
 LorentzVector.h:185
 LorentzVector.h:186
 LorentzVector.h:187
 LorentzVector.h:188
 LorentzVector.h:189
 LorentzVector.h:190
 LorentzVector.h:191
 LorentzVector.h:192
 LorentzVector.h:193
 LorentzVector.h:194
 LorentzVector.h:195
 LorentzVector.h:196
 LorentzVector.h:197
 LorentzVector.h:198
 LorentzVector.h:199
 LorentzVector.h:200
 LorentzVector.h:201
 LorentzVector.h:202
 LorentzVector.h:203
 LorentzVector.h:204
 LorentzVector.h:205
 LorentzVector.h:206
 LorentzVector.h:207
 LorentzVector.h:208
 LorentzVector.h:209
 LorentzVector.h:210
 LorentzVector.h:211
 LorentzVector.h:212
 LorentzVector.h:213
 LorentzVector.h:214
 LorentzVector.h:215
 LorentzVector.h:216
 LorentzVector.h:217
 LorentzVector.h:218
 LorentzVector.h:219
 LorentzVector.h:220
 LorentzVector.h:221
 LorentzVector.h:222
 LorentzVector.h:223
 LorentzVector.h:224
 LorentzVector.h:225
 LorentzVector.h:226
 LorentzVector.h:227
 LorentzVector.h:228
 LorentzVector.h:229
 LorentzVector.h:230
 LorentzVector.h:231
 LorentzVector.h:232
 LorentzVector.h:233
 LorentzVector.h:234
 LorentzVector.h:235
 LorentzVector.h:236
 LorentzVector.h:237
 LorentzVector.h:238
 LorentzVector.h:239
 LorentzVector.h:240
 LorentzVector.h:241
 LorentzVector.h:242
 LorentzVector.h:243
 LorentzVector.h:244
 LorentzVector.h:245
 LorentzVector.h:246
 LorentzVector.h:247
 LorentzVector.h:248
 LorentzVector.h:249
 LorentzVector.h:250
 LorentzVector.h:251
 LorentzVector.h:252
 LorentzVector.h:253
 LorentzVector.h:254
 LorentzVector.h:255
 LorentzVector.h:256
 LorentzVector.h:257
 LorentzVector.h:258
 LorentzVector.h:259
 LorentzVector.h:260
 LorentzVector.h:261
 LorentzVector.h:262
 LorentzVector.h:263
 LorentzVector.h:264
 LorentzVector.h:265
 LorentzVector.h:266
 LorentzVector.h:267
 LorentzVector.h:268
 LorentzVector.h:269
 LorentzVector.h:270
 LorentzVector.h:271
 LorentzVector.h:272
 LorentzVector.h:273
 LorentzVector.h:274
 LorentzVector.h:275
 LorentzVector.h:276
 LorentzVector.h:277
 LorentzVector.h:278
 LorentzVector.h:279
 LorentzVector.h:280
 LorentzVector.h:281
 LorentzVector.h:282
 LorentzVector.h:283
 LorentzVector.h:284
 LorentzVector.h:285
 LorentzVector.h:286
 LorentzVector.h:287
 LorentzVector.h:288
 LorentzVector.h:289
 LorentzVector.h:290
 LorentzVector.h:291
 LorentzVector.h:292
 LorentzVector.h:293
 LorentzVector.h:294
 LorentzVector.h:295
 LorentzVector.h:296
 LorentzVector.h:297
 LorentzVector.h:298
 LorentzVector.h:299
 LorentzVector.h:300
 LorentzVector.h:301
 LorentzVector.h:302
 LorentzVector.h:303
 LorentzVector.h:304
 LorentzVector.h:305
 LorentzVector.h:306
 LorentzVector.h:307
 LorentzVector.h:308
 LorentzVector.h:309
 LorentzVector.h:310
 LorentzVector.h:311
 LorentzVector.h:312
 LorentzVector.h:313
 LorentzVector.h:314
 LorentzVector.h:315
 LorentzVector.h:316
 LorentzVector.h:317
 LorentzVector.h:318
 LorentzVector.h:319
 LorentzVector.h:320
 LorentzVector.h:321
 LorentzVector.h:322
 LorentzVector.h:323
 LorentzVector.h:324
 LorentzVector.h:325
 LorentzVector.h:326
 LorentzVector.h:327
 LorentzVector.h:328
 LorentzVector.h:329
 LorentzVector.h:330
 LorentzVector.h:331
 LorentzVector.h:332
 LorentzVector.h:333
 LorentzVector.h:334
 LorentzVector.h:335
 LorentzVector.h:336
 LorentzVector.h:337
 LorentzVector.h:338
 LorentzVector.h:339
 LorentzVector.h:340
 LorentzVector.h:341
 LorentzVector.h:342
 LorentzVector.h:343
 LorentzVector.h:344
 LorentzVector.h:345
 LorentzVector.h:346
 LorentzVector.h:347
 LorentzVector.h:348
 LorentzVector.h:349
 LorentzVector.h:350
 LorentzVector.h:351
 LorentzVector.h:352
 LorentzVector.h:353
 LorentzVector.h:354
 LorentzVector.h:355
 LorentzVector.h:356
 LorentzVector.h:357
 LorentzVector.h:358
 LorentzVector.h:359
 LorentzVector.h:360
 LorentzVector.h:361
 LorentzVector.h:362
 LorentzVector.h:363
 LorentzVector.h:364
 LorentzVector.h:365
 LorentzVector.h:366
 LorentzVector.h:367
 LorentzVector.h:368
 LorentzVector.h:369
 LorentzVector.h:370
 LorentzVector.h:371
 LorentzVector.h:372
 LorentzVector.h:373
 LorentzVector.h:374
 LorentzVector.h:375
 LorentzVector.h:376
 LorentzVector.h:377
 LorentzVector.h:378
 LorentzVector.h:379
 LorentzVector.h:380
 LorentzVector.h:381
 LorentzVector.h:382
 LorentzVector.h:383
 LorentzVector.h:384
 LorentzVector.h:385
 LorentzVector.h:386
 LorentzVector.h:387
 LorentzVector.h:388
 LorentzVector.h:389
 LorentzVector.h:390
 LorentzVector.h:391
 LorentzVector.h:392
 LorentzVector.h:393
 LorentzVector.h:394
 LorentzVector.h:395
 LorentzVector.h:396
 LorentzVector.h:397
 LorentzVector.h:398
 LorentzVector.h:399
 LorentzVector.h:400
 LorentzVector.h:401
 LorentzVector.h:402
 LorentzVector.h:403
 LorentzVector.h:404
 LorentzVector.h:405
 LorentzVector.h:406
 LorentzVector.h:407
 LorentzVector.h:408
 LorentzVector.h:409
 LorentzVector.h:410
 LorentzVector.h:411
 LorentzVector.h:412
 LorentzVector.h:413
 LorentzVector.h:414
 LorentzVector.h:415
 LorentzVector.h:416
 LorentzVector.h:417
 LorentzVector.h:418
 LorentzVector.h:419
 LorentzVector.h:420
 LorentzVector.h:421
 LorentzVector.h:422
 LorentzVector.h:423
 LorentzVector.h:424
 LorentzVector.h:425
 LorentzVector.h:426
 LorentzVector.h:427
 LorentzVector.h:428
 LorentzVector.h:429
 LorentzVector.h:430
 LorentzVector.h:431
 LorentzVector.h:432
 LorentzVector.h:433
 LorentzVector.h:434
 LorentzVector.h:435
 LorentzVector.h:436
 LorentzVector.h:437
 LorentzVector.h:438
 LorentzVector.h:439
 LorentzVector.h:440
 LorentzVector.h:441
 LorentzVector.h:442
 LorentzVector.h:443
 LorentzVector.h:444
 LorentzVector.h:445
 LorentzVector.h:446
 LorentzVector.h:447
 LorentzVector.h:448
 LorentzVector.h:449
 LorentzVector.h:450
 LorentzVector.h:451
 LorentzVector.h:452
 LorentzVector.h:453
 LorentzVector.h:454
 LorentzVector.h:455
 LorentzVector.h:456
 LorentzVector.h:457
 LorentzVector.h:458
 LorentzVector.h:459
 LorentzVector.h:460
 LorentzVector.h:461
 LorentzVector.h:462
 LorentzVector.h:463
 LorentzVector.h:464
 LorentzVector.h:465
 LorentzVector.h:466
 LorentzVector.h:467
 LorentzVector.h:468
 LorentzVector.h:469
 LorentzVector.h:470
 LorentzVector.h:471
 LorentzVector.h:472
 LorentzVector.h:473
 LorentzVector.h:474
 LorentzVector.h:475
 LorentzVector.h:476
 LorentzVector.h:477
 LorentzVector.h:478
 LorentzVector.h:479
 LorentzVector.h:480
 LorentzVector.h:481
 LorentzVector.h:482
 LorentzVector.h:483
 LorentzVector.h:484
 LorentzVector.h:485
 LorentzVector.h:486
 LorentzVector.h:487
 LorentzVector.h:488
 LorentzVector.h:489
 LorentzVector.h:490
 LorentzVector.h:491
 LorentzVector.h:492
 LorentzVector.h:493
 LorentzVector.h:494
 LorentzVector.h:495
 LorentzVector.h:496
 LorentzVector.h:497
 LorentzVector.h:498
 LorentzVector.h:499
 LorentzVector.h:500
 LorentzVector.h:501
 LorentzVector.h:502
 LorentzVector.h:503
 LorentzVector.h:504
 LorentzVector.h:505
 LorentzVector.h:506
 LorentzVector.h:507
 LorentzVector.h:508
 LorentzVector.h:509
 LorentzVector.h:510
 LorentzVector.h:511
 LorentzVector.h:512
 LorentzVector.h:513
 LorentzVector.h:514
 LorentzVector.h:515
 LorentzVector.h:516
 LorentzVector.h:517
 LorentzVector.h:518
 LorentzVector.h:519
 LorentzVector.h:520
 LorentzVector.h:521
 LorentzVector.h:522
 LorentzVector.h:523
 LorentzVector.h:524
 LorentzVector.h:525
 LorentzVector.h:526
 LorentzVector.h:527
 LorentzVector.h:528
 LorentzVector.h:529
 LorentzVector.h:530
 LorentzVector.h:531
 LorentzVector.h:532
 LorentzVector.h:533
 LorentzVector.h:534
 LorentzVector.h:535
 LorentzVector.h:536
 LorentzVector.h:537
 LorentzVector.h:538
 LorentzVector.h:539
 LorentzVector.h:540
 LorentzVector.h:541
 LorentzVector.h:542
 LorentzVector.h:543
 LorentzVector.h:544
 LorentzVector.h:545
 LorentzVector.h:546
 LorentzVector.h:547
 LorentzVector.h:548
 LorentzVector.h:549
 LorentzVector.h:550
 LorentzVector.h:551
 LorentzVector.h:552
 LorentzVector.h:553
 LorentzVector.h:554
 LorentzVector.h:555
 LorentzVector.h:556
 LorentzVector.h:557
 LorentzVector.h:558
 LorentzVector.h:559
 LorentzVector.h:560
 LorentzVector.h:561
 LorentzVector.h:562
 LorentzVector.h:563
 LorentzVector.h:564
 LorentzVector.h:565
 LorentzVector.h:566
 LorentzVector.h:567
 LorentzVector.h:568
 LorentzVector.h:569
 LorentzVector.h:570
 LorentzVector.h:571
 LorentzVector.h:572
 LorentzVector.h:573
 LorentzVector.h:574
 LorentzVector.h:575
 LorentzVector.h:576
 LorentzVector.h:577
 LorentzVector.h:578
 LorentzVector.h:579
 LorentzVector.h:580
 LorentzVector.h:581
 LorentzVector.h:582
 LorentzVector.h:583
 LorentzVector.h:584
 LorentzVector.h:585
 LorentzVector.h:586
 LorentzVector.h:587
 LorentzVector.h:588
 LorentzVector.h:589
 LorentzVector.h:590
 LorentzVector.h:591
 LorentzVector.h:592
 LorentzVector.h:593
 LorentzVector.h:594
 LorentzVector.h:595
 LorentzVector.h:596
 LorentzVector.h:597
 LorentzVector.h:598
 LorentzVector.h:599
 LorentzVector.h:600
 LorentzVector.h:601
 LorentzVector.h:602
 LorentzVector.h:603
 LorentzVector.h:604
 LorentzVector.h:605
 LorentzVector.h:606
 LorentzVector.h:607
 LorentzVector.h:608
 LorentzVector.h:609
 LorentzVector.h:610
 LorentzVector.h:611
 LorentzVector.h:612
 LorentzVector.h:613
 LorentzVector.h:614
 LorentzVector.h:615
 LorentzVector.h:616
 LorentzVector.h:617
 LorentzVector.h:618
 LorentzVector.h:619
 LorentzVector.h:620
 LorentzVector.h:621
 LorentzVector.h:622
 LorentzVector.h:623
 LorentzVector.h:624
 LorentzVector.h:625
 LorentzVector.h:626
 LorentzVector.h:627
 LorentzVector.h:628
 LorentzVector.h:629
 LorentzVector.h:630
 LorentzVector.h:631
 LorentzVector.h:632
 LorentzVector.h:633
 LorentzVector.h:634
 LorentzVector.h:635
 LorentzVector.h:636
 LorentzVector.h:637
 LorentzVector.h:638
 LorentzVector.h:639
 LorentzVector.h:640
 LorentzVector.h:641
 LorentzVector.h:642
 LorentzVector.h:643
 LorentzVector.h:644
 LorentzVector.h:645
 LorentzVector.h:646
 LorentzVector.h:647
 LorentzVector.h:648
 LorentzVector.h:649
 LorentzVector.h:650
 LorentzVector.h:651
 LorentzVector.h:652
 LorentzVector.h:653
 LorentzVector.h:654
 LorentzVector.h:655
 LorentzVector.h:656
 LorentzVector.h:657
 LorentzVector.h:658
 LorentzVector.h:659
 LorentzVector.h:660
 LorentzVector.h:661
 LorentzVector.h:662
 LorentzVector.h:663
 LorentzVector.h:664
 LorentzVector.h:665
 LorentzVector.h:666
 LorentzVector.h:667
 LorentzVector.h:668
 LorentzVector.h:669
 LorentzVector.h:670
 LorentzVector.h:671
 LorentzVector.h:672
 LorentzVector.h:673
 LorentzVector.h:674
 LorentzVector.h:675
 LorentzVector.h:676
 LorentzVector.h:677
 LorentzVector.h:678
 LorentzVector.h:679
 LorentzVector.h:680
 LorentzVector.h:681
 LorentzVector.h:682
 LorentzVector.h:683
 LorentzVector.h:684
 LorentzVector.h:685
 LorentzVector.h:686
 LorentzVector.h:687
 LorentzVector.h:688
 LorentzVector.h:689
 LorentzVector.h:690
 LorentzVector.h:691
 LorentzVector.h:692
 LorentzVector.h:693
 LorentzVector.h:694
 LorentzVector.h:695
 LorentzVector.h:696
 LorentzVector.h:697
 LorentzVector.h:698
 LorentzVector.h:699
 LorentzVector.h:700
 LorentzVector.h:701
 LorentzVector.h:702
 LorentzVector.h:703
 LorentzVector.h:704
 LorentzVector.h:705
 LorentzVector.h:706
 LorentzVector.h:707
 LorentzVector.h:708
 LorentzVector.h:709
 LorentzVector.h:710
 LorentzVector.h:711
 LorentzVector.h:712
 LorentzVector.h:713
 LorentzVector.h:714
 LorentzVector.h:715
 LorentzVector.h:716
 LorentzVector.h:717
 LorentzVector.h:718
 LorentzVector.h:719
 LorentzVector.h:720
 LorentzVector.h:721
 LorentzVector.h:722
 LorentzVector.h:723
 LorentzVector.h:724
 LorentzVector.h:725
 LorentzVector.h:726
 LorentzVector.h:727
 LorentzVector.h:728
 LorentzVector.h:729
 LorentzVector.h:730
 LorentzVector.h:731
 LorentzVector.h:732
 LorentzVector.h:733
 LorentzVector.h:734
 LorentzVector.h:735
 LorentzVector.h:736
 LorentzVector.h:737
 LorentzVector.h:738
 LorentzVector.h:739
 LorentzVector.h:740
 LorentzVector.h:741
 LorentzVector.h:742
 LorentzVector.h:743
 LorentzVector.h:744
 LorentzVector.h:745
 LorentzVector.h:746
 LorentzVector.h:747
 LorentzVector.h:748
 LorentzVector.h:749
 LorentzVector.h:750
 LorentzVector.h:751
 LorentzVector.h:752
 LorentzVector.h:753
 LorentzVector.h:754
 LorentzVector.h:755
 LorentzVector.h:756
 LorentzVector.h:757
 LorentzVector.h:758
 LorentzVector.h:759
 LorentzVector.h:760
 LorentzVector.h:761
 LorentzVector.h:762
 LorentzVector.h:763
 LorentzVector.h:764
 LorentzVector.h:765
 LorentzVector.h:766
 LorentzVector.h:767
 LorentzVector.h:768
 LorentzVector.h:769