rivet is hosted by Hepforge, IPPP Durham
Rivet 4.0.2
Vector4.hh
1#ifndef RIVET_MATH_VECTOR4
2#define RIVET_MATH_VECTOR4
3
4#include "Rivet/Tools/TypeTraits.hh"
5#include "Rivet/Math/MathConstants.hh"
6#include "Rivet/Math/MathUtils.hh"
7#include "Rivet/Math/VectorN.hh"
8#include "Rivet/Math/Vector3.hh"
9
10// Forward declaration
11namespace fastjet { class PseudoJet; }
12
13namespace Rivet {
14
15
16 class FourVector;
17 typedef FourVector Vector4;
18 typedef FourVector V4;
19
20 class FourMomentum;
21 typedef FourMomentum P4;
22
23 class LorentzTransform;
24 FourVector transform(const LorentzTransform& lt, const FourVector& v4);
25
26
30 class FourVector : public Vector<4> {
31 friend FourVector multiply(const double a, const FourVector& v);
32 friend FourVector multiply(const FourVector& v, const double a);
33 friend FourVector add(const FourVector& a, const FourVector& b);
34 friend FourVector transform(const LorentzTransform& lt, const FourVector& v4);
35
36 public:
37
38 FourVector() : Vector<4>() { }
39
41 FourVector(const V4TYPE& other) {
42 this->setT(other.t());
43 this->setX(other.x());
44 this->setY(other.y());
45 this->setZ(other.z());
46 }
47
49 : Vector<4>(other) { }
50
51 FourVector(const double t, const double x, const double y, const double z) {
52 this->setT(t);
53 this->setX(x);
54 this->setY(y);
55 this->setZ(z);
56 }
57
58 virtual ~FourVector() { }
59
64 operator fastjet::PseudoJet () const;
65
66
67 public:
68
69 double t() const { return get(0); }
70 double t2() const { return sqr(t()); }
71 FourVector& setT(const double t) { set(0, t); return *this; }
72
73 double x() const { return get(1); }
74 double x2() const { return sqr(x()); }
75 FourVector& setX(const double x) { set(1, x); return *this; }
76
77 double y() const { return get(2); }
78 double y2() const { return sqr(y()); }
79 FourVector& setY(const double y) { set(2, y); return *this; }
80
81 double z() const { return get(3); }
82 double z2() const { return sqr(z()); }
83 FourVector& setZ(const double z) { set(3, z); return *this; }
84
85 double invariant() const {
86 // Done this way for numerical precision
87 return (t() + z())*(t() - z()) - x()*x() - y()*y();
88 }
89
90 bool isNull() const {
91 return Rivet::isZero(invariant());
92 }
93
95 double angle(const FourVector& v) const {
96 return vector3().angle( v.vector3() );
97 }
99 double angle(const Vector3& v3) const {
100 return vector3().angle(v3);
101 }
102
106 double polarRadius2() const {
107 return vector3().polarRadius2();
108 }
110 double perp2() const {
111 return vector3().perp2();
112 }
114 double rho2() const {
115 return vector3().rho2();
116 }
117
119 double polarRadius() const {
120 return vector3().polarRadius();
121 }
123 double perp() const {
124 return vector3().perp();
125 }
127 double rho() const {
128 return vector3().rho();
129 }
130
133 return vector3().polarVec();
134 }
136 Vector3 perpVec() const {
137 return vector3().perpVec();
138 }
140 Vector3 rhoVec() const {
141 return vector3().rhoVec();
142 }
143
145 double azimuthalAngle(const PhiMapping mapping=ZERO_2PI) const {
147 }
149 double phi(const PhiMapping mapping=ZERO_2PI) const {
150 return vector3().phi(mapping);
151 }
152
154 double polarAngle() const {
155 return vector3().polarAngle();
156 }
158 double theta() const {
159 return vector3().theta();
160 }
161
163 double pseudorapidity() const {
164 return vector3().pseudorapidity();
165 }
167 double eta() const {
168 return vector3().eta();
169 }
170
172 double abspseudorapidity() const { return fabs(eta()); }
174 double abseta() const { return fabs(eta()); }
175
177 Vector3 vector3() const {
178 return Vector3(get(1), get(2), get(3));
179 }
180
182 operator Vector3 () const { return vector3(); }
183
184
185 public:
186
188 double contract(const FourVector& v) const {
189 const double result = t()*v.t() - x()*v.x() - y()*v.y() - z()*v.z();
190 return result;
191 }
192
194 double dot(const FourVector& v) const {
195 return contract(v);
196 }
197
199 double operator * (const FourVector& v) const {
200 return contract(v);
201 }
202
205 _vec = multiply(a, *this)._vec;
206 return *this;
207 }
208
211 _vec = multiply(1.0/a, *this)._vec;
212 return *this;
213 }
214
217 _vec = add(*this, v)._vec;
218 return *this;
219 }
220
223 _vec = add(*this, -v)._vec;
224 return *this;
225 }
226
230 result._vec = -_vec;
231 return result;
232 }
233
236 FourVector result = -*this;
237 result.setT(-result.t());
238 return result;
239 }
240
241 };
242
243
245 inline double contract(const FourVector& a, const FourVector& b) {
246 return a.contract(b);
247 }
248
250 inline double dot(const FourVector& a, const FourVector& b) {
251 return contract(a, b);
252 }
253
254 inline FourVector multiply(const double a, const FourVector& v) {
255 FourVector result;
256 result._vec = a * v._vec;
257 return result;
258 }
259
260 inline FourVector multiply(const FourVector& v, const double a) {
261 return multiply(a, v);
262 }
263
264 inline FourVector operator * (const double a, const FourVector& v) {
265 return multiply(a, v);
266 }
267
268 inline FourVector operator * (const FourVector& v, const double a) {
269 return multiply(a, v);
270 }
271
272 inline FourVector operator / (const FourVector& v, const double a) {
273 return multiply(1.0/a, v);
274 }
275
276 inline FourVector add(const FourVector& a, const FourVector& b) {
277 FourVector result;
278 result._vec = a._vec + b._vec;
279 return result;
280 }
281
282 inline FourVector operator+(const FourVector& a, const FourVector& b) {
283 return add(a, b);
284 }
285
286 inline FourVector operator-(const FourVector& a, const FourVector& b) {
287 return add(a, -b);
288 }
289
292 inline double invariant(const FourVector& lv) {
293 return lv.invariant();
294 }
295
297 inline double angle(const FourVector& a, const FourVector& b) {
298 return a.angle(b);
299 }
300
302 inline double angle(const Vector3& a, const FourVector& b) {
303 return angle( a, b.vector3() );
304 }
305
307 inline double angle(const FourVector& a, const Vector3& b) {
308 return a.angle(b);
309 }
310
311
313
314
316 class FourMomentum : public FourVector {
317 friend FourMomentum multiply(const double a, const FourMomentum& v);
318 friend FourMomentum multiply(const FourMomentum& v, const double a);
319 friend FourMomentum add(const FourMomentum& a, const FourMomentum& b);
320 friend FourMomentum transform(const LorentzTransform& lt, const FourMomentum& v4);
321
322 public:
323 FourMomentum() { }
324
326 FourMomentum(const V4TYPE& other) {
327 this->setE(other.t());
328 this->setPx(other.x());
329 this->setPy(other.y());
330 this->setPz(other.z());
331 }
332
334 : FourVector(other) { }
335
336 FourMomentum(const double E, const double px, const double py, const double pz) {
337 this->setE(E);
338 this->setPx(px);
339 this->setPy(py);
340 this->setPz(pz);
341 }
342
343 ~FourMomentum() {}
344
345 public:
346
347
350
352 FourMomentum& setE(double E) {
353 setT(E);
354 return *this;
355 }
356
359 setX(px);
360 return *this;
361 }
362
365 setY(py);
366 return *this;
367 }
368
371 setZ(pz);
372 return *this;
373 }
374
375
377 FourMomentum& setPE(double px, double py, double pz, double E) {
378 if (E < 0)
379 throw std::invalid_argument("Negative energy given as argument: " + to_str(E));
380 setPx(px); setPy(py); setPz(pz); setE(E);
381 return *this;
382 }
384 FourMomentum& setXYZE(double px, double py, double pz, double E) {
385 return setPE(px, py, pz, E);
386 }
387 // /// Near-alias with switched arg order
388 // FourMomentum& setEP(double E, double px, double py, double pz) {
389 // return setPE(px, py, pz, E);
390 // }
391 // /// Alias for setEP
392 // FourMomentum& setEXYZ(double E, double px, double py, double pz) {
393 // return setEP(E, px, py, pz);
394 // }
395
396
398 FourMomentum& setPM(double px, double py, double pz, double mass) {
399 if (mass < 0)
400 throw std::invalid_argument("Negative mass given as argument: " + to_str(mass));
401 const double E = sqrt( sqr(mass) + sqr(px) + sqr(py) + sqr(pz) );
402 // setPx(px); setPy(py); setPz(pz); setE(E);
403 return setPE(px, py, pz, E);
404 }
406 FourMomentum& setXYZM(double px, double py, double pz, double mass) {
407 return setPM(px, py, pz, mass);
408 }
409
410
415 FourMomentum& setEtaPhiME(double eta, double phi, double mass, double E) {
416 if (mass < 0)
417 throw std::invalid_argument("Negative mass given as argument");
418 if (E < 0)
419 throw std::invalid_argument("Negative energy given as argument");
420 const double theta = 2 * atan(exp(-eta));
422 throw std::domain_error("Polar angle outside 0..pi in calculation");
424 return *this;
425 }
426
431 FourMomentum& setEtaPhiMPt(double eta, double phi, double mass, double pt) {
432 if (mass < 0)
433 throw std::invalid_argument("Negative mass given as argument");
434 if (pt < 0)
435 throw std::invalid_argument("Negative transverse momentum given as argument");
436 const double theta = 2 * atan(exp(-eta));
438 throw std::domain_error("Polar angle outside 0..pi in calculation");
439 const double p = pt / sin(theta);
440 const double E = sqrt( sqr(p) + sqr(mass) );
442 return *this;
443 }
444
453 FourMomentum& setRapPhiME(double y, double phi, double mass, double E) {
454 if (mass < 0)
455 throw std::invalid_argument("Negative mass given as argument");
456 if (E < 0)
457 throw std::invalid_argument("Negative energy given as argument");
458 const double sqrt_pt2_m2 = E / cosh(y);
459 const double pt = sqrt( sqr(sqrt_pt2_m2) - sqr(mass) );
460 if (pt < 0)
461 throw std::domain_error("Negative transverse momentum in calculation");
462 const double pz = sqrt_pt2_m2 * sinh(y);
463 const double px = pt * cos(phi);
464 const double py = pt * sin(phi);
465 setPE(px, py, pz, E);
466 return *this;
467 }
468
473 FourMomentum& setRapPhiMPt(double y, double phi, double mass, double pt) {
474 if (mass < 0)
475 throw std::invalid_argument("Negative mass given as argument");
476 if (pt < 0)
477 throw std::invalid_argument("Negative transverse mass given as argument");
478 const double E = sqrt( sqr(pt) + sqr(mass) ) * cosh(y);
479 if (E < 0)
480 throw std::domain_error("Negative energy in calculation");
481 setRapPhiME(y, phi, mass, E);
482 return *this;
483 }
484
490 FourMomentum& setThetaPhiME(double theta, double phi, double mass, double E) {
492 throw std::invalid_argument("Polar angle outside 0..pi given as argument");
493 if (mass < 0)
494 throw std::invalid_argument("Negative mass given as argument");
495 if (E < 0)
496 throw std::invalid_argument("Negative energy given as argument");
497 const double p = sqrt( sqr(E) - sqr(mass) );
498 const double pz = p * cos(theta);
499 const double pt = p * sin(theta);
500 if (pt < 0)
501 throw std::invalid_argument("Negative transverse momentum in calculation");
502 const double px = pt * cos(phi);
503 const double py = pt * sin(phi);
504 setPE(px, py, pz, E);
505 return *this;
506 }
507
513 FourMomentum& setThetaPhiMPt(double theta, double phi, double mass, double pt) {
515 throw std::invalid_argument("Polar angle outside 0..pi given as argument");
516 if (mass < 0)
517 throw std::invalid_argument("Negative mass given as argument");
518 if (pt < 0)
519 throw std::invalid_argument("Negative transverse momentum given as argument");
520 const double p = pt / sin(theta);
521 const double px = pt * cos(phi);
522 const double py = pt * sin(phi);
523 const double pz = p * cos(theta);
524 const double E = sqrt( sqr(p) + sqr(mass) );
525 setPE(px, py, pz, E);
526 return *this;
527 }
528
532 FourMomentum& setPtPhiME(double pt, double phi, double mass, double E) {
533 if (pt < 0)
534 throw std::invalid_argument("Negative transverse momentum given as argument");
535 if (mass < 0)
536 throw std::invalid_argument("Negative mass given as argument");
537 if (E < 0)
538 throw std::invalid_argument("Negative energy given as argument");
539 const double px = pt * cos(phi);
540 const double py = pt * sin(phi);
541 const double pz = sqrt(sqr(E) - sqr(mass) - sqr(pt));
542 setPE(px, py, pz, E);
543 return *this;
544 }
545
547
548
551
553 double E() const { return t(); }
555 double E2() const { return t2(); }
556
558 double px() const { return x(); }
560 double px2() const { return x2(); }
561
563 double py() const { return y(); }
565 double py2() const { return y2(); }
566
568 double pz() const { return z(); }
570 double pz2() const { return z2(); }
571
572
576 double mass() const {
577 // assert(Rivet::isZero(mass2()) || mass2() > 0);
578 // if (Rivet::isZero(mass2())) {
579 // return 0.0;
580 // } else {
581 // return sqrt(mass2());
582 // }
583 return sign(mass2()) * sqrt(fabs(mass2()));
584 }
585
587 double mass2() const {
588 return invariant();
589 }
590
591
593 Vector3 p3() const { return vector3(); }
594
596 double p() const {
597 return p3().mod();
598 }
599
601 double p2() const {
602 return p3().mod2();
603 }
604
605
607 double rapidity() const {
608 if (E() == 0.0) return 0.0;
609 if (E() == fabs(pz())) return std::copysign(INF, pz());
610 return 0.5 * std::log( (E() + pz()) / (E() - pz()) );
611 }
613 double rap() const {
614 return rapidity();
615 }
616
618 double absrapidity() const {
619 return fabs(rapidity());
620 }
622 double absrap() const {
623 return fabs(rap());
624 }
625
627 Vector3 pTvec() const {
628 return p3().polarVec();
629 }
631 Vector3 ptvec() const {
632 return pTvec();
633 }
634
636 double pT2() const {
637 return vector3().polarRadius2();
638 }
640 double pt2() const {
641 return vector3().polarRadius2();
642 }
643
645 double pT() const {
646 return sqrt(pT2());
647 }
649 double pt() const {
650 return sqrt(pT2());
651 }
652
654 double Et2() const {
655 return Et() * Et();
656 }
658 double Et() const {
659 return E() * sin(polarAngle());
660 }
661
663
664
667
670 double gamma() const {
671 return sqrt(E2()/mass2());
672 }
673
677 return gamma() * p3().unit();
678 }
679
682 double beta() const {
683 return p()/E();
684 }
685
688 Vector3 betaVec() const {
689 // return Vector3(px()/E(), py()/E(), pz()/E());
690 return p3()/E();
691 }
692
694
695
697
700
703 _vec = multiply(a, *this)._vec;
704 return *this;
705 }
706
709 _vec = multiply(1.0/a, *this)._vec;
710 return *this;
711 }
712
715 _vec = add(*this, v)._vec;
716 return *this;
717 }
718
721 _vec = add(*this, -v)._vec;
722 return *this;
723 }
724
728 result._vec = -_vec;
729 return result;
730 }
731
734 FourMomentum result = -*this;
735 result.setE(-result.E());
736 return result;
737 }
738
740
741
743
744
747
749 static FourMomentum mkXYZE(double px, double py, double pz, double E) {
750 return FourMomentum().setPE(px, py, pz, E);
751 }
752
754 static FourMomentum mkXYZM(double px, double py, double pz, double mass) {
755 return FourMomentum().setPM(px, py, pz, mass);
756 }
757
759 static FourMomentum mkEtaPhiME(double eta, double phi, double mass, double E) {
760 return FourMomentum().setEtaPhiME(eta, phi, mass, E);
761 }
762
764 static FourMomentum mkEtaPhiMPt(double eta, double phi, double mass, double pt) {
765 return FourMomentum().setEtaPhiMPt(eta, phi, mass, pt);
766 }
767
769 static FourMomentum mkRapPhiME(double y, double phi, double mass, double E) {
770 return FourMomentum().setRapPhiME(y, phi, mass, E);
771 }
772
774 static FourMomentum mkRapPhiMPt(double y, double phi, double mass, double pt) {
775 return FourMomentum().setRapPhiMPt(y, phi, mass, pt);
776 }
777
779 static FourMomentum mkThetaPhiME(double theta, double phi, double mass, double E) {
781 }
782
784 static FourMomentum mkThetaPhiMPt(double theta, double phi, double mass, double pt) {
786 }
787
789 static FourMomentum mkPtPhiME(double pt, double phi, double mass, double E) {
790 return FourMomentum().setPtPhiME(pt, phi, mass, E);
791 }
792
794
795
796 };
797
798
799 inline FourMomentum multiply(const double a, const FourMomentum& v) {
800 FourMomentum result;
801 result._vec = a * v._vec;
802 return result;
803 }
804
805 inline FourMomentum multiply(const FourMomentum& v, const double a) {
806 return multiply(a, v);
807 }
808
809 inline FourMomentum operator*(const double a, const FourMomentum& v) {
810 return multiply(a, v);
811 }
812
813 inline FourMomentum operator*(const FourMomentum& v, const double a) {
814 return multiply(a, v);
815 }
816
817 inline FourMomentum operator/(const FourMomentum& v, const double a) {
818 return multiply(1.0/a, v);
819 }
820
821 inline FourMomentum add(const FourMomentum& a, const FourMomentum& b) {
822 FourMomentum result;
823 result._vec = a._vec + b._vec;
824 return result;
825 }
826
827 inline FourMomentum operator+(const FourMomentum& a, const FourMomentum& b) {
828 return add(a, b);
829 }
830
831 inline FourMomentum operator-(const FourMomentum& a, const FourMomentum& b) {
832 return add(a, -b);
833 }
834
835
837
838
841
850 inline double deltaR2(const FourVector& a, const FourVector& b,
851 RapScheme scheme=PSEUDORAPIDITY) {
852 switch (scheme) {
853 case PSEUDORAPIDITY :
854 return deltaR2(a.vector3(), b.vector3());
855 case RAPIDITY:
856 {
857 const FourMomentum* ma = dynamic_cast<const FourMomentum*>(&a);
858 const FourMomentum* mb = dynamic_cast<const FourMomentum*>(&b);
859 if (!ma || !mb) {
860 string err = "deltaR with scheme RAPIDITY can only be called with FourMomentum objects, not FourVectors";
861 throw std::runtime_error(err);
862 }
863 return deltaR2(*ma, *mb, scheme);
864 }
865 default:
866 throw std::runtime_error("The specified deltaR scheme is not yet implemented");
867 }
868 }
869
878 inline double deltaR(const FourVector& a, const FourVector& b,
879 RapScheme scheme=PSEUDORAPIDITY) {
880 return sqrt(deltaR2(a, b, scheme));
881 }
882
883
884
891 inline double deltaR2(const FourVector& v,
892 double eta2, double phi2,
893 RapScheme scheme=PSEUDORAPIDITY) {
894 switch (scheme) {
895 case PSEUDORAPIDITY :
896 return deltaR2(v.vector3(), eta2, phi2);
897 case RAPIDITY:
898 {
899 const FourMomentum* mv = dynamic_cast<const FourMomentum*>(&v);
900 if (!mv) {
901 string err = "deltaR with scheme RAPIDITY can only be called with FourMomentum objects, not FourVectors";
902 throw std::runtime_error(err);
903 }
904 return deltaR2(*mv, eta2, phi2, scheme);
905 }
906 default:
907 throw std::runtime_error("The specified deltaR scheme is not yet implemented");
908 }
909 }
910
917 inline double deltaR(const FourVector& v,
918 double eta2, double phi2,
919 RapScheme scheme=PSEUDORAPIDITY) {
920 return sqrt(deltaR2(v, eta2, phi2, scheme));
921 }
922
923
930 inline double deltaR2(double eta1, double phi1,
931 const FourVector& v,
932 RapScheme scheme=PSEUDORAPIDITY) {
933 switch (scheme) {
934 case PSEUDORAPIDITY :
935 return deltaR2(eta1, phi1, v.vector3());
936 case RAPIDITY:
937 {
938 const FourMomentum* mv = dynamic_cast<const FourMomentum*>(&v);
939 if (!mv) {
940 string err = "deltaR with scheme RAPIDITY can only be called with FourMomentum objects, not FourVectors";
941 throw std::runtime_error(err);
942 }
943 return deltaR2(eta1, phi1, *mv, scheme);
944 }
945 default:
946 throw std::runtime_error("The specified deltaR scheme is not yet implemented");
947 }
948 }
949
956 inline double deltaR(double eta1, double phi1,
957 const FourVector& v,
958 RapScheme scheme=PSEUDORAPIDITY) {
959 return sqrt(deltaR2(eta1, phi1, v, scheme));
960 }
961
962
969 inline double deltaR2(const FourMomentum& a, const FourMomentum& b,
970 RapScheme scheme=PSEUDORAPIDITY) {
971 switch (scheme) {
972 case PSEUDORAPIDITY:
973 return deltaR2(a.vector3(), b.vector3());
974 case RAPIDITY:
975 return deltaR2(a.rapidity(), a.azimuthalAngle(), b.rapidity(), b.azimuthalAngle());
976 default:
977 throw std::runtime_error("The specified deltaR scheme is not yet implemented");
978 }
979 }
980
987 inline double deltaR(const FourMomentum& a, const FourMomentum& b,
988 RapScheme scheme=PSEUDORAPIDITY) {
989 return sqrt(deltaR2(a, b, scheme));
990 }
991
992
998 inline double deltaR2(const FourMomentum& v,
999 double eta2, double phi2,
1000 RapScheme scheme=PSEUDORAPIDITY) {
1001 switch (scheme) {
1002 case PSEUDORAPIDITY:
1003 return deltaR2(v.vector3(), eta2, phi2);
1004 case RAPIDITY:
1005 return deltaR2(v.rapidity(), v.azimuthalAngle(), eta2, phi2);
1006 default:
1007 throw std::runtime_error("The specified deltaR scheme is not yet implemented");
1008 }
1009 }
1010
1016 inline double deltaR(const FourMomentum& v,
1017 double eta2, double phi2,
1018 RapScheme scheme=PSEUDORAPIDITY) {
1019 return sqrt(deltaR2(v, eta2, phi2, scheme));
1020 }
1021
1022
1028 inline double deltaR2(double eta1, double phi1,
1029 const FourMomentum& v,
1030 RapScheme scheme=PSEUDORAPIDITY) {
1031 switch (scheme) {
1032 case PSEUDORAPIDITY:
1033 return deltaR2(eta1, phi1, v.vector3());
1034 case RAPIDITY:
1035 return deltaR2(eta1, phi1, v.rapidity(), v.azimuthalAngle());
1036 default:
1037 throw std::runtime_error("The specified deltaR scheme is not yet implemented");
1038 }
1039 }
1040
1046 inline double deltaR(double eta1, double phi1,
1047 const FourMomentum& v,
1048 RapScheme scheme=PSEUDORAPIDITY) {
1049 return sqrt(deltaR2(eta1, phi1, v, scheme));
1050 }
1051
1052
1058 inline double deltaR2(const FourMomentum& a, const FourVector& b,
1059 RapScheme scheme=PSEUDORAPIDITY) {
1060 switch (scheme) {
1061 case PSEUDORAPIDITY:
1062 return deltaR2(a.vector3(), b.vector3());
1063 case RAPIDITY:
1065 default:
1066 throw std::runtime_error("The specified deltaR scheme is not yet implemented");
1067 }
1068 }
1069
1075 inline double deltaR(const FourMomentum& a, const FourVector& b,
1076 RapScheme scheme=PSEUDORAPIDITY) {
1077 return sqrt(deltaR2(a, b, scheme));
1078 }
1079
1080
1086 inline double deltaR2(const FourVector& a, const FourMomentum& b,
1087 RapScheme scheme=PSEUDORAPIDITY) {
1088 return deltaR2(b, a, scheme); //< note reversed args
1089 }
1090
1096 inline double deltaR(const FourVector& a, const FourMomentum& b,
1097 RapScheme scheme=PSEUDORAPIDITY) {
1098 return deltaR(b, a, scheme); //< note reversed args
1099 }
1100
1101
1104 inline double deltaR2(const FourMomentum& a, const Vector3& b) {
1105 return deltaR2(a.vector3(), b);
1106 }
1107
1110 inline double deltaR(const FourMomentum& a, const Vector3& b) {
1111 return deltaR(a.vector3(), b);
1112 }
1113
1116 inline double deltaR2(const Vector3& a, const FourMomentum& b) {
1117 return deltaR2(a, b.vector3());
1118 }
1119
1122 inline double deltaR(const Vector3& a, const FourMomentum& b) {
1123 return deltaR(a, b.vector3());
1124 }
1125
1128 inline double deltaR2(const FourVector& a, const Vector3& b) {
1129 return deltaR2(a.vector3(), b);
1130 }
1131
1134 inline double deltaR(const FourVector& a, const Vector3& b) {
1135 return deltaR(a.vector3(), b);
1136 }
1137
1140 inline double deltaR2(const Vector3& a, const FourVector& b) {
1141 return deltaR2(a, b.vector3());
1142 }
1143
1146 inline double deltaR(const Vector3& a, const FourVector& b) {
1147 return deltaR(a, b.vector3());
1148 }
1149
1151
1152
1154
1155
1158
1160 inline double deltaPhi(const FourMomentum& a, const FourMomentum& b, bool sign=false) {
1161 return deltaPhi(a.vector3(), b.vector3(), sign);
1162 }
1163
1165 inline double deltaPhi(const FourMomentum& v, double phi2, bool sign=false) {
1166 return deltaPhi(v.vector3(), phi2, sign);
1167 }
1168
1170 inline double deltaPhi(double phi1, const FourMomentum& v, bool sign=false) {
1171 return deltaPhi(phi1, v.vector3(), sign);
1172 }
1173
1175 inline double deltaPhi(const FourVector& a, const FourVector& b, bool sign=false) {
1176 return deltaPhi(a.vector3(), b.vector3(), sign);
1177 }
1178
1180 inline double deltaPhi(const FourVector& v, double phi2, bool sign=false) {
1181 return deltaPhi(v.vector3(), phi2, sign);
1182 }
1183
1185 inline double deltaPhi(double phi1, const FourVector& v, bool sign=false) {
1186 return deltaPhi(phi1, v.vector3(), sign);
1187 }
1188
1190 inline double deltaPhi(const FourVector& a, const FourMomentum& b, bool sign=false) {
1191 return deltaPhi(a.vector3(), b.vector3(), sign);
1192 }
1193
1195 inline double deltaPhi(const FourMomentum& a, const FourVector& b, bool sign=false) {
1196 return deltaPhi(a.vector3(), b.vector3(), sign);
1197 }
1198
1200 inline double deltaPhi(const FourVector& a, const Vector3& b, bool sign=false) {
1201 return deltaPhi(a.vector3(), b, sign);
1202 }
1203
1205 inline double deltaPhi(const Vector3& a, const FourVector& b, bool sign=false) {
1206 return deltaPhi(a, b.vector3(), sign);
1207 }
1208
1210 inline double deltaPhi(const FourMomentum& a, const Vector3& b, bool sign=false) {
1211 return deltaPhi(a.vector3(), b, sign);
1212 }
1213
1215 inline double deltaPhi(const Vector3& a, const FourMomentum& b, bool sign=false) {
1216 return deltaPhi(a, b.vector3(), sign);
1217 }
1218
1220
1221
1223
1224
1227
1229 inline double deltaEta(const FourMomentum& a, const FourMomentum& b, bool sign=false) {
1230 return deltaEta(a.vector3(), b.vector3(), sign);
1231 }
1232
1234 inline double deltaEta(const FourMomentum& v, double eta2, bool sign=false) {
1235 return deltaEta(v.vector3(), eta2, sign);
1236 }
1237
1239 inline double deltaEta(double eta1, const FourMomentum& v, bool sign=false) {
1240 return deltaEta(eta1, v.vector3(), sign);
1241 }
1242
1244 inline double deltaEta(const FourVector& a, const FourVector& b, bool sign=false) {
1245 return deltaEta(a.vector3(), b.vector3(), sign);
1246 }
1247
1249 inline double deltaEta(const FourVector& v, double eta2, bool sign=false) {
1250 return deltaEta(v.vector3(), eta2, sign);
1251 }
1252
1254 inline double deltaEta(double eta1, const FourVector& v, bool sign=false) {
1255 return deltaEta(eta1, v.vector3(), sign);
1256 }
1257
1259 inline double deltaEta(const FourVector& a, const FourMomentum& b, bool sign=false) {
1260 return deltaEta(a.vector3(), b.vector3(), sign);
1261 }
1262
1264 inline double deltaEta(const FourMomentum& a, const FourVector& b, bool sign=false) {
1265 return deltaEta(a.vector3(), b.vector3(), sign);
1266 }
1267
1269 inline double deltaEta(const FourVector& a, const Vector3& b, bool sign=false) {
1270 return deltaEta(a.vector3(), b, sign);
1271 }
1272
1274 inline double deltaEta(const Vector3& a, const FourVector& b, bool sign=false) {
1275 return deltaEta(a, b.vector3(), sign);
1276 }
1277
1279 inline double deltaEta(const FourMomentum& a, const Vector3& b, bool sign=false) {
1280 return deltaEta(a.vector3(), b, sign);
1281 }
1282
1284 inline double deltaEta(const Vector3& a, const FourMomentum& b, bool sign=false) {
1285 return deltaEta(a, b.vector3(), sign);
1286 }
1287
1289
1290
1293
1295 inline double deltaRap(const FourMomentum& a, const FourMomentum& b, bool sign=false) {
1296 return deltaRap(a.rapidity(), b.rapidity(), sign);
1297 }
1298
1300 inline double deltaRap(const FourMomentum& v, double y2, bool sign=false) {
1301 return deltaRap(v.rapidity(), y2, sign);
1302 }
1303
1305 inline double deltaRap(double y1, const FourMomentum& v, bool sign=false) {
1306 return deltaRap(y1, v.rapidity(), sign);
1307 }
1308
1310
1311
1313
1314
1317
1320
1322 inline bool cmpMomByPt(const FourMomentum& a, const FourMomentum& b) {
1323 return a.pt() > b.pt();
1324 }
1326 inline bool cmpMomByAscPt(const FourMomentum& a, const FourMomentum& b) {
1327 return a.pt() < b.pt();
1328 }
1329
1331 inline bool cmpMomByP(const FourMomentum& a, const FourMomentum& b) {
1332 return a.vector3().mod() > b.vector3().mod();
1333 }
1335 inline bool cmpMomByAscP(const FourMomentum& a, const FourMomentum& b) {
1336 return a.vector3().mod() < b.vector3().mod();
1337 }
1338
1340 inline bool cmpMomByEt(const FourMomentum& a, const FourMomentum& b) {
1341 return a.Et() > b.Et();
1342 }
1344 inline bool cmpMomByAscEt(const FourMomentum& a, const FourMomentum& b) {
1345 return a.Et() < b.Et();
1346 }
1347
1349 inline bool cmpMomByE(const FourMomentum& a, const FourMomentum& b) {
1350 return a.E() > b.E();
1351 }
1353 inline bool cmpMomByAscE(const FourMomentum& a, const FourMomentum& b) {
1354 return a.E() < b.E();
1355 }
1356
1358 inline bool cmpMomByMass(const FourMomentum& a, const FourMomentum& b) {
1359 return a.mass() > b.mass();
1360 }
1362 inline bool cmpMomByAscMass(const FourMomentum& a, const FourMomentum& b) {
1363 return a.mass() < b.mass();
1364 }
1365
1367 inline bool cmpMomByEta(const FourMomentum& a, const FourMomentum& b) {
1368 return a.eta() < b.eta();
1369 }
1370
1372 inline bool cmpMomByDescEta(const FourMomentum& a, const FourMomentum& b) {
1373 return a.pseudorapidity() > b.pseudorapidity();
1374 }
1375
1377 inline bool cmpMomByAbsEta(const FourMomentum& a, const FourMomentum& b) {
1378 return fabs(a.eta()) < fabs(b.eta());
1379 }
1380
1382 inline bool cmpMomByDescAbsEta(const FourMomentum& a, const FourMomentum& b) {
1383 return fabs(a.eta()) > fabs(b.eta());
1384 }
1385
1387 inline bool cmpMomByRap(const FourMomentum& a, const FourMomentum& b) {
1388 return a.rapidity() < b.rapidity();
1389 }
1390
1392 inline bool cmpMomByDescRap(const FourMomentum& a, const FourMomentum& b) {
1393 return a.rapidity() > b.rapidity();
1394 }
1395
1397 inline bool cmpMomByAbsRap(const FourMomentum& a, const FourMomentum& b) {
1398 return fabs(a.rapidity()) < fabs(b.rapidity());
1399 }
1400
1402 inline bool cmpMomByDescAbsRap(const FourMomentum& a, const FourMomentum& b) {
1403 return fabs(a.rapidity()) > fabs(b.rapidity());
1404 }
1405
1407
1408
1410 template<typename MOMS, typename CMP>
1411 inline MOMS& isortBy(MOMS& pbs, const CMP& cmp) {
1412 std::sort(pbs.begin(), pbs.end(), cmp);
1413 return pbs;
1414 }
1416 template<typename MOMS, typename CMP>
1417 inline MOMS sortBy(const MOMS& pbs, const CMP& cmp) {
1418 MOMS rtn = pbs;
1419 std::sort(rtn.begin(), rtn.end(), cmp);
1420 return rtn;
1421 }
1422
1424 template<typename MOMS>
1425 inline MOMS& isortByPt(MOMS& pbs) {
1426 return isortBy(pbs, cmpMomByPt);
1427 }
1429 template<typename MOMS>
1430 inline MOMS sortByPt(const MOMS& pbs) {
1431 return sortBy(pbs, cmpMomByPt);
1432 }
1433
1435 template<typename MOMS>
1436 inline MOMS& isortByE(MOMS& pbs) {
1437 return isortBy(pbs, cmpMomByE);
1438 }
1440 template<typename MOMS>
1441 inline MOMS sortByE(const MOMS& pbs) {
1442 return sortBy(pbs, cmpMomByE);
1443 }
1444
1446 template<typename MOMS>
1447 inline MOMS& isortByEt(MOMS& pbs) {
1448 return isortBy(pbs, cmpMomByEt);
1449 }
1451 template<typename MOMS>
1452 inline MOMS sortByEt(const MOMS& pbs) {
1453 return sortBy(pbs, cmpMomByEt);
1454 }
1455
1457
1458
1461
1463 inline double mass(const FourMomentum& a, const FourMomentum& b) {
1464 return (a + b).mass();
1465 }
1466
1468 inline double mass2(const FourMomentum& a, const FourMomentum& b) {
1469 return (a + b).mass2();
1470 }
1471
1476 inline double mT(const FourMomentum& vis, const FourMomentum& invis) {
1477 return mT(vis.p3(), invis.p3());
1478 }
1479
1484 inline double mT(const FourMomentum& vis, const Vector3& invis) {
1485 return mT(vis.p3(), invis);
1486 }
1487
1492 inline double mT(const Vector3& vis, const FourMomentum& invis) {
1493 return mT(vis, invis.p3());
1494 }
1495
1497 inline double pT(const FourMomentum& vis, const FourMomentum& invis) {
1498 return pT(vis.p3(), invis.p3());
1499 }
1500
1502 inline double pT(const FourMomentum& vis, const Vector3& invis) {
1503 return pT(vis.p3(), invis);
1504 }
1505
1507 inline double pT(const Vector3& vis, const FourMomentum& invis) {
1508 return pT(vis, invis.p3());
1509 }
1510
1512
1513
1515
1516
1519
1521 inline std::string toString(const FourVector& lv) {
1522 std::ostringstream out;
1523 out << "(" << (fabs(lv.t()) < 1E-30 ? 0.0 : lv.t())
1524 << "; " << (fabs(lv.x()) < 1E-30 ? 0.0 : lv.x())
1525 << ", " << (fabs(lv.y()) < 1E-30 ? 0.0 : lv.y())
1526 << ", " << (fabs(lv.z()) < 1E-30 ? 0.0 : lv.z())
1527 << ")";
1528 return out.str();
1529 }
1530
1532 inline std::ostream& operator<<(std::ostream& out, const FourVector& lv) {
1533 out << toString(lv);
1534 return out;
1535 }
1536
1538
1541 typedef std::vector<FourVector> FourVectors;
1542 typedef std::vector<FourMomentum> FourMomenta;
1544
1546
1547
1548}
1549
1550#endif
Specialized version of the FourVector with momentum/energy functionality.
Definition Vector4.hh:316
double p2() const
Get the modulus-squared of the 3-momentum.
Definition Vector4.hh:601
double E2() const
Get energy-squared .
Definition Vector4.hh:555
double mass() const
Get the mass (the Lorentz self-invariant).
Definition Vector4.hh:576
double pz2() const
Get z-squared .
Definition Vector4.hh:570
double rapidity() const
Calculate the rapidity.
Definition Vector4.hh:607
FourMomentum reverse() const
Multiply space components only by -1.
Definition Vector4.hh:733
double px() const
Get x-component of momentum .
Definition Vector4.hh:558
FourMomentum & setPtPhiME(double pt, double phi, double mass, double E)
Definition Vector4.hh:532
FourMomentum & setRapPhiMPt(double y, double phi, double mass, double pt)
Definition Vector4.hh:473
Vector3 p3() const
Get 3-momentum part, .
Definition Vector4.hh:593
FourMomentum & setPz(double pz)
Set z-component of momentum .
Definition Vector4.hh:370
double p() const
Get the modulus of the 3-momentum.
Definition Vector4.hh:596
FourMomentum & setXYZE(double px, double py, double pz, double E)
Alias for setPE.
Definition Vector4.hh:384
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:779
Vector3 gammaVec() const
Definition Vector4.hh:676
double pt() const
Calculate the transverse momentum .
Definition Vector4.hh:649
FourMomentum & setPE(double px, double py, double pz, double E)
Set the p coordinates and energy simultaneously.
Definition Vector4.hh:377
double py() const
Get y-component of momentum .
Definition Vector4.hh:563
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:769
double beta() const
Definition Vector4.hh:682
double Et() const
Calculate the transverse energy .
Definition Vector4.hh:658
double Et2() const
Calculate the transverse energy .
Definition Vector4.hh:654
static FourMomentum mkXYZE(double px, double py, double pz, double E)
Make a vector from (px,py,pz,E) coordinates.
Definition Vector4.hh:749
FourMomentum & setPy(double py)
Set y-component of momentum .
Definition Vector4.hh:364
double rap() const
Alias for rapidity.
Definition Vector4.hh:613
FourMomentum & setPx(double px)
Set x-component of momentum .
Definition Vector4.hh:358
double pt2() const
Calculate the squared transverse momentum .
Definition Vector4.hh:640
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:789
Vector3 pTvec() const
Calculate the transverse momentum vector .
Definition Vector4.hh:627
FourMomentum & setRapPhiME(double y, double phi, double mass, double E)
Definition Vector4.hh:453
Vector3 ptvec() const
Synonym for pTvec.
Definition Vector4.hh:631
double mass2() const
Get the squared mass (the Lorentz self-invariant).
Definition Vector4.hh:587
FourMomentum & setEtaPhiMPt(double eta, double phi, double mass, double pt)
Definition Vector4.hh:431
double pz() const
Get z-component of momentum .
Definition Vector4.hh:568
double py2() const
Get y-squared .
Definition Vector4.hh:565
double absrap() const
Absolute rapidity.
Definition Vector4.hh:622
double pT2() const
Calculate the squared transverse momentum .
Definition Vector4.hh:636
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:759
double px2() const
Get x-squared .
Definition Vector4.hh:560
Vector3 betaVec() const
Definition Vector4.hh:688
FourMomentum & operator-=(const FourMomentum &v)
Subtract from this 4-vector. NB time as well as space components are subtracted.
Definition Vector4.hh:720
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:754
FourMomentum & setThetaPhiME(double theta, double phi, double mass, double E)
Definition Vector4.hh:490
double E() const
Get energy (time component of momentum).
Definition Vector4.hh:553
FourMomentum & setXYZM(double px, double py, double pz, double mass)
Alias for setPM.
Definition Vector4.hh:406
FourMomentum & operator+=(const FourMomentum &v)
Add to this 4-vector. NB time as well as space components are added.
Definition Vector4.hh:714
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:764
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:774
double pT() const
Calculate the transverse momentum .
Definition Vector4.hh:645
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:784
FourMomentum & setThetaPhiMPt(double theta, double phi, double mass, double pt)
Definition Vector4.hh:513
FourMomentum operator-() const
Multiply all components (time and space) by -1.
Definition Vector4.hh:726
double absrapidity() const
Absolute rapidity.
Definition Vector4.hh:618
FourMomentum & setEtaPhiME(double eta, double phi, double mass, double E)
Definition Vector4.hh:415
FourMomentum & operator/=(double a)
Divide by a scalar.
Definition Vector4.hh:708
FourMomentum & setE(double E)
Set energy (time component of momentum).
Definition Vector4.hh:352
double gamma() const
Definition Vector4.hh:670
FourMomentum & operator*=(double a)
Multiply by a scalar.
Definition Vector4.hh:702
FourMomentum & setPM(double px, double py, double pz, double mass)
Set the p coordinates and mass simultaneously.
Definition Vector4.hh:398
Specialisation of VectorN to a general (non-momentum) Lorentz 4-vector.
Definition Vector4.hh:30
FourVector operator-() const
Multiply all components (space and time) by -1.
Definition Vector4.hh:228
double perp2() const
Synonym for polarRadius2.
Definition Vector4.hh:110
FourVector & operator/=(double a)
Divide by a scalar.
Definition Vector4.hh:210
double rho() const
Synonym for polarRadius.
Definition Vector4.hh:127
Vector3 perpVec() const
Synonym for polarVec.
Definition Vector4.hh:136
Vector3 vector3() const
Get the spatial part of the 4-vector as a 3-vector.
Definition Vector4.hh:177
double angle(const FourVector &v) const
Angle between this vector and another.
Definition Vector4.hh:95
double phi(const PhiMapping mapping=ZERO_2PI) const
Synonym for azimuthalAngle.
Definition Vector4.hh:149
double contract(const FourVector &v) const
Contract two 4-vectors, with metric signature (+ - - -).
Definition Vector4.hh:188
double abseta() const
Get the directly (alias).
Definition Vector4.hh:174
double azimuthalAngle(const PhiMapping mapping=ZERO_2PI) const
Angle subtended by the 3-vector's projection in x-y and the x-axis.
Definition Vector4.hh:145
double dot(const FourVector &v) const
Contract two 4-vectors, with metric signature (+ - - -).
Definition Vector4.hh:194
double perp() const
Synonym for polarRadius.
Definition Vector4.hh:123
FourVector & operator*=(double a)
Multiply by a scalar.
Definition Vector4.hh:204
double operator*(const FourVector &v) const
Contract two 4-vectors, with metric signature (+ - - -).
Definition Vector4.hh:199
double eta() const
Synonym for pseudorapidity.
Definition Vector4.hh:167
double theta() const
Synonym for polarAngle.
Definition Vector4.hh:158
double pseudorapidity() const
Pseudorapidity (defined purely by the 3-vector components)
Definition Vector4.hh:163
double polarRadius() const
Magnitude of projection of 3-vector on to the plane.
Definition Vector4.hh:119
Vector3 rhoVec() const
Synonym for polarVec.
Definition Vector4.hh:140
double angle(const Vector3 &v3) const
Angle between this vector and another (3-vector)
Definition Vector4.hh:99
Vector3 polarVec() const
Projection of 3-vector on to the plane.
Definition Vector4.hh:132
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:106
FourVector reverse() const
Multiply space components only by -1.
Definition Vector4.hh:235
double polarAngle() const
Angle subtended by the 3-vector and the z-axis.
Definition Vector4.hh:154
FourVector & operator-=(const FourVector &v)
Subtract from this 4-vector. NB time as well as space components are subtracted.
Definition Vector4.hh:222
double abspseudorapidity() const
Get the directly.
Definition Vector4.hh:172
FourVector & operator+=(const FourVector &v)
Add to this 4-vector.
Definition Vector4.hh:216
double rho2() const
Synonym for polarRadius2.
Definition Vector4.hh:114
Object implementing Lorentz transform calculations and boosts.
Definition LorentzTrans.hh:21
Three-dimensional specialisation of Vector.
Definition Vector3.hh:40
Vector3 rhoVec() const
Synonym for polarVec.
Definition Vector3.hh:140
double polarRadius() const
Polar radius.
Definition Vector3.hh:158
double rho() const
Synonym for polarRadius.
Definition Vector3.hh:166
double eta() const
Synonym for pseudorapidity.
Definition Vector3.hh:220
double perp() const
Synonym for polarRadius.
Definition Vector3.hh:162
double perp2() const
Synonym for polarRadius2.
Definition Vector3.hh:149
double pseudorapidity() const
Purely geometric approximation to rapidity.
Definition Vector3.hh:212
Vector3 unit() const
Synonym for unitVec.
Definition Vector3.hh:124
double theta() const
Synonym for polarAngle.
Definition Vector3.hh:200
double phi(const PhiMapping mapping=ZERO_2PI) const
Synonym for azimuthalAngle.
Definition Vector3.hh:183
double polarRadius2() const
Square of the polar radius (.
Definition Vector3.hh:145
Vector3 perpVec() const
Synonym for polarVec.
Definition Vector3.hh:136
double azimuthalAngle(const PhiMapping mapping=ZERO_2PI) const
Angle subtended by the vector's projection in x-y and the x-axis.
Definition Vector3.hh:174
double rho2() const
Synonym for polarRadius2.
Definition Vector3.hh:153
double angle(const Vector3 &v) const
Angle in radians to another vector.
Definition Vector3.hh:108
double polarAngle() const
Angle subtended by the vector and the z-axis.
Definition Vector3.hh:193
Vector3 polarVec() const
Polar projection of this vector into the x-y plane.
Definition Vector3.hh:130
A minimal base class for -dimensional vectors.
Definition VectorN.hh:23
double mod() const
Calculate the modulus of a vector. .
Definition VectorN.hh:95
double mod2() const
Calculate the modulus-squared of a vector. .
Definition VectorN.hh:84
Vector< N > & set(const size_t index, const double value)
Set indexed value.
Definition VectorN.hh:60
MOMS & isortBy(MOMS &pbs, const CMP &cmp)
Sort a container of momenta by cmp and return by reference for non-const inputs.
Definition Vector4.hh:1411
MOMS & isortByE(MOMS &pbs)
Sort a container of momenta by E (decreasing) and return by reference for non-const inputs.
Definition Vector4.hh:1436
MOMS & isortByEt(MOMS &pbs)
Sort a container of momenta by Et (decreasing) and return by reference for non-const inputs.
Definition Vector4.hh:1447
bool cmpMomByDescEta(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by decreasing eta (pseudorapidity)
Definition Vector4.hh:1372
bool cmpMomByMass(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by decreasing mass.
Definition Vector4.hh:1358
bool cmpMomByP(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by decreasing 3-momentum magnitude |p|.
Definition Vector4.hh:1331
bool cmpMomByAscEt(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by increasing transverse energy.
Definition Vector4.hh:1344
bool cmpMomByRap(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by increasing rapidity.
Definition Vector4.hh:1387
bool cmpMomByE(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by decreasing energy.
Definition Vector4.hh:1349
MOMS sortBy(const MOMS &pbs, const CMP &cmp)
Sort a container of momenta by cmp and return by value for const inputs.
Definition Vector4.hh:1417
bool cmpMomByAbsRap(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by increasing absolute rapidity.
Definition Vector4.hh:1397
MOMS sortByPt(const MOMS &pbs)
Sort a container of momenta by pT (decreasing) and return by value for const inputs.
Definition Vector4.hh:1430
MOMS sortByEt(const MOMS &pbs)
Sort a container of momenta by Et (decreasing) and return by value for const inputs.
Definition Vector4.hh:1452
bool cmpMomByEt(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by decreasing transverse energy.
Definition Vector4.hh:1340
bool cmpMomByAscPt(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by increasing pT.
Definition Vector4.hh:1326
bool cmpMomByAscP(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by increasing 3-momentum magnitude |p|.
Definition Vector4.hh:1335
bool cmpMomByAscMass(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by increasing mass.
Definition Vector4.hh:1362
MOMS sortByE(const MOMS &pbs)
Sort a container of momenta by E (decreasing) and return by value for const inputs.
Definition Vector4.hh:1441
bool cmpMomByPt(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by decreasing pT.
Definition Vector4.hh:1322
bool cmpMomByAbsEta(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by increasing absolute eta (pseudorapidity)
Definition Vector4.hh:1377
bool cmpMomByDescAbsEta(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by increasing absolute eta (pseudorapidity)
Definition Vector4.hh:1382
bool cmpMomByEta(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by increasing eta (pseudorapidity)
Definition Vector4.hh:1367
bool cmpMomByAscE(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by increasing energy.
Definition Vector4.hh:1353
bool cmpMomByDescAbsRap(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by decreasing absolute rapidity.
Definition Vector4.hh:1402
MOMS & isortByPt(MOMS &pbs)
Sort a container of momenta by pT (decreasing) and return by reference for non-const inputs.
Definition Vector4.hh:1425
bool cmpMomByDescRap(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by decreasing rapidity.
Definition Vector4.hh:1392
double mass(const FourMomentum &a, const FourMomentum &b)
Calculate mass of two 4-vectors.
Definition Vector4.hh:1463
double mass2(const FourMomentum &a, const FourMomentum &b)
Calculate mass^2 of two 4-vectors.
Definition Vector4.hh:1468
double pT(const Vector3 &a, const Vector3 &b)
Calculate transverse momentum of pair of 3-vectors.
Definition Vector3.hh:640
std::vector< FourVector > FourVectors
Definition Vector4.hh:1541
string to_str(const T &x)
Convert any object to a string.
Definition Utils.hh:67
Definition MC_CENT_PPB_Projections.hh:10
constexpr std::enable_if_t< std::is_arithmetic_v< NUM >, int > sign(NUM val)
Find the sign of a number.
Definition MathUtils.hh:264
double deltaR(double rap1, double phi1, double rap2, double phi2)
Definition MathUtils.hh:697
double deltaPhi(double phi1, double phi2, bool sign=false)
Calculate the difference between two angles in radians.
Definition MathUtils.hh:667
double deltaEta(double eta1, double eta2, bool sign=false)
Definition MathUtils.hh:675
PhiMapping
Enum for range of to be mapped into.
Definition MathConstants.hh:49
double deltaR2(double rap1, double phi1, double rap2, double phi2)
Definition MathUtils.hh:690
double mT(double pT1, double pT2, double dphi)
Definition MathUtils.hh:719
std::ostream & operator<<(std::ostream &os, const AnalysisInfo &ai)
Stream an AnalysisInfo as a text description.
Definition AnalysisInfo.hh:362
std::enable_if_t< std::is_floating_point_v< NUM >, bool > isZero(NUM val, double tolerance=1e-8)
Compare a number to zero.
Definition MathUtils.hh:24
std::enable_if_t< std::is_arithmetic_v< NUM >, NUM > sqr(NUM a)
Named number-type squaring operation.
Definition MathUtils.hh:218
double invariant(const FourVector &lv)
Definition Vector4.hh:292
double contract(const FourVector &a, const FourVector &b)
Contract two 4-vectors, with metric signature (+ - - -).
Definition Vector4.hh:245
RapScheme
Enum for rapidity variable to be used in calculating , applying rapidity cuts, etc.
Definition MathConstants.hh:46
double deltaRap(double y1, double y2, bool sign=false)
Definition MathUtils.hh:683
Cmp< T > cmp(const T &t1, const T &t2)
Global helper function for easy creation of Cmp objects.
Definition Cmp.hh:255
std::string toString(const AnalysisInfo &ai)
String representation.
double angle(const Vector2 &a, const Vector2 &b)
Angle (in radians) between two 2-vectors.
Definition Vector2.hh:177
double rapidity(double E, double pz)
Calculate a rapidity value from the supplied energy E and longitudinal momentum pz.
Definition MathUtils.hh:702