rivet is hosted by Hepforge, IPPP Durham
Rivet  2.7.0
Vector4.hh
1 #ifndef RIVET_MATH_VECTOR4
2 #define RIVET_MATH_VECTOR4
3 
4 #include "Rivet/Math/MathHeader.hh"
5 #include "Rivet/Math/MathUtils.hh"
6 #include "Rivet/Math/VectorN.hh"
7 #include "Rivet/Math/Vector3.hh"
8 
9 namespace Rivet {
10 
11 
12  class FourVector;
13  class FourMomentum;
14  class LorentzTransform;
15  typedef FourVector Vector4;
16  FourVector transform(const LorentzTransform& lt, const FourVector& v4);
17 
18 
22  class FourVector : public Vector<4> {
23  friend FourVector multiply(const double a, const FourVector& v);
24  friend FourVector multiply(const FourVector& v, const double a);
25  friend FourVector add(const FourVector& a, const FourVector& b);
26  friend FourVector transform(const LorentzTransform& lt, const FourVector& v4);
27 
28  public:
29 
30  FourVector() : Vector<4>() { }
31 
32  template<typename V4>
33  FourVector(const V4& other) {
34  this->setT(other.t());
35  this->setX(other.x());
36  this->setY(other.y());
37  this->setZ(other.z());
38  }
39 
40  FourVector(const Vector<4>& other)
41  : Vector<4>(other) { }
42 
43  FourVector(const double t, const double x, const double y, const double z) {
44  this->setT(t);
45  this->setX(x);
46  this->setY(y);
47  this->setZ(z);
48  }
49 
50  virtual ~FourVector() { }
51 
52  public:
53 
54  double t() const { return get(0); }
55  double t2() const { return sqr(t()); }
56  FourVector& setT(const double t) { set(0, t); return *this; }
57 
58  double x() const { return get(1); }
59  double x2() const { return sqr(x()); }
60  FourVector& setX(const double x) { set(1, x); return *this; }
61 
62  double y() const { return get(2); }
63  double y2() const { return sqr(y()); }
64  FourVector& setY(const double y) { set(2, y); return *this; }
65 
66  double z() const { return get(3); }
67  double z2() const { return sqr(z()); }
68  FourVector& setZ(const double z) { set(3, z); return *this; }
69 
70  double invariant() const {
71  // Done this way for numerical precision
72  return (t() + z())*(t() - z()) - x()*x() - y()*y();
73  }
74 
75  bool isNull() const {
76  return Rivet::isZero(invariant());
77  }
78 
80  double angle(const FourVector& v) const {
81  return vector3().angle( v.vector3() );
82  }
84  double angle(const Vector3& v3) const {
85  return vector3().angle(v3);
86  }
87 
91  double polarRadius2() const {
92  return vector3().polarRadius2();
93  }
95  double perp2() const {
96  return vector3().perp2();
97  }
99  double rho2() const {
100  return vector3().rho2();
101  }
102 
104  double polarRadius() const {
105  return vector3().polarRadius();
106  }
108  double perp() const {
109  return vector3().perp();
110  }
112  double rho() const {
113  return vector3().rho();
114  }
115 
117  Vector3 polarVec() const {
118  return vector3().polarVec();
119  }
121  Vector3 perpVec() const {
122  return vector3().perpVec();
123  }
125  Vector3 rhoVec() const {
126  return vector3().rhoVec();
127  }
128 
130  double azimuthalAngle(const PhiMapping mapping=ZERO_2PI) const {
131  return vector3().azimuthalAngle(mapping);
132  }
134  double phi(const PhiMapping mapping=ZERO_2PI) const {
135  return vector3().phi(mapping);
136  }
137 
139  double polarAngle() const {
140  return vector3().polarAngle();
141  }
143  double theta() const {
144  return vector3().theta();
145  }
146 
148  double pseudorapidity() const {
149  return vector3().pseudorapidity();
150  }
152  double eta() const {
153  return vector3().eta();
154  }
155 
157  double abspseudorapidity() const { return fabs(eta()); }
159  double abseta() const { return fabs(eta()); }
160 
162  Vector3 vector3() const {
163  return Vector3(get(1), get(2), get(3));
164  }
165 
167  operator Vector3 () const { return vector3(); }
168 
169 
170  public:
171 
173  double contract(const FourVector& v) const {
174  const double result = t()*v.t() - x()*v.x() - y()*v.y() - z()*v.z();
175  return result;
176  }
177 
179  double dot(const FourVector& v) const {
180  return contract(v);
181  }
182 
184  double operator*(const FourVector& v) const {
185  return contract(v);
186  }
187 
189  FourVector& operator*=(double a) {
190  _vec = multiply(a, *this)._vec;
191  return *this;
192  }
193 
195  FourVector& operator/=(double a) {
196  _vec = multiply(1.0/a, *this)._vec;
197  return *this;
198  }
199 
201  FourVector& operator+=(const FourVector& v) {
202  _vec = add(*this, v)._vec;
203  return *this;
204  }
205 
207  FourVector& operator-=(const FourVector& v) {
208  _vec = add(*this, -v)._vec;
209  return *this;
210  }
211 
213  FourVector operator-() const {
214  FourVector result;
215  result._vec = -_vec;
216  return result;
217  }
218 
220  FourVector reverse() const {
221  FourVector result = -*this;
222  result.setT(-result.t());
223  return result;
224  }
225 
226  };
227 
228 
230  inline double contract(const FourVector& a, const FourVector& b) {
231  return a.contract(b);
232  }
233 
235  inline double dot(const FourVector& a, const FourVector& b) {
236  return contract(a, b);
237  }
238 
239  inline FourVector multiply(const double a, const FourVector& v) {
240  FourVector result;
241  result._vec = a * v._vec;
242  return result;
243  }
244 
245  inline FourVector multiply(const FourVector& v, const double a) {
246  return multiply(a, v);
247  }
248 
249  inline FourVector operator*(const double a, const FourVector& v) {
250  return multiply(a, v);
251  }
252 
253  inline FourVector operator*(const FourVector& v, const double a) {
254  return multiply(a, v);
255  }
256 
257  inline FourVector operator/(const FourVector& v, const double a) {
258  return multiply(1.0/a, v);
259  }
260 
261  inline FourVector add(const FourVector& a, const FourVector& b) {
262  FourVector result;
263  result._vec = a._vec + b._vec;
264  return result;
265  }
266 
267  inline FourVector operator+(const FourVector& a, const FourVector& b) {
268  return add(a, b);
269  }
270 
271  inline FourVector operator-(const FourVector& a, const FourVector& b) {
272  return add(a, -b);
273  }
274 
277  inline double invariant(const FourVector& lv) {
278  return lv.invariant();
279  }
280 
282  inline double angle(const FourVector& a, const FourVector& b) {
283  return a.angle(b);
284  }
285 
287  inline double angle(const Vector3& a, const FourVector& b) {
288  return angle( a, b.vector3() );
289  }
290 
292  inline double angle(const FourVector& a, const Vector3& b) {
293  return a.angle(b);
294  }
295 
296 
298 
299 
301  class FourMomentum : public FourVector {
302  friend FourMomentum multiply(const double a, const FourMomentum& v);
303  friend FourMomentum multiply(const FourMomentum& v, const double a);
304  friend FourMomentum add(const FourMomentum& a, const FourMomentum& b);
305  friend FourMomentum transform(const LorentzTransform& lt, const FourMomentum& v4);
306 
307  public:
308  FourMomentum() { }
309 
310  template<typename V4>
311  FourMomentum(const V4& other) {
312  this->setE(other.t());
313  this->setPx(other.x());
314  this->setPy(other.y());
315  this->setPz(other.z());
316  }
317 
318  FourMomentum(const Vector<4>& other)
319  : FourVector(other) { }
320 
321  FourMomentum(const double E, const double px, const double py, const double pz) {
322  this->setE(E);
323  this->setPx(px);
324  this->setPy(py);
325  this->setPz(pz);
326  }
327 
328  ~FourMomentum() {}
329 
330  public:
331 
332 
334 
335 
337  FourMomentum& setE(double E) {
338  setT(E);
339  return *this;
340  }
341 
343  FourMomentum& setPx(double px) {
344  setX(px);
345  return *this;
346  }
347 
349  FourMomentum& setPy(double py) {
350  setY(py);
351  return *this;
352  }
353 
355  FourMomentum& setPz(double pz) {
356  setZ(pz);
357  return *this;
358  }
359 
360 
362  FourMomentum& setPE(double px, double py, double pz, double E) {
363  if (E < 0)
364  throw std::invalid_argument("Negative energy given as argument: " + to_str(E));
365  setPx(px); setPy(py); setPz(pz); setE(E);
366  return *this;
367  }
369  FourMomentum& setXYZE(double px, double py, double pz, double E) {
370  return setPE(px, py, pz, E);
371  }
372  // /// Near-alias with switched arg order
373  // FourMomentum& setEP(double E, double px, double py, double pz) {
374  // return setPE(px, py, pz, E);
375  // }
376  // /// Alias for setEP
377  // FourMomentum& setEXYZ(double E, double px, double py, double pz) {
378  // return setEP(E, px, py, pz);
379  // }
380 
381 
383  FourMomentum& setPM(double px, double py, double pz, double mass) {
384  if (mass < 0)
385  throw std::invalid_argument("Negative mass given as argument: " + to_str(mass));
386  const double E = sqrt( sqr(mass) + sqr(px) + sqr(py) + sqr(pz) );
387  // setPx(px); setPy(py); setPz(pz); setE(E);
388  return setPE(px, py, pz, E);
389  }
391  FourMomentum& setXYZM(double px, double py, double pz, double mass) {
392  return setPM(px, py, pz, mass);
393  }
394 
395 
400  FourMomentum& setEtaPhiME(double eta, double phi, double mass, double E) {
401  if (mass < 0)
402  throw std::invalid_argument("Negative mass given as argument");
403  if (E < 0)
404  throw std::invalid_argument("Negative energy given as argument");
405  const double theta = 2 * atan(exp(-eta));
406  if (theta < 0 || theta > M_PI)
407  throw std::domain_error("Polar angle outside 0..pi in calculation");
408  setThetaPhiME(theta, phi, mass, E);
409  return *this;
410  }
411 
416  FourMomentum& setEtaPhiMPt(double eta, double phi, double mass, double pt) {
417  if (mass < 0)
418  throw std::invalid_argument("Negative mass given as argument");
419  if (pt < 0)
420  throw std::invalid_argument("Negative transverse momentum given as argument");
421  const double theta = 2 * atan(exp(-eta));
422  if (theta < 0 || theta > M_PI)
423  throw std::domain_error("Polar angle outside 0..pi in calculation");
424  const double p = pt / sin(theta);
425  const double E = sqrt( sqr(p) + sqr(mass) );
426  setThetaPhiME(theta, phi, mass, E);
427  return *this;
428  }
429 
438  FourMomentum& setRapPhiME(double y, double phi, double mass, double E) {
439  if (mass < 0)
440  throw std::invalid_argument("Negative mass given as argument");
441  if (E < 0)
442  throw std::invalid_argument("Negative energy given as argument");
443  const double sqrt_pt2_m2 = E / cosh(y);
444  const double pt = sqrt( sqr(sqrt_pt2_m2) - sqr(mass) );
445  if (pt < 0)
446  throw std::domain_error("Negative transverse momentum in calculation");
447  const double pz = sqrt_pt2_m2 * sinh(y);
448  const double px = pt * cos(phi);
449  const double py = pt * sin(phi);
450  setPE(px, py, pz, E);
451  return *this;
452  }
453 
458  FourMomentum& setRapPhiMPt(double y, double phi, double mass, double pt) {
459  if (mass < 0)
460  throw std::invalid_argument("Negative mass given as argument");
461  if (pt < 0)
462  throw std::invalid_argument("Negative transverse mass given as argument");
463  const double E = sqrt( sqr(pt) + sqr(mass) ) * cosh(y);
464  if (E < 0)
465  throw std::domain_error("Negative energy in calculation");
466  setRapPhiME(y, phi, mass, E);
467  return *this;
468  }
469 
475  FourMomentum& setThetaPhiME(double theta, double phi, double mass, double E) {
476  if (theta < 0 || theta > M_PI)
477  throw std::invalid_argument("Polar angle outside 0..pi given as argument");
478  if (mass < 0)
479  throw std::invalid_argument("Negative mass given as argument");
480  if (E < 0)
481  throw std::invalid_argument("Negative energy given as argument");
482  const double p = sqrt( sqr(E) - sqr(mass) );
483  const double pz = p * cos(theta);
484  const double pt = p * sin(theta);
485  if (pt < 0)
486  throw std::invalid_argument("Negative transverse momentum in calculation");
487  const double px = pt * cos(phi);
488  const double py = pt * sin(phi);
489  setPE(px, py, pz, E);
490  return *this;
491  }
492 
498  FourMomentum& setThetaPhiMPt(double theta, double phi, double mass, double pt) {
499  if (theta < 0 || theta > M_PI)
500  throw std::invalid_argument("Polar angle outside 0..pi given as argument");
501  if (mass < 0)
502  throw std::invalid_argument("Negative mass given as argument");
503  if (pt < 0)
504  throw std::invalid_argument("Negative transverse momentum given as argument");
505  const double p = pt / sin(theta);
506  const double px = pt * cos(phi);
507  const double py = pt * sin(phi);
508  const double pz = p * cos(theta);
509  const double E = sqrt( sqr(p) + sqr(mass) );
510  setPE(px, py, pz, E);
511  return *this;
512  }
513 
517  FourMomentum& setPtPhiME(double pt, double phi, double mass, double E) {
518  if (pt < 0)
519  throw std::invalid_argument("Negative transverse momentum given as argument");
520  if (mass < 0)
521  throw std::invalid_argument("Negative mass given as argument");
522  if (E < 0)
523  throw std::invalid_argument("Negative energy given as argument");
524  const double px = pt * cos(phi);
525  const double py = pt * sin(phi);
526  const double pz = sqrt(sqr(E) - sqr(mass) - sqr(pt));
527  setPE(px, py, pz, E);
528  return *this;
529  }
530 
532 
533 
535 
536 
538  double E() const { return t(); }
540  double E2() const { return t2(); }
541 
543  double px() const { return x(); }
545  double px2() const { return x2(); }
546 
548  double py() const { return y(); }
550  double py2() const { return y2(); }
551 
553  double pz() const { return z(); }
555  double pz2() const { return z2(); }
556 
557 
561  double mass() const {
562  // assert(Rivet::isZero(mass2()) || mass2() > 0);
563  // if (Rivet::isZero(mass2())) {
564  // return 0.0;
565  // } else {
566  // return sqrt(mass2());
567  // }
568  return sign(mass2()) * sqrt(fabs(mass2()));
569  }
570 
572  double mass2() const {
573  return invariant();
574  }
575 
576 
578  Vector3 p3() const { return vector3(); }
579 
581  double p() const {
582  return p3().mod();
583  }
584 
586  double p2() const {
587  return p3().mod2();
588  }
589 
590 
592  double rapidity() const {
593  return 0.5 * std::log( (E() + pz()) / (E() - pz()) );
594  }
596  double rap() const {
597  return rapidity();
598  }
599 
601  double absrapidity() const {
602  return fabs(rapidity());
603  }
605  double absrap() const {
606  return fabs(rap());
607  }
608 
610  Vector3 pTvec() const {
611  return p3().polarVec();
612  }
614  Vector3 ptvec() const {
615  return pTvec();
616  }
617 
619  double pT2() const {
620  return vector3().polarRadius2();
621  }
623  double pt2() const {
624  return vector3().polarRadius2();
625  }
626 
628  double pT() const {
629  return sqrt(pT2());
630  }
632  double pt() const {
633  return sqrt(pT2());
634  }
635 
637  double Et2() const {
638  return Et() * Et();
639  }
641  double Et() const {
642  return E() * sin(polarAngle());
643  }
644 
646 
647 
649 
650 
653  double gamma() const {
654  return sqrt(E2()/mass2());
655  }
656 
659  Vector3 gammaVec() const {
660  return gamma() * p3().unit();
661  }
662 
665  double beta() const {
666  return p()/E();
667  }
668 
671  Vector3 betaVec() const {
672  // return Vector3(px()/E(), py()/E(), pz()/E());
673  return p3()/E();
674  }
675 
678  Vector3 boostVector() const { return betaVec(); }
679 
681 
682 
684 
685 
687 
688 
690  struct byEAscending {
691  bool operator()(const FourMomentum& left, const FourMomentum& right) const{
692  const double pt2left = left.E();
693  const double pt2right = right.E();
694  return pt2left < pt2right;
695  }
696 
697  bool operator()(const FourMomentum* left, const FourMomentum* right) const{
698  return (*this)(*left, *right);
699  }
700  };
701 
702 
704  struct byEDescending {
705  bool operator()(const FourMomentum& left, const FourMomentum& right) const{
706  return byEAscending()(right, left);
707  }
708 
709  bool operator()(const FourMomentum* left, const FourVector* right) const{
710  return (*this)(*left, *right);
711  }
712  };
713 
715 
716 
718 
719 
721 
722 
724  FourMomentum& operator*=(double a) {
725  _vec = multiply(a, *this)._vec;
726  return *this;
727  }
728 
730  FourMomentum& operator/=(double a) {
731  _vec = multiply(1.0/a, *this)._vec;
732  return *this;
733  }
734 
736  FourMomentum& operator+=(const FourMomentum& v) {
737  _vec = add(*this, v)._vec;
738  return *this;
739  }
740 
742  FourMomentum& operator-=(const FourMomentum& v) {
743  _vec = add(*this, -v)._vec;
744  return *this;
745  }
746 
748  FourMomentum operator-() const {
749  FourMomentum result;
750  result._vec = -_vec;
751  return result;
752  }
753 
755  FourMomentum reverse() const {
756  FourMomentum result = -*this;
757  result.setE(-result.E());
758  return result;
759  }
760 
762 
763 
765 
766 
768 
769 
771  static FourMomentum mkXYZE(double px, double py, double pz, double E) {
772  return FourMomentum().setPE(px, py, pz, E);
773  }
774 
776  static FourMomentum mkXYZM(double px, double py, double pz, double mass) {
777  return FourMomentum().setPM(px, py, pz, mass);
778  }
779 
781  static FourMomentum mkEtaPhiME(double eta, double phi, double mass, double E) {
782  return FourMomentum().setEtaPhiME(eta, phi, mass, E);
783  }
784 
786  static FourMomentum mkEtaPhiMPt(double eta, double phi, double mass, double pt) {
787  return FourMomentum().setEtaPhiMPt(eta, phi, mass, pt);
788  }
789 
791  static FourMomentum mkRapPhiME(double y, double phi, double mass, double E) {
792  return FourMomentum().setRapPhiME(y, phi, mass, E);
793  }
794 
796  static FourMomentum mkRapPhiMPt(double y, double phi, double mass, double pt) {
797  return FourMomentum().setRapPhiMPt(y, phi, mass, pt);
798  }
799 
801  static FourMomentum mkThetaPhiME(double theta, double phi, double mass, double E) {
802  return FourMomentum().setThetaPhiME(theta, phi, mass, E);
803  }
804 
806  static FourMomentum mkThetaPhiMPt(double theta, double phi, double mass, double pt) {
807  return FourMomentum().setThetaPhiMPt(theta, phi, mass, pt);
808  }
809 
811  static FourMomentum mkPtPhiME(double pt, double phi, double mass, double E) {
812  return FourMomentum().setPtPhiME(pt, phi, mass, E);
813  }
814 
816 
817 
818  };
819 
820 
821 
822  inline FourMomentum multiply(const double a, const FourMomentum& v) {
823  FourMomentum result;
824  result._vec = a * v._vec;
825  return result;
826  }
827 
828  inline FourMomentum multiply(const FourMomentum& v, const double a) {
829  return multiply(a, v);
830  }
831 
832  inline FourMomentum operator*(const double a, const FourMomentum& v) {
833  return multiply(a, v);
834  }
835 
836  inline FourMomentum operator*(const FourMomentum& v, const double a) {
837  return multiply(a, v);
838  }
839 
840  inline FourMomentum operator/(const FourMomentum& v, const double a) {
841  return multiply(1.0/a, v);
842  }
843 
844  inline FourMomentum add(const FourMomentum& a, const FourMomentum& b) {
845  FourMomentum result;
846  result._vec = a._vec + b._vec;
847  return result;
848  }
849 
850  inline FourMomentum operator+(const FourMomentum& a, const FourMomentum& b) {
851  return add(a, b);
852  }
853 
854  inline FourMomentum operator-(const FourMomentum& a, const FourMomentum& b) {
855  return add(a, -b);
856  }
857 
858 
860 
861 
863 
864 
873  inline double deltaR2(const FourVector& a, const FourVector& b,
874  RapScheme scheme=PSEUDORAPIDITY) {
875  switch (scheme) {
876  case PSEUDORAPIDITY :
877  return deltaR2(a.vector3(), b.vector3());
878  case RAPIDITY:
879  {
880  const FourMomentum* ma = dynamic_cast<const FourMomentum*>(&a);
881  const FourMomentum* mb = dynamic_cast<const FourMomentum*>(&b);
882  if (!ma || !mb) {
883  string err = "deltaR with scheme RAPIDITY can only be called with FourMomentum objects, not FourVectors";
884  throw std::runtime_error(err);
885  }
886  return deltaR2(*ma, *mb, scheme);
887  }
888  default:
889  throw std::runtime_error("The specified deltaR scheme is not yet implemented");
890  }
891  }
892 
901  inline double deltaR(const FourVector& a, const FourVector& b,
902  RapScheme scheme=PSEUDORAPIDITY) {
903  return sqrt(deltaR2(a, b, scheme));
904  }
905 
906 
907 
914  inline double deltaR2(const FourVector& v,
915  double eta2, double phi2,
916  RapScheme scheme=PSEUDORAPIDITY) {
917  switch (scheme) {
918  case PSEUDORAPIDITY :
919  return deltaR2(v.vector3(), eta2, phi2);
920  case RAPIDITY:
921  {
922  const FourMomentum* mv = dynamic_cast<const FourMomentum*>(&v);
923  if (!mv) {
924  string err = "deltaR with scheme RAPIDITY can only be called with FourMomentum objects, not FourVectors";
925  throw std::runtime_error(err);
926  }
927  return deltaR2(*mv, eta2, phi2, scheme);
928  }
929  default:
930  throw std::runtime_error("The specified deltaR scheme is not yet implemented");
931  }
932  }
933 
940  inline double deltaR(const FourVector& v,
941  double eta2, double phi2,
942  RapScheme scheme=PSEUDORAPIDITY) {
943  return sqrt(deltaR2(v, eta2, phi2, scheme));
944  }
945 
946 
953  inline double deltaR2(double eta1, double phi1,
954  const FourVector& v,
955  RapScheme scheme=PSEUDORAPIDITY) {
956  switch (scheme) {
957  case PSEUDORAPIDITY :
958  return deltaR2(eta1, phi1, v.vector3());
959  case RAPIDITY:
960  {
961  const FourMomentum* mv = dynamic_cast<const FourMomentum*>(&v);
962  if (!mv) {
963  string err = "deltaR with scheme RAPIDITY can only be called with FourMomentum objects, not FourVectors";
964  throw std::runtime_error(err);
965  }
966  return deltaR2(eta1, phi1, *mv, scheme);
967  }
968  default:
969  throw std::runtime_error("The specified deltaR scheme is not yet implemented");
970  }
971  }
972 
979  inline double deltaR(double eta1, double phi1,
980  const FourVector& v,
981  RapScheme scheme=PSEUDORAPIDITY) {
982  return sqrt(deltaR2(eta1, phi1, v, scheme));
983  }
984 
985 
992  inline double deltaR2(const FourMomentum& a, const FourMomentum& b,
993  RapScheme scheme=PSEUDORAPIDITY) {
994  switch (scheme) {
995  case PSEUDORAPIDITY:
996  return deltaR2(a.vector3(), b.vector3());
997  case RAPIDITY:
998  return deltaR2(a.rapidity(), a.azimuthalAngle(), b.rapidity(), b.azimuthalAngle());
999  default:
1000  throw std::runtime_error("The specified deltaR scheme is not yet implemented");
1001  }
1002  }
1003 
1010  inline double deltaR(const FourMomentum& a, const FourMomentum& b,
1011  RapScheme scheme=PSEUDORAPIDITY) {
1012  return sqrt(deltaR2(a, b, scheme));
1013  }
1014 
1015 
1021  inline double deltaR2(const FourMomentum& v,
1022  double eta2, double phi2,
1023  RapScheme scheme=PSEUDORAPIDITY) {
1024  switch (scheme) {
1025  case PSEUDORAPIDITY:
1026  return deltaR2(v.vector3(), eta2, phi2);
1027  case RAPIDITY:
1028  return deltaR2(v.rapidity(), v.azimuthalAngle(), eta2, phi2);
1029  default:
1030  throw std::runtime_error("The specified deltaR scheme is not yet implemented");
1031  }
1032  }
1033 
1039  inline double deltaR(const FourMomentum& v,
1040  double eta2, double phi2,
1041  RapScheme scheme=PSEUDORAPIDITY) {
1042  return sqrt(deltaR2(v, eta2, phi2, scheme));
1043  }
1044 
1045 
1051  inline double deltaR2(double eta1, double phi1,
1052  const FourMomentum& v,
1053  RapScheme scheme=PSEUDORAPIDITY) {
1054  switch (scheme) {
1055  case PSEUDORAPIDITY:
1056  return deltaR2(eta1, phi1, v.vector3());
1057  case RAPIDITY:
1058  return deltaR2(eta1, phi1, v.rapidity(), v.azimuthalAngle());
1059  default:
1060  throw std::runtime_error("The specified deltaR scheme is not yet implemented");
1061  }
1062  }
1063 
1069  inline double deltaR(double eta1, double phi1,
1070  const FourMomentum& v,
1071  RapScheme scheme=PSEUDORAPIDITY) {
1072  return sqrt(deltaR2(eta1, phi1, v, scheme));
1073  }
1074 
1075 
1081  inline double deltaR2(const FourMomentum& a, const FourVector& b,
1082  RapScheme scheme=PSEUDORAPIDITY) {
1083  switch (scheme) {
1084  case PSEUDORAPIDITY:
1085  return deltaR2(a.vector3(), b.vector3());
1086  case RAPIDITY:
1087  return deltaR2(a.rapidity(), a.azimuthalAngle(), FourMomentum(b).rapidity(), b.azimuthalAngle());
1088  default:
1089  throw std::runtime_error("The specified deltaR scheme is not yet implemented");
1090  }
1091  }
1092 
1098  inline double deltaR(const FourMomentum& a, const FourVector& b,
1099  RapScheme scheme=PSEUDORAPIDITY) {
1100  return sqrt(deltaR2(a, b, scheme));
1101  }
1102 
1103 
1109  inline double deltaR2(const FourVector& a, const FourMomentum& b,
1110  RapScheme scheme=PSEUDORAPIDITY) {
1111  return deltaR2(b, a, scheme); //< note reversed args
1112  }
1113 
1119  inline double deltaR(const FourVector& a, const FourMomentum& b,
1120  RapScheme scheme=PSEUDORAPIDITY) {
1121  return deltaR(b, a, scheme); //< note reversed args
1122  }
1123 
1124 
1127  inline double deltaR2(const FourMomentum& a, const Vector3& b) {
1128  return deltaR2(a.vector3(), b);
1129  }
1130 
1133  inline double deltaR(const FourMomentum& a, const Vector3& b) {
1134  return deltaR(a.vector3(), b);
1135  }
1136 
1139  inline double deltaR2(const Vector3& a, const FourMomentum& b) {
1140  return deltaR2(a, b.vector3());
1141  }
1142 
1145  inline double deltaR(const Vector3& a, const FourMomentum& b) {
1146  return deltaR(a, b.vector3());
1147  }
1148 
1151  inline double deltaR2(const FourVector& a, const Vector3& b) {
1152  return deltaR2(a.vector3(), b);
1153  }
1154 
1157  inline double deltaR(const FourVector& a, const Vector3& b) {
1158  return deltaR(a.vector3(), b);
1159  }
1160 
1163  inline double deltaR2(const Vector3& a, const FourVector& b) {
1164  return deltaR2(a, b.vector3());
1165  }
1166 
1169  inline double deltaR(const Vector3& a, const FourVector& b) {
1170  return deltaR(a, b.vector3());
1171  }
1172 
1174 
1175 
1177 
1178 
1180 
1181 
1183  inline double deltaPhi(const FourMomentum& a, const FourMomentum& b, bool sign=false) {
1184  return deltaPhi(a.vector3(), b.vector3(), sign);
1185  }
1186 
1188  inline double deltaPhi(const FourMomentum& v, double phi2, bool sign=false) {
1189  return deltaPhi(v.vector3(), phi2, sign);
1190  }
1191 
1193  inline double deltaPhi(double phi1, const FourMomentum& v, bool sign=false) {
1194  return deltaPhi(phi1, v.vector3(), sign);
1195  }
1196 
1198  inline double deltaPhi(const FourVector& a, const FourVector& b, bool sign=false) {
1199  return deltaPhi(a.vector3(), b.vector3(), sign);
1200  }
1201 
1203  inline double deltaPhi(const FourVector& v, double phi2, bool sign=false) {
1204  return deltaPhi(v.vector3(), phi2, sign);
1205  }
1206 
1208  inline double deltaPhi(double phi1, const FourVector& v, bool sign=false) {
1209  return deltaPhi(phi1, v.vector3(), sign);
1210  }
1211 
1213  inline double deltaPhi(const FourVector& a, const FourMomentum& b, bool sign=false) {
1214  return deltaPhi(a.vector3(), b.vector3(), sign);
1215  }
1216 
1218  inline double deltaPhi(const FourMomentum& a, const FourVector& b, bool sign=false) {
1219  return deltaPhi(a.vector3(), b.vector3(), sign);
1220  }
1221 
1223  inline double deltaPhi(const FourVector& a, const Vector3& b, bool sign=false) {
1224  return deltaPhi(a.vector3(), b, sign);
1225  }
1226 
1228  inline double deltaPhi(const Vector3& a, const FourVector& b, bool sign=false) {
1229  return deltaPhi(a, b.vector3(), sign);
1230  }
1231 
1233  inline double deltaPhi(const FourMomentum& a, const Vector3& b, bool sign=false) {
1234  return deltaPhi(a.vector3(), b, sign);
1235  }
1236 
1238  inline double deltaPhi(const Vector3& a, const FourMomentum& b, bool sign=false) {
1239  return deltaPhi(a, b.vector3(), sign);
1240  }
1241 
1243 
1244 
1246 
1247 
1249 
1250 
1252  inline double deltaEta(const FourMomentum& a, const FourMomentum& b) {
1253  return deltaEta(a.vector3(), b.vector3());
1254  }
1255 
1257  inline double deltaEta(const FourMomentum& v, double eta2) {
1258  return deltaEta(v.vector3(), eta2);
1259  }
1260 
1262  inline double deltaEta(double eta1, const FourMomentum& v) {
1263  return deltaEta(eta1, v.vector3());
1264  }
1265 
1267  inline double deltaEta(const FourVector& a, const FourVector& b) {
1268  return deltaEta(a.vector3(), b.vector3());
1269  }
1270 
1272  inline double deltaEta(const FourVector& v, double eta2) {
1273  return deltaEta(v.vector3(), eta2);
1274  }
1275 
1277  inline double deltaEta(double eta1, const FourVector& v) {
1278  return deltaEta(eta1, v.vector3());
1279  }
1280 
1282  inline double deltaEta(const FourVector& a, const FourMomentum& b) {
1283  return deltaEta(a.vector3(), b.vector3());
1284  }
1285 
1287  inline double deltaEta(const FourMomentum& a, const FourVector& b) {
1288  return deltaEta(a.vector3(), b.vector3());
1289  }
1290 
1292  inline double deltaEta(const FourVector& a, const Vector3& b) {
1293  return deltaEta(a.vector3(), b);
1294  }
1295 
1297  inline double deltaEta(const Vector3& a, const FourVector& b) {
1298  return deltaEta(a, b.vector3());
1299  }
1300 
1302  inline double deltaEta(const FourMomentum& a, const Vector3& b) {
1303  return deltaEta(a.vector3(), b);
1304  }
1305 
1307  inline double deltaEta(const Vector3& a, const FourMomentum& b) {
1308  return deltaEta(a, b.vector3());
1309  }
1310 
1312 
1313 
1315 
1316 
1318  inline double deltaRap(const FourMomentum& a, const FourMomentum& b) {
1319  return deltaRap(a.rapidity(), b.rapidity());
1320  }
1321 
1323  inline double deltaRap(const FourMomentum& v, double y2) {
1324  return deltaRap(v.rapidity(), y2);
1325  }
1326 
1328  inline double deltaRap(double y1, const FourMomentum& v) {
1329  return deltaRap(y1, v.rapidity());
1330  }
1331 
1333 
1334 
1336 
1337 
1339 
1340 
1342  inline bool cmpMomByPt(const FourMomentum& a, const FourMomentum& b) {
1343  return a.pt() > b.pt();
1344  }
1346  inline bool cmpMomByAscPt(const FourMomentum& a, const FourMomentum& b) {
1347  return a.pt() < b.pt();
1348  }
1349 
1351  inline bool cmpMomByP(const FourMomentum& a, const FourMomentum& b) {
1352  return a.vector3().mod() > b.vector3().mod();
1353  }
1355  inline bool cmpMomByAscP(const FourMomentum& a, const FourMomentum& b) {
1356  return a.vector3().mod() < b.vector3().mod();
1357  }
1358 
1360  inline bool cmpMomByEt(const FourMomentum& a, const FourMomentum& b) {
1361  return a.Et() > b.Et();
1362  }
1364  inline bool cmpMomByAscEt(const FourMomentum& a, const FourMomentum& b) {
1365  return a.Et() < b.Et();
1366  }
1367 
1369  inline bool cmpMomByE(const FourMomentum& a, const FourMomentum& b) {
1370  return a.E() > b.E();
1371  }
1373  inline bool cmpMomByAscE(const FourMomentum& a, const FourMomentum& b) {
1374  return a.E() < b.E();
1375  }
1376 
1378  inline bool cmpMomByMass(const FourMomentum& a, const FourMomentum& b) {
1379  return a.mass() > b.mass();
1380  }
1382  inline bool cmpMomByAscMass(const FourMomentum& a, const FourMomentum& b) {
1383  return a.mass() < b.mass();
1384  }
1385 
1387  inline bool cmpMomByEta(const FourMomentum& a, const FourMomentum& b) {
1388  return a.eta() < b.eta();
1389  }
1390 
1392  inline bool cmpMomByDescEta(const FourMomentum& a, const FourMomentum& b) {
1393  return a.pseudorapidity() > b.pseudorapidity();
1394  }
1395 
1397  inline bool cmpMomByAbsEta(const FourMomentum& a, const FourMomentum& b) {
1398  return fabs(a.eta()) < fabs(b.eta());
1399  }
1400 
1402  inline bool cmpMomByDescAbsEta(const FourMomentum& a, const FourMomentum& b) {
1403  return fabs(a.eta()) > fabs(b.eta());
1404  }
1405 
1407  inline bool cmpMomByRap(const FourMomentum& a, const FourMomentum& b) {
1408  return a.rapidity() < b.rapidity();
1409  }
1410 
1412  inline bool cmpMomByDescRap(const FourMomentum& a, const FourMomentum& b) {
1413  return a.rapidity() > b.rapidity();
1414  }
1415 
1417  inline bool cmpMomByAbsRap(const FourMomentum& a, const FourMomentum& b) {
1418  return fabs(a.rapidity()) < fabs(b.rapidity());
1419  }
1420 
1422  inline bool cmpMomByDescAbsRap(const FourMomentum& a, const FourMomentum& b) {
1423  return fabs(a.rapidity()) > fabs(b.rapidity());
1424  }
1425 
1427 
1428 
1430  template<typename MOMS, typename CMP>
1431  inline MOMS& sortBy(MOMS& pbs, const CMP& cmp) {
1432  std::sort(pbs.begin(), pbs.end(), cmp);
1433  return pbs;
1434  }
1436  template<typename MOMS, typename CMP>
1437  inline MOMS sortBy(const MOMS& pbs, const CMP& cmp) {
1438  MOMS rtn = pbs;
1439  std::sort(rtn.begin(), rtn.end(), cmp);
1440  return rtn;
1441  }
1442 
1444  template<typename MOMS>
1445  inline MOMS& sortByPt(MOMS& pbs) {
1446  return sortBy(pbs, cmpMomByPt);
1447  }
1449  template<typename MOMS>
1450  inline MOMS sortByPt(const MOMS& pbs) {
1451  return sortBy(pbs, cmpMomByPt);
1452  }
1453 
1455  template<typename MOMS>
1456  inline MOMS& sortByE(MOMS& pbs) {
1457  return sortBy(pbs, cmpMomByE);
1458  }
1460  template<typename MOMS>
1461  inline MOMS sortByE(const MOMS& pbs) {
1462  return sortBy(pbs, cmpMomByE);
1463  }
1464 
1466  template<typename MOMS>
1467  inline MOMS& sortByEt(MOMS& pbs) {
1468  return sortBy(pbs, cmpMomByEt);
1469  }
1471  template<typename MOMS>
1472  inline MOMS sortByEt(const MOMS& pbs) {
1473  return sortBy(pbs, cmpMomByEt);
1474  }
1475 
1477 
1478 
1480 
1481 
1483  inline double mT(const Vector3& vis, const Vector3& invis) {
1484  // return sqrt(2*vis.perp()*invis.perp() * (1 - cos(deltaPhi(vis, invis))) );
1485  return mT(vis.perp(), invis.perp(), deltaPhi(vis, invis));
1486  }
1487 
1489  inline double mT(const FourMomentum& vis, const FourMomentum& invis) {
1490  return mT(vis.p3(), invis.p3());
1491  }
1492 
1494  inline double mT(const FourMomentum& vis, const Vector3& invis) {
1495  return mT(vis.p3(), invis);
1496  }
1497 
1499  inline double mT(const Vector3& vis, const FourMomentum& invis) {
1500  return mT(vis, invis.p3());
1501  }
1502 
1504 
1505 
1507 
1508 
1510 
1511 
1513  inline std::string toString(const FourVector& lv) {
1514  ostringstream out;
1515  out << "(" << (fabs(lv.t()) < 1E-30 ? 0.0 : lv.t())
1516  << "; " << (fabs(lv.x()) < 1E-30 ? 0.0 : lv.x())
1517  << ", " << (fabs(lv.y()) < 1E-30 ? 0.0 : lv.y())
1518  << ", " << (fabs(lv.z()) < 1E-30 ? 0.0 : lv.z())
1519  << ")";
1520  return out.str();
1521  }
1522 
1524  inline std::ostream& operator<<(std::ostream& out, const FourVector& lv) {
1525  out << toString(lv);
1526  return out;
1527  }
1528 
1530 
1531 
1534 
1535  //typedef FourVector V4; //< generic
1536  typedef FourVector X4; //< spatial
1537  typedef FourMomentum P4; //< momentum
1539 
1541 
1542  typedef std::vector<FourVector> FourVectors;
1543  typedef std::vector<FourMomentum> FourMomenta;
1544  typedef std::vector<X4> X4s;
1545  typedef std::vector<P4> P4a;
1547 
1548 
1549 }
1550 
1551 #endif
FourMomentum & setPy(double py)
Set y-component of momentum .
Definition: Vector4.hh:349
FourVector operator-() const
Multiply all components (space and time) by -1.
Definition: Vector4.hh:213
Definition: ALICE_2010_I880049.cc:13
double abseta() const
Get the directly (alias).
Definition: Vector4.hh:159
double py2() const
Get y-squared .
Definition: Vector4.hh:550
MOMS & sortByPt(MOMS &pbs)
Sort a container of momenta by pT (decreasing) and return by reference for non-const inputs...
Definition: Vector4.hh:1445
FourMomentum & setXYZM(double px, double py, double pz, double mass)
Alias for setPM.
Definition: Vector4.hh:391
Vector3 rhoVec() const
Synonym for polarVec.
Definition: Vector3.hh:121
double pT() const
Calculate the transverse momentum .
Definition: Vector4.hh:628
FourMomentum & operator+=(const FourMomentum &v)
Add to this 4-vector. NB time as well as space components are added.
Definition: Vector4.hh:736
Vector3 betaVec() const
Definition: Vector4.hh:671
double perp() const
Synonym for polarRadius.
Definition: Vector4.hh:108
double rapidity(double E, double pz)
Calculate a rapidity value from the supplied energy E and longitudinal momentum pz.
Definition: MathUtils.hh:602
double px() const
Get x-component of momentum .
Definition: Vector4.hh:543
double mod2() const
Calculate the modulus-squared of a vector. .
Definition: VectorN.hh:84
double eta() const
Synonym for pseudorapidity.
Definition: Vector4.hh:152
static FourMomentum mkRapPhiMPt(double y, double phi, double mass, double pt)
Make a vector from (y,phi,pT) coordinates and the mass.
Definition: Vector4.hh:796
FourMomentum & setPz(double pz)
Set z-component of momentum .
Definition: Vector4.hh:355
FourVector & operator*=(double a)
Multiply by a scalar.
Definition: Vector4.hh:189
double polarRadius() const
Polar radius.
Definition: Vector3.hh:139
double absrap() const
Absolute rapidity.
Definition: Vector4.hh:605
bool cmpMomByEt(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by decreasing transverse energy.
Definition: Vector4.hh:1360
double pt2() const
Calculate the squared transverse momentum .
Definition: Vector4.hh:623
Vector3 polarVec() const
Projection of 3-vector on to the plane.
Definition: Vector4.hh:117
double azimuthalAngle(const PhiMapping mapping=ZERO_2PI) const
Angle subtended by the 3-vector&#39;s projection in x-y and the x-axis.
Definition: Vector4.hh:130
bool cmpMomByMass(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by decreasing mass.
Definition: Vector4.hh:1378
double Et() const
Calculate the transverse energy .
Definition: Vector4.hh:641
Vector3 ptvec() const
Synonym for pTvec.
Definition: Vector4.hh:614
double rho() const
Synonym for polarRadius.
Definition: Vector3.hh:147
double deltaRap(double y1, double y2)
Definition: MathUtils.hh:584
Vector3 gammaVec() const
Definition: Vector4.hh:659
double absrapidity() const
Absolute rapidity.
Definition: Vector4.hh:601
double p2() const
Get the modulus-squared of the 3-momentum.
Definition: Vector4.hh:586
FourMomentum & operator*=(double a)
Multiply by a scalar.
Definition: Vector4.hh:724
double perp() const
Synonym for polarRadius.
Definition: Vector3.hh:143
static FourMomentum mkThetaPhiME(double theta, double phi, double mass, double E)
Make a vector from (theta,phi,energy) coordinates and the mass.
Definition: Vector4.hh:801
FourMomentum & setE(double E)
Set energy (time component of momentum).
Definition: Vector4.hh:337
double contract(const FourVector &v) const
Contract two 4-vectors, with metric signature (+ - - -).
Definition: Vector4.hh:173
static FourMomentum mkEtaPhiMPt(double eta, double phi, double mass, double pt)
Make a vector from (eta,phi,pT) coordinates and the mass.
Definition: Vector4.hh:786
static FourMomentum mkEtaPhiME(double eta, double phi, double mass, double E)
Make a vector from (eta,phi,energy) coordinates and the mass.
Definition: Vector4.hh:781
double mass2() const
Get the squared mass (the Lorentz self-invariant).
Definition: Vector4.hh:572
bool cmpMomByDescRap(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by decreasing rapidity.
Definition: Vector4.hh:1412
double abspseudorapidity() const
Get the directly.
Definition: Vector4.hh:157
double rapidity() const
Calculate the rapidity.
Definition: Vector4.hh:592
double pz2() const
Get z-squared .
Definition: Vector4.hh:555
bool cmpMomByAscE(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by increasing energy.
Definition: Vector4.hh:1373
FourMomentum & operator-=(const FourMomentum &v)
Subtract from this 4-vector. NB time as well as space components are subtracted.
Definition: Vector4.hh:742
Vector3 p3() const
Get 3-momentum part, .
Definition: Vector4.hh:578
double pseudorapidity() const
Pseudorapidity (defined purely by the 3-vector components)
Definition: Vector4.hh:148
bool cmpMomByAscEt(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by increasing transverse energy.
Definition: Vector4.hh:1364
double polarRadius() const
Magnitude of projection of 3-vector on to the plane.
Definition: Vector4.hh:104
double E2() const
Get energy-squared .
Definition: Vector4.hh:540
bool cmpMomByEta(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by increasing eta (pseudorapidity)
Definition: Vector4.hh:1387
FourMomentum & setThetaPhiME(double theta, double phi, double mass, double E)
Definition: Vector4.hh:475
double angle(const Vector3 &v) const
Angle in radians to another vector.
Definition: Vector3.hh:89
FourMomentum & setEtaPhiMPt(double eta, double phi, double mass, double pt)
Definition: Vector4.hh:416
bool cmpMomByE(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by decreasing energy.
Definition: Vector4.hh:1369
Vector3 polarVec() const
Polar projection of this vector into the x-y plane.
Definition: Vector3.hh:111
FourMomentum & setPtPhiME(double pt, double phi, double mass, double E)
Definition: Vector4.hh:517
FourMomentum & setThetaPhiMPt(double theta, double phi, double mass, double pt)
Definition: Vector4.hh:498
static FourMomentum mkPtPhiME(double pt, double phi, double mass, double E)
Make a vector from (pT,phi,energy) coordinates and the mass.
Definition: Vector4.hh:811
FourMomentum & setPx(double px)
Set x-component of momentum .
Definition: Vector4.hh:343
Struct for sorting by increasing energy.
Definition: Vector4.hh:690
double polarAngle() const
Angle subtended by the vector and the z-axis.
Definition: Vector3.hh:166
PhiMapping
Enum for range of to be mapped into.
Definition: MathHeader.hh:31
double gamma() const
Definition: Vector4.hh:653
double dot(const FourVector &v) const
Contract two 4-vectors, with metric signature (+ - - -).
Definition: Vector4.hh:179
double pseudorapidity() const
Definition: Vector3.hh:182
double deltaPhi(double phi1, double phi2, bool sign=false)
Calculate the difference between two angles in radians.
Definition: MathUtils.hh:569
bool cmpMomByAscMass(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by increasing mass.
Definition: Vector4.hh:1382
double angle(const Vector3 &v3) const
Angle between this vector and another (3-vector)
Definition: Vector4.hh:84
bool cmpMomByAscPt(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by increasing pT.
Definition: Vector4.hh:1346
double rho2() const
Synonym for polarRadius2.
Definition: Vector3.hh:134
double mass() const
Get the mass (the Lorentz self-invariant).
Definition: Vector4.hh:561
double beta() const
Definition: Vector4.hh:665
double rho2() const
Synonym for polarRadius2.
Definition: Vector4.hh:99
Object implementing Lorentz transform calculations and boosts.
Definition: LorentzTrans.hh:21
FourMomentum operator-() const
Multiply all components (time and space) by -1.
Definition: Vector4.hh:748
bool cmpMomByDescEta(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by decreasing eta (pseudorapidity)
Definition: Vector4.hh:1392
bool cmpMomByP(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by decreasing 3-momentum magnitude |p|.
Definition: Vector4.hh:1351
Specialisation of VectorN to a general (non-momentum) Lorentz 4-vector.
Definition: Vector4.hh:22
double mT(double pT1, double pT2, double dphi)
Definition: MathUtils.hh:619
double p() const
Get the modulus of the 3-momentum.
Definition: Vector4.hh:581
double invariant(const FourVector &lv)
Definition: Vector4.hh:277
static FourMomentum mkXYZE(double px, double py, double pz, double E)
Make a vector from (px,py,pz,E) coordinates.
Definition: Vector4.hh:771
double deltaEta(double eta1, double eta2)
Definition: MathUtils.hh:577
A minimal base class for -dimensional vectors.
Definition: VectorN.hh:13
double deltaR2(double rap1, double phi1, double rap2, double phi2)
Definition: MathUtils.hh:590
double polarRadius2() const
Mod-square of the projection of the 3-vector on to the plane This is a more efficient function than ...
Definition: Vector4.hh:91
double deltaR(double rap1, double phi1, double rap2, double phi2)
Definition: MathUtils.hh:597
double E() const
Get energy (time component of momentum).
Definition: Vector4.hh:538
bool cmpMomByPt(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by decreasing pT.
Definition: Vector4.hh:1342
Vector3 perpVec() const
Synonym for polarVec.
Definition: Vector4.hh:121
MOMS & sortByEt(MOMS &pbs)
Sort a container of momenta by Et (decreasing) and return by reference for non-const inputs...
Definition: Vector4.hh:1467
RapScheme
Enum for rapidity variable to be used in calculating , applying rapidity cuts, etc.
Definition: MathHeader.hh:28
bool cmpMomByAscP(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by increasing 3-momentum magnitude |p|.
Definition: Vector4.hh:1355
FourMomentum & operator/=(double a)
Divide by a scalar.
Definition: Vector4.hh:730
Vector3 perpVec() const
Synonym for polarVec.
Definition: Vector3.hh:117
double operator*(const FourVector &v) const
Contract two 4-vectors, with metric signature (+ - - -).
Definition: Vector4.hh:184
bool cmpMomByRap(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by increasing rapidity.
Definition: Vector4.hh:1407
FourMomentum & setRapPhiMPt(double y, double phi, double mass, double pt)
Definition: Vector4.hh:458
MOMS & sortBy(MOMS &pbs, const CMP &cmp)
Sort a container of momenta by cmp and return by reference for non-const inputs.
Definition: Vector4.hh:1431
double perp2() const
Synonym for polarRadius2.
Definition: Vector4.hh:95
double pt() const
Calculate the transverse momentum .
Definition: Vector4.hh:632
Vector3 vector3() const
Get the spatial part of the 4-vector as a 3-vector.
Definition: Vector4.hh:162
Vector3 unit() const
Synonym for unitVec.
Definition: Vector3.hh:105
FourMomentum & setXYZE(double px, double py, double pz, double E)
Alias for setPE.
Definition: Vector4.hh:369
string to_str(const T &x)
Convert any object to a string.
Definition: Utils.hh:82
static FourMomentum mkThetaPhiMPt(double theta, double phi, double mass, double pt)
Make a vector from (theta,phi,pT) coordinates and the mass.
Definition: Vector4.hh:806
FourMomentum & setRapPhiME(double y, double phi, double mass, double E)
Definition: Vector4.hh:438
double angle(const FourVector &v) const
Angle between this vector and another.
Definition: Vector4.hh:80
double theta() const
Synonym for polarAngle.
Definition: Vector4.hh:143
std::enable_if< std::is_floating_point< NUM >::value, bool >::type isZero(NUM val, double tolerance=1e-8)
Compare a number to zero.
Definition: MathUtils.hh:21
double polarAngle() const
Angle subtended by the 3-vector and the z-axis.
Definition: Vector4.hh:139
Vector3 pTvec() const
Calculate the transverse momentum vector .
Definition: Vector4.hh:610
FourMomentum & setEtaPhiME(double eta, double phi, double mass, double E)
Definition: Vector4.hh:400
double pz() const
Get z-component of momentum .
Definition: Vector4.hh:553
double phi(const PhiMapping mapping=ZERO_2PI) const
Synonym for azimuthalAngle.
Definition: Vector4.hh:134
double py() const
Get y-component of momentum .
Definition: Vector4.hh:548
double azimuthalAngle(const PhiMapping mapping=ZERO_2PI) const
Angle subtended by the vector&#39;s projection in x-y and the x-axis.
Definition: Vector3.hh:152
std::string toString(const AnalysisInfo &ai)
String representation.
Definition: AnalysisInfo.cc:144
FourVector & operator-=(const FourVector &v)
Subtract from this 4-vector. NB time as well as space components are subtracted.
Definition: Vector4.hh:207
double px2() const
Get x-squared .
Definition: Vector4.hh:545
FourMomentum & setPE(double px, double py, double pz, double E)
Set the p coordinates and energy simultaneously.
Definition: Vector4.hh:362
bool cmpMomByAbsRap(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by increasing absolute rapidity.
Definition: Vector4.hh:1417
double theta() const
Synonym for polarAngle.
Definition: Vector3.hh:173
double phi(const PhiMapping mapping=ZERO_2PI) const
Synonym for azimuthalAngle.
Definition: Vector3.hh:161
Three-dimensional specialisation of Vector.
Definition: Vector3.hh:26
MOMS & sortByE(MOMS &pbs)
Sort a container of momenta by E (decreasing) and return by reference for non-const inputs...
Definition: Vector4.hh:1456
FourVector & operator+=(const FourVector &v)
Add to this 4-vector.
Definition: Vector4.hh:201
bool cmpMomByDescAbsEta(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by increasing absolute eta (pseudorapidity)
Definition: Vector4.hh:1402
Struct for sorting by decreasing energy.
Definition: Vector4.hh:704
Vector3 rhoVec() const
Synonym for polarVec.
Definition: Vector4.hh:125
FourVector & operator/=(double a)
Divide by a scalar.
Definition: Vector4.hh:195
Specialized version of the FourVector with momentum/energy functionality.
Definition: Vector4.hh:301
double eta() const
Synonym for pseudorapidity.
Definition: Vector3.hh:192
double Et2() const
Calculate the transverse energy .
Definition: Vector4.hh:637
static FourMomentum mkRapPhiME(double y, double phi, double mass, double E)
Make a vector from (y,phi,energy) coordinates and the mass.
Definition: Vector4.hh:791
static FourMomentum mkXYZM(double px, double py, double pz, double mass)
Make a vector from (px,py,pz) coordinates and the mass.
Definition: Vector4.hh:776
FourMomentum reverse() const
Multiply space components only by -1.
Definition: Vector4.hh:755
bool cmpMomByDescAbsRap(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by decreasing absolute rapidity.
Definition: Vector4.hh:1422
FourMomentum & setPM(double px, double py, double pz, double mass)
Set the p coordinates and mass simultaneously.
Definition: Vector4.hh:383
std::enable_if< std::is_arithmetic< NUM >::value, NUM >::type sqr(NUM a)
Named number-type squaring operation.
Definition: MathUtils.hh:189
double rap() const
Alias for rapidity.
Definition: Vector4.hh:596
double polarRadius2() const
Square of the polar radius (.
Definition: Vector3.hh:126
double mod() const
Calculate the modulus of a vector. .
Definition: VectorN.hh:95
double rho() const
Synonym for polarRadius.
Definition: Vector4.hh:112
double perp2() const
Synonym for polarRadius2.
Definition: Vector3.hh:130
bool cmpMomByAbsEta(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by increasing absolute eta (pseudorapidity)
Definition: Vector4.hh:1397
FourVector reverse() const
Multiply space components only by -1.
Definition: Vector4.hh:220
std::enable_if< std::is_arithmetic< NUM >::value, int >::type sign(NUM val)
Find the sign of a number.
Definition: MathUtils.hh:236
double pT2() const
Calculate the squared transverse momentum .
Definition: Vector4.hh:619
Vector3 boostVector() const
Deprecated alias for betaVec.
Definition: Vector4.hh:678
Cmp< T > cmp(const T &t1, const T &t2)
Global helper function for easy creation of Cmp objects.
Definition: Cmp.hh:285