Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
LorentzVector.h
Go to the documentation of this file.
1 // @(#)root/mathcore:$Id$
2 // Authors: W. Brown, M. Fischler, L. Moneta 2005
3 
4 /**********************************************************************
5  * *
6  * Copyright (c) 2005 , LCG ROOT MathLib Team *
7  * *
8  * *
9  **********************************************************************/
10 
11 // Header file for class LorentzVector
12 //
13 // Created by: moneta at Tue May 31 17:06:09 2005
14 // Major mods by: fischler at Wed Jul 20 2005
15 //
16 // Last update: $Id$
17 //
18 #ifndef ROOT_Math_GenVector_LorentzVector
19 #define ROOT_Math_GenVector_LorentzVector 1
20 
22 
24 
26 
27 #include <cmath>
28 
29 namespace ROOT {
30 
31  namespace Math {
32 
33 //__________________________________________________________________________________________
34  /**
35  Class describing a generic LorentzVector in the 4D space-time,
36  using the specified coordinate system for the spatial vector part.
37  The metric used for the LorentzVector is (-,-,-,+).
38  In the case of LorentzVector we don't distinguish the concepts
39  of points and displacement vectors as in the 3D case,
40  since the main use case for 4D Vectors is to describe the kinematics of
41  relativistic particles. A LorentzVector behaves like a
42  DisplacementVector in 4D. The Minkowski components could be viewed as
43  v and t, or for kinematic 4-vectors, as p and E.
44 
45  ROOT provides specialisations and aliases to them of the ROOT::Math::LorentzVector template:
46  - ROOT::Math::PtEtaPhiMVector based on pt (rho),eta,phi and M (t) coordinates in double precision
47  - ROOT::Math::PtEtaPhiEVector based on pt (rho),eta,phi and E (t) coordinates in double precision
48  - ROOT::Math::PxPyPzMVector based on px,py,pz and M (mass) coordinates in double precision
49  - ROOT::Math::PxPyPzEVector based on px,py,pz and E (energy) coordinates in double precision
50  - ROOT::Math::XYZTVector based on x,y,z,t coordinates (cartesian) in double precision (same as PxPyPzEVector)
51  - ROOT::Math::XYZTVectorF based on x,y,z,t coordinates (cartesian) in float precision (same as PxPyPzEVector but float)
52 
53 More details about the GenVector package can be found [here](Vector.html).
54 
55 
56  @ingroup GenVector
57  */
58  template< class CoordSystem >
59  class LorentzVector {
60 
61  public:
62 
63  // ------ ctors ------
64 
65  typedef typename CoordSystem::Scalar Scalar;
66  typedef CoordSystem CoordinateType;
67 
68  /**
69  default constructor of an empty vector (Px = Py = Pz = E = 0 )
70  */
71  LorentzVector ( ) : fCoordinates() { }
72 
73  /**
74  generic constructors from four scalar values.
75  The association between values and coordinate depends on the
76  coordinate system. For PxPyPzE4D,
77  \param a scalar value (Px)
78  \param b scalar value (Py)
79  \param c scalar value (Pz)
80  \param d scalar value (E)
81  */
82  LorentzVector(const Scalar & a,
83  const Scalar & b,
84  const Scalar & c,
85  const Scalar & d) :
86  fCoordinates(a , b, c, d) { }
87 
88  /**
89  constructor from a LorentzVector expressed in different
90  coordinates, or using a different Scalar type
91  */
92  template< class Coords >
93  explicit LorentzVector(const LorentzVector<Coords> & v ) :
94  fCoordinates( v.Coordinates() ) { }
95 
96  /**
97  Construct from a foreign 4D vector type, for example, HepLorentzVector
98  Precondition: v must implement methods x(), y(), z(), and t()
99  */
100  template<class ForeignLorentzVector>
101  explicit LorentzVector( const ForeignLorentzVector & v) :
102  fCoordinates(PxPyPzE4D<Scalar>( v.x(), v.y(), v.z(), v.t() ) ) { }
103 
104 #ifdef LATER
105  /**
106  construct from a generic linear algebra vector implementing operator []
107  and with a size of at least 4. This could be also a C array
108  In this case v[0] is the first data member
109  ( Px for a PxPyPzE4D base)
110  \param v LA vector
111  \param index0 index of first vector element (Px)
112  */
113  template< class LAVector >
114  explicit LorentzVector(const LAVector & v, size_t index0 ) {
115  fCoordinates = CoordSystem ( v[index0], v[index0+1], v[index0+2], v[index0+3] );
116  }
117 #endif
118 
119 
120  // ------ assignment ------
121 
122  /**
123  Assignment operator from a lorentz vector of arbitrary type
124  */
125  template< class OtherCoords >
126  LorentzVector & operator= ( const LorentzVector<OtherCoords> & v) {
127  fCoordinates = v.Coordinates();
128  return *this;
129  }
130 
131  /**
132  assignment from any other Lorentz vector implementing
133  x(), y(), z() and t()
134  */
135  template<class ForeignLorentzVector>
136  LorentzVector & operator = ( const ForeignLorentzVector & v) {
137  SetXYZT( v.x(), v.y(), v.z(), v.t() );
138  return *this;
139  }
140 
141 #ifdef LATER
142  /**
143  assign from a generic linear algebra vector implementing operator []
144  and with a size of at least 4
145  In this case v[0] is the first data member
146  ( Px for a PxPyPzE4D base)
147  \param v LA vector
148  \param index0 index of first vector element (Px)
149  */
150  template< class LAVector >
151  LorentzVector & AssignFrom(const LAVector & v, size_t index0=0 ) {
152  fCoordinates.SetCoordinates( v[index0], v[index0+1], v[index0+2], v[index0+3] );
153  return *this;
154  }
155 #endif
156 
157  // ------ Set, Get, and access coordinate data ------
158 
159  /**
160  Retrieve a const reference to the coordinates object
161  */
162  const CoordSystem & Coordinates() const {
163  return fCoordinates;
164  }
165 
166  /**
167  Set internal data based on an array of 4 Scalar numbers
168  */
169  LorentzVector<CoordSystem>& SetCoordinates( const Scalar src[] ) {
170  fCoordinates.SetCoordinates(src);
171  return *this;
172  }
173 
174  /**
175  Set internal data based on 4 Scalar numbers
176  */
177  LorentzVector<CoordSystem>& SetCoordinates( Scalar a, Scalar b, Scalar c, Scalar d ) {
178  fCoordinates.SetCoordinates(a, b, c, d);
179  return *this;
180  }
181 
182  /**
183  Set internal data based on 4 Scalars at *begin to *end
184  */
185 //#ifdef NDEBUG
186  //this does not compile in CINT
187 // template< class IT >
188 // LorentzVector<CoordSystem>& SetCoordinates( IT begin, IT /* end */ ) {
189 // #endif
190  template< class IT >
191 #ifndef NDEBUG
192  LorentzVector<CoordSystem>& SetCoordinates( IT begin, IT end ) {
193 #else
194  LorentzVector<CoordSystem>& SetCoordinates( IT begin, IT /* end */ ) {
195 #endif
196  IT a = begin; IT b = ++begin; IT c = ++begin; IT d = ++begin;
197  assert (++begin==end);
198  SetCoordinates (*a,*b,*c,*d);
199  return *this;
200  }
201 
202  /**
203  get internal data into 4 Scalar numbers
204  */
205  void GetCoordinates( Scalar& a, Scalar& b, Scalar& c, Scalar & d ) const
206  { fCoordinates.GetCoordinates(a, b, c, d); }
207 
208  /**
209  get internal data into an array of 4 Scalar numbers
210  */
211  void GetCoordinates( Scalar dest[] ) const
212  { fCoordinates.GetCoordinates(dest); }
213 
214  /**
215  get internal data into 4 Scalars at *begin to *end
216  */
217  template <class IT>
218 #ifndef NDEBUG
219  void GetCoordinates( IT begin, IT end ) const
220 #else
221  void GetCoordinates( IT begin, IT /* end */ ) const
222 #endif
223  { IT a = begin; IT b = ++begin; IT c = ++begin; IT d = ++begin;
224  assert (++begin==end);
225  GetCoordinates (*a,*b,*c,*d);
226  }
227 
228  /**
229  get internal data into 4 Scalars at *begin
230  */
231  template <class IT>
232  void GetCoordinates( IT begin ) const {
233  Scalar a,b,c,d = 0;
234  GetCoordinates (a,b,c,d);
235  *begin++ = a;
236  *begin++ = b;
237  *begin++ = c;
238  *begin = d;
239  }
240 
241  /**
242  set the values of the vector from the cartesian components (x,y,z,t)
243  (if the vector is held in another coordinates, like (Pt,eta,phi,m)
244  then (x, y, z, t) are converted to that form)
245  */
246  LorentzVector<CoordSystem>& SetXYZT (Scalar xx, Scalar yy, Scalar zz, Scalar tt) {
247  fCoordinates.SetPxPyPzE(xx,yy,zz,tt);
248  return *this;
249  }
250  LorentzVector<CoordSystem>& SetPxPyPzE (Scalar xx, Scalar yy, Scalar zz, Scalar ee) {
251  fCoordinates.SetPxPyPzE(xx,yy,zz,ee);
252  return *this;
253  }
254 
255  // ------------------- Equality -----------------
256 
257  /**
258  Exact equality
259  */
260  bool operator==(const LorentzVector & rhs) const {
261  return fCoordinates==rhs.fCoordinates;
262  }
263  bool operator!= (const LorentzVector & rhs) const {
264  return !(operator==(rhs));
265  }
266 
267  // ------ Individual element access, in various coordinate systems ------
268 
269  // individual coordinate accessors in various coordinate systems
270 
271  /**
272  spatial X component
273  */
274  Scalar Px() const { return fCoordinates.Px(); }
275  Scalar X() const { return fCoordinates.Px(); }
276  /**
277  spatial Y component
278  */
279  Scalar Py() const { return fCoordinates.Py(); }
280  Scalar Y() const { return fCoordinates.Py(); }
281  /**
282  spatial Z component
283  */
284  Scalar Pz() const { return fCoordinates.Pz(); }
285  Scalar Z() const { return fCoordinates.Pz(); }
286  /**
287  return 4-th component (time, or energy for a 4-momentum vector)
288  */
289  Scalar E() const { return fCoordinates.E(); }
290  Scalar T() const { return fCoordinates.E(); }
291  /**
292  return magnitude (mass) squared M2 = T**2 - X**2 - Y**2 - Z**2
293  (we use -,-,-,+ metric)
294  */
295  Scalar M2() const { return fCoordinates.M2(); }
296  /**
297  return magnitude (mass) using the (-,-,-,+) metric.
298  If M2 is negative (space-like vector) a GenVector_exception
299  is suggested and if continuing, - sqrt( -M2) is returned
300  */
301  Scalar M() const { return fCoordinates.M();}
302  /**
303  return the spatial (3D) magnitude ( sqrt(X**2 + Y**2 + Z**2) )
304  */
305  Scalar R() const { return fCoordinates.R(); }
306  Scalar P() const { return fCoordinates.R(); }
307  /**
308  return the square of the spatial (3D) magnitude ( X**2 + Y**2 + Z**2 )
309  */
310  Scalar P2() const { return P() * P(); }
311  /**
312  return the square of the transverse spatial component ( X**2 + Y**2 )
313  */
314  Scalar Perp2( ) const { return fCoordinates.Perp2();}
315 
316  /**
317  return the transverse spatial component sqrt ( X**2 + Y**2 )
318  */
319  Scalar Pt() const { return fCoordinates.Pt(); }
320  Scalar Rho() const { return fCoordinates.Pt(); }
321 
322  /**
323  return the transverse mass squared
324  \f[ m_t^2 = E^2 - p{_z}^2 \f]
325  */
326  Scalar Mt2() const { return fCoordinates.Mt2(); }
327 
328  /**
329  return the transverse mass
330  \f[ \sqrt{ m_t^2 = E^2 - p{_z}^2} X sign(E^ - p{_z}^2) \f]
331  */
332  Scalar Mt() const { return fCoordinates.Mt(); }
333 
334  /**
335  return the transverse energy squared
336  \f[ e_t = \frac{E^2 p_{\perp}^2 }{ |p|^2 } \f]
337  */
338  Scalar Et2() const { return fCoordinates.Et2(); }
339 
340  /**
341  return the transverse energy
342  \f[ e_t = \sqrt{ \frac{E^2 p_{\perp}^2 }{ |p|^2 } } X sign(E) \f]
343  */
344  Scalar Et() const { return fCoordinates.Et(); }
345 
346  /**
347  azimuthal Angle
348  */
349  Scalar Phi() const { return fCoordinates.Phi();}
350 
351  /**
352  polar Angle
353  */
354  Scalar Theta() const { return fCoordinates.Theta(); }
355 
356  /**
357  pseudorapidity
358  \f[ \eta = - \ln { \tan { \frac { \theta} {2} } } \f]
359  */
360  Scalar Eta() const { return fCoordinates.Eta(); }
361 
362  /**
363  get the spatial components of the Vector in a
364  DisplacementVector based on Cartesian Coordinates
365  */
366  ::ROOT::Math::DisplacementVector3D<Cartesian3D<Scalar> > Vect() const {
367  return ::ROOT::Math::DisplacementVector3D<Cartesian3D<Scalar> >( X(), Y(), Z() );
368  }
369 
370  // ------ Operations combining two Lorentz vectors ------
371 
372  /**
373  scalar (Dot) product of two LorentzVector vectors (metric is -,-,-,+)
374  Enable the product using any other LorentzVector implementing
375  the x(), y() , y() and t() member functions
376  \param q any LorentzVector implementing the x(), y() , z() and t()
377  member functions
378  \return the result of v.q of type according to the base scalar type of v
379  */
380 
381  template< class OtherLorentzVector >
382  Scalar Dot(const OtherLorentzVector & q) const {
383  return t()*q.t() - x()*q.x() - y()*q.y() - z()*q.z();
384  }
385 
386  /**
387  Self addition with another Vector ( v+= q )
388  Enable the addition with any other LorentzVector
389  \param q any LorentzVector implementing the x(), y() , z() and t()
390  member functions
391  */
392  template< class OtherLorentzVector >
393  inline LorentzVector & operator += ( const OtherLorentzVector & q)
394  {
395  SetXYZT( x() + q.x(), y() + q.y(), z() + q.z(), t() + q.t() );
396  return *this;
397  }
398 
399  /**
400  Self subtraction of another Vector from this ( v-= q )
401  Enable the addition with any other LorentzVector
402  \param q any LorentzVector implementing the x(), y() , z() and t()
403  member functions
404  */
405  template< class OtherLorentzVector >
406  LorentzVector & operator -= ( const OtherLorentzVector & q) {
407  SetXYZT( x() - q.x(), y() - q.y(), z() - q.z(), t() - q.t() );
408  return *this;
409  }
410 
411  /**
412  addition of two LorentzVectors (v3 = v1 + v2)
413  Enable the addition with any other LorentzVector
414  \param v2 any LorentzVector implementing the x(), y() , z() and t()
415  member functions
416  \return a new LorentzVector of the same type as v1
417  */
418  template<class OtherLorentzVector>
419  LorentzVector operator + ( const OtherLorentzVector & v2) const
420  {
421  LorentzVector<CoordinateType> v3(*this);
422  v3 += v2;
423  return v3;
424  }
425 
426  /**
427  subtraction of two LorentzVectors (v3 = v1 - v2)
428  Enable the subtraction of any other LorentzVector
429  \param v2 any LorentzVector implementing the x(), y() , z() and t()
430  member functions
431  \return a new LorentzVector of the same type as v1
432  */
433  template<class OtherLorentzVector>
434  LorentzVector operator - ( const OtherLorentzVector & v2) const {
435  LorentzVector<CoordinateType> v3(*this);
436  v3 -= v2;
437  return v3;
438  }
439 
440  //--- scaling operations ------
441 
442  /**
443  multiplication by a scalar quantity v *= a
444  */
445  LorentzVector & operator *= (Scalar a) {
446  fCoordinates.Scale(a);
447  return *this;
448  }
449 
450  /**
451  division by a scalar quantity v /= a
452  */
453  LorentzVector & operator /= (Scalar a) {
454  fCoordinates.Scale(1/a);
455  return *this;
456  }
457 
458  /**
459  product of a LorentzVector by a scalar quantity
460  \param a scalar quantity of type a
461  \return a new mathcoreLorentzVector q = v * a same type as v
462  */
463  LorentzVector operator * ( const Scalar & a) const {
464  LorentzVector tmp(*this);
465  tmp *= a;
466  return tmp;
467  }
468 
469  /**
470  Divide a LorentzVector by a scalar quantity
471  \param a scalar quantity of type a
472  \return a new mathcoreLorentzVector q = v / a same type as v
473  */
474  LorentzVector<CoordSystem> operator / ( const Scalar & a) const {
475  LorentzVector<CoordSystem> tmp(*this);
476  tmp /= a;
477  return tmp;
478  }
479 
480  /**
481  Negative of a LorentzVector (q = - v )
482  \return a new LorentzVector with opposite direction and time
483  */
484  LorentzVector operator - () const {
485  //LorentzVector<CoordinateType> v(*this);
486  //v.Negate();
487  return operator*( Scalar(-1) );
488  }
489  LorentzVector operator + () const {
490  return *this;
491  }
492 
493  // ---- Relativistic Properties ----
494 
495  /**
496  Rapidity relative to the Z axis: .5 log [(E+Pz)/(E-Pz)]
497  */
498  Scalar Rapidity() const {
499  // TODO - It would be good to check that E > Pz and use the Throw()
500  // mechanism or at least load a NAN if not.
501  // We should then move the code to a .cpp file.
502  const Scalar ee = E();
503  const Scalar ppz = Pz();
504  return Scalar(0.5) * log((ee + ppz) / (ee - ppz));
505  }
506 
507  /**
508  Rapidity in the direction of travel: atanh (|P|/E)=.5 log[(E+P)/(E-P)]
509  */
510  Scalar ColinearRapidity() const {
511  // TODO - It would be good to check that E > P and use the Throw()
512  // mechanism or at least load a NAN if not.
513  const Scalar ee = E();
514  const Scalar pp = P();
515  return Scalar(0.5) * log((ee + pp) / (ee - pp));
516  }
517 
518  /**
519  Determine if momentum-energy can represent a physical massive particle
520  */
521  bool isTimelike( ) const {
522  Scalar ee = E(); Scalar pp = P(); return ee*ee > pp*pp;
523  }
524 
525  /**
526  Determine if momentum-energy can represent a massless particle
527  */
528  bool isLightlike( Scalar tolerance
529  = 100*std::numeric_limits<Scalar>::epsilon() ) const {
530  Scalar ee = E(); Scalar pp = P(); Scalar delta = ee-pp;
531  if ( ee==0 ) return pp==0;
532  return delta*delta < tolerance * ee*ee;
533  }
534 
535  /**
536  Determine if momentum-energy is spacelike, and represents a tachyon
537  */
538  bool isSpacelike( ) const {
539  Scalar ee = E(); Scalar pp = P(); return ee*ee < pp*pp;
540  }
541 
542  typedef DisplacementVector3D< Cartesian3D<Scalar> > BetaVector;
543 
544  /**
545  The beta vector for the boost that would bring this vector into
546  its center of mass frame (zero momentum)
547  */
548  BetaVector BoostToCM( ) const {
549  if (E() == 0) {
550  if (P() == 0) {
551  return BetaVector();
552  } else {
553  // TODO - should attempt to Throw with msg about
554  // boostVector computed for LorentzVector with t=0
555  return -Vect()/E();
556  }
557  }
558  if (M2() <= 0) {
559  // TODO - should attempt to Throw with msg about
560  // boostVector computed for a non-timelike LorentzVector
561  }
562  return -Vect()/E();
563  }
564 
565  /**
566  The beta vector for the boost that would bring this vector into
567  its center of mass frame (zero momentum)
568  */
569  template <class Other4Vector>
570  BetaVector BoostToCM(const Other4Vector& v ) const {
571  Scalar eSum = E() + v.E();
572  DisplacementVector3D< Cartesian3D<Scalar> > vecSum = Vect() + v.Vect();
573  if (eSum == 0) {
574  if (vecSum.Mag2() == 0) {
575  return BetaVector();
576  } else {
577  // TODO - should attempt to Throw with msg about
578  // boostToCM computed for two 4-vectors with combined t=0
579  return BetaVector(vecSum/eSum);
580  }
581  // TODO - should attempt to Throw with msg about
582  // boostToCM computed for two 4-vectors with combined e=0
583  }
584  return BetaVector (vecSum * (-1./eSum));
585  }
586 
587  //beta and gamma
588 
589  /**
590  Return beta scalar value
591  */
592  Scalar Beta() const {
593  if ( E() == 0 ) {
594  if ( P2() == 0)
595  // to avoid Nan
596  return 0;
597  else {
598  GenVector::Throw ("LorentzVector::Beta() - beta computed for LorentzVector with t = 0. Return an Infinite result");
599  return 1./E();
600  }
601  }
602  if ( M2() <= 0 ) {
603  GenVector::Throw ("LorentzVector::Beta() - beta computed for non-timelike LorentzVector . Result is physically meaningless" );
604  }
605  return P() / E();
606  }
607  /**
608  Return Gamma scalar value
609  */
610  Scalar Gamma() const {
611  const Scalar v2 = P2();
612  const Scalar t2 = E() * E();
613  if (E() == 0) {
614  if ( P2() == 0) {
615  return 1;
616  } else {
617  GenVector::Throw ("LorentzVector::Gamma() - gamma computed for LorentzVector with t = 0. Return a zero result");
618 
619  }
620  }
621  if ( t2 < v2 ) {
622  GenVector::Throw ("LorentzVector::Gamma() - gamma computed for a spacelike LorentzVector. Imaginary result");
623  return 0;
624  }
625  else if ( t2 == v2 ) {
626  GenVector::Throw ("LorentzVector::Gamma() - gamma computed for a lightlike LorentzVector. Infinite result");
627  }
628  return Scalar(1) / sqrt(Scalar(1) - v2 / t2);
629  } /* gamma */
630 
631 
632  // Method providing limited backward name compatibility with CLHEP ----
633 
634  Scalar x() const { return fCoordinates.Px(); }
635  Scalar y() const { return fCoordinates.Py(); }
636  Scalar z() const { return fCoordinates.Pz(); }
637  Scalar t() const { return fCoordinates.E(); }
638  Scalar px() const { return fCoordinates.Px(); }
639  Scalar py() const { return fCoordinates.Py(); }
640  Scalar pz() const { return fCoordinates.Pz(); }
641  Scalar e() const { return fCoordinates.E(); }
642  Scalar r() const { return fCoordinates.R(); }
643  Scalar theta() const { return fCoordinates.Theta(); }
644  Scalar phi() const { return fCoordinates.Phi(); }
645  Scalar rho() const { return fCoordinates.Rho(); }
646  Scalar eta() const { return fCoordinates.Eta(); }
647  Scalar pt() const { return fCoordinates.Pt(); }
648  Scalar perp2() const { return fCoordinates.Perp2(); }
649  Scalar mag2() const { return fCoordinates.M2(); }
650  Scalar mag() const { return fCoordinates.M(); }
651  Scalar mt() const { return fCoordinates.Mt(); }
652  Scalar mt2() const { return fCoordinates.Mt2(); }
653 
654 
655  // Methods requested by CMS ---
656  Scalar energy() const { return fCoordinates.E(); }
657  Scalar mass() const { return fCoordinates.M(); }
658  Scalar mass2() const { return fCoordinates.M2(); }
659 
660 
661  /**
662  Methods setting a Single-component
663  Work only if the component is one of which the vector is represented.
664  For example SetE will work for a PxPyPzE Vector but not for a PxPyPzM Vector.
665  */
666  LorentzVector<CoordSystem>& SetE ( Scalar a ) { fCoordinates.SetE (a); return *this; }
667  LorentzVector<CoordSystem>& SetEta( Scalar a ) { fCoordinates.SetEta(a); return *this; }
668  LorentzVector<CoordSystem>& SetM ( Scalar a ) { fCoordinates.SetM (a); return *this; }
669  LorentzVector<CoordSystem>& SetPhi( Scalar a ) { fCoordinates.SetPhi(a); return *this; }
670  LorentzVector<CoordSystem>& SetPt ( Scalar a ) { fCoordinates.SetPt (a); return *this; }
671  LorentzVector<CoordSystem>& SetPx ( Scalar a ) { fCoordinates.SetPx (a); return *this; }
672  LorentzVector<CoordSystem>& SetPy ( Scalar a ) { fCoordinates.SetPy (a); return *this; }
673  LorentzVector<CoordSystem>& SetPz ( Scalar a ) { fCoordinates.SetPz (a); return *this; }
674 
675  private:
676 
677  CoordSystem fCoordinates; // internal coordinate system
678 
679 
680  }; // LorentzVector<>
681 
682 
683 
684  // global nethods
685 
686  /**
687  Scale of a LorentzVector with a scalar quantity a
688  \param a scalar quantity of typpe a
689  \param v mathcore::LorentzVector based on any coordinate system
690  \return a new mathcoreLorentzVector q = v * a same type as v
691  */
692  template< class CoordSystem >
693  inline LorentzVector<CoordSystem> operator *
694  ( const typename LorentzVector<CoordSystem>::Scalar & a,
695  const LorentzVector<CoordSystem>& v) {
696  LorentzVector<CoordSystem> tmp(v);
697  tmp *= a;
698  return tmp;
699  }
700 
701  // ------------- I/O to/from streams -------------
702 
703  template< class char_t, class traits_t, class Coords >
704  inline
705  std::basic_ostream<char_t,traits_t> &
706  operator << ( std::basic_ostream<char_t,traits_t> & os
707  , LorentzVector<Coords> const & v
708  )
709  {
710  if( !os ) return os;
711 
712  typename Coords::Scalar a, b, c, d;
713  v.GetCoordinates(a, b, c, d);
714 
715  if( detail::get_manip( os, detail::bitforbit ) ) {
716  detail::set_manip( os, detail::bitforbit, '\00' );
717  // TODO: call MF's bitwise-accurate functions on each of a, b, c, d
718  }
719  else {
720  os << detail::get_manip( os, detail::open ) << a
721  << detail::get_manip( os, detail::sep ) << b
722  << detail::get_manip( os, detail::sep ) << c
723  << detail::get_manip( os, detail::sep ) << d
724  << detail::get_manip( os, detail::close );
725  }
726 
727  return os;
728 
729  } // op<< <>()
730 
731 
732  template< class char_t, class traits_t, class Coords >
733  inline
734  std::basic_istream<char_t,traits_t> &
735  operator >> ( std::basic_istream<char_t,traits_t> & is
736  , LorentzVector<Coords> & v
737  )
738  {
739  if( !is ) return is;
740 
741  typename Coords::Scalar a, b, c, d;
742 
743  if( detail::get_manip( is, detail::bitforbit ) ) {
744  detail::set_manip( is, detail::bitforbit, '\00' );
745  // TODO: call MF's bitwise-accurate functions on each of a, b, c
746  }
747  else {
748  detail::require_delim( is, detail::open ); is >> a;
749  detail::require_delim( is, detail::sep ); is >> b;
750  detail::require_delim( is, detail::sep ); is >> c;
751  detail::require_delim( is, detail::sep ); is >> d;
752  detail::require_delim( is, detail::close );
753  }
754 
755  if( is )
756  v.SetCoordinates(a, b, c, d);
757  return is;
758 
759  } // op>> <>()
760 
761 
762 
763  } // end namespace Math
764 
765 } // end namespace ROOT
766 
767 #include <sstream>
768 namespace cling
769 {
770 template<typename CoordSystem>
771 std::string printValue(const ROOT::Math::LorentzVector<CoordSystem> *v)
772 {
773  std::stringstream s;
774  s << *v;
775  return s.str();
776 }
777 
778 } // end namespace cling
779 
780 #endif
781 
782 //#include "Math/GenVector/LorentzVectorOperations.h"
783 
784 
785