rivet is hosted by Hepforge, IPPP Durham
Rivet 4.0.0
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 return 0.5 * std::log( (E() + pz()) / (E() - pz()) );
609 }
611 double rap() const {
612 return rapidity();
613 }
614
616 double absrapidity() const {
617 return fabs(rapidity());
618 }
620 double absrap() const {
621 return fabs(rap());
622 }
623
625 Vector3 pTvec() const {
626 return p3().polarVec();
627 }
629 Vector3 ptvec() const {
630 return pTvec();
631 }
632
634 double pT2() const {
635 return vector3().polarRadius2();
636 }
638 double pt2() const {
639 return vector3().polarRadius2();
640 }
641
643 double pT() const {
644 return sqrt(pT2());
645 }
647 double pt() const {
648 return sqrt(pT2());
649 }
650
652 double Et2() const {
653 return Et() * Et();
654 }
656 double Et() const {
657 return E() * sin(polarAngle());
658 }
659
661
662
665
668 double gamma() const {
669 return sqrt(E2()/mass2());
670 }
671
675 return gamma() * p3().unit();
676 }
677
680 double beta() const {
681 return p()/E();
682 }
683
686 Vector3 betaVec() const {
687 // return Vector3(px()/E(), py()/E(), pz()/E());
688 return p3()/E();
689 }
690
692
693
695
698
701 _vec = multiply(a, *this)._vec;
702 return *this;
703 }
704
707 _vec = multiply(1.0/a, *this)._vec;
708 return *this;
709 }
710
713 _vec = add(*this, v)._vec;
714 return *this;
715 }
716
719 _vec = add(*this, -v)._vec;
720 return *this;
721 }
722
726 result._vec = -_vec;
727 return result;
728 }
729
732 FourMomentum result = -*this;
733 result.setE(-result.E());
734 return result;
735 }
736
738
739
741
742
745
747 static FourMomentum mkXYZE(double px, double py, double pz, double E) {
748 return FourMomentum().setPE(px, py, pz, E);
749 }
750
752 static FourMomentum mkXYZM(double px, double py, double pz, double mass) {
753 return FourMomentum().setPM(px, py, pz, mass);
754 }
755
757 static FourMomentum mkEtaPhiME(double eta, double phi, double mass, double E) {
758 return FourMomentum().setEtaPhiME(eta, phi, mass, E);
759 }
760
762 static FourMomentum mkEtaPhiMPt(double eta, double phi, double mass, double pt) {
763 return FourMomentum().setEtaPhiMPt(eta, phi, mass, pt);
764 }
765
767 static FourMomentum mkRapPhiME(double y, double phi, double mass, double E) {
768 return FourMomentum().setRapPhiME(y, phi, mass, E);
769 }
770
772 static FourMomentum mkRapPhiMPt(double y, double phi, double mass, double pt) {
773 return FourMomentum().setRapPhiMPt(y, phi, mass, pt);
774 }
775
777 static FourMomentum mkThetaPhiME(double theta, double phi, double mass, double E) {
779 }
780
782 static FourMomentum mkThetaPhiMPt(double theta, double phi, double mass, double pt) {
784 }
785
787 static FourMomentum mkPtPhiME(double pt, double phi, double mass, double E) {
788 return FourMomentum().setPtPhiME(pt, phi, mass, E);
789 }
790
792
793
794 };
795
796
797 inline FourMomentum multiply(const double a, const FourMomentum& v) {
798 FourMomentum result;
799 result._vec = a * v._vec;
800 return result;
801 }
802
803 inline FourMomentum multiply(const FourMomentum& v, const double a) {
804 return multiply(a, v);
805 }
806
807 inline FourMomentum operator*(const double a, const FourMomentum& v) {
808 return multiply(a, v);
809 }
810
811 inline FourMomentum operator*(const FourMomentum& v, const double a) {
812 return multiply(a, v);
813 }
814
815 inline FourMomentum operator/(const FourMomentum& v, const double a) {
816 return multiply(1.0/a, v);
817 }
818
819 inline FourMomentum add(const FourMomentum& a, const FourMomentum& b) {
820 FourMomentum result;
821 result._vec = a._vec + b._vec;
822 return result;
823 }
824
825 inline FourMomentum operator+(const FourMomentum& a, const FourMomentum& b) {
826 return add(a, b);
827 }
828
829 inline FourMomentum operator-(const FourMomentum& a, const FourMomentum& b) {
830 return add(a, -b);
831 }
832
833
835
836
839
848 inline double deltaR2(const FourVector& a, const FourVector& b,
849 RapScheme scheme=PSEUDORAPIDITY) {
850 switch (scheme) {
851 case PSEUDORAPIDITY :
852 return deltaR2(a.vector3(), b.vector3());
853 case RAPIDITY:
854 {
855 const FourMomentum* ma = dynamic_cast<const FourMomentum*>(&a);
856 const FourMomentum* mb = dynamic_cast<const FourMomentum*>(&b);
857 if (!ma || !mb) {
858 string err = "deltaR with scheme RAPIDITY can only be called with FourMomentum objects, not FourVectors";
859 throw std::runtime_error(err);
860 }
861 return deltaR2(*ma, *mb, scheme);
862 }
863 default:
864 throw std::runtime_error("The specified deltaR scheme is not yet implemented");
865 }
866 }
867
876 inline double deltaR(const FourVector& a, const FourVector& b,
877 RapScheme scheme=PSEUDORAPIDITY) {
878 return sqrt(deltaR2(a, b, scheme));
879 }
880
881
882
889 inline double deltaR2(const FourVector& v,
890 double eta2, double phi2,
891 RapScheme scheme=PSEUDORAPIDITY) {
892 switch (scheme) {
893 case PSEUDORAPIDITY :
894 return deltaR2(v.vector3(), eta2, phi2);
895 case RAPIDITY:
896 {
897 const FourMomentum* mv = dynamic_cast<const FourMomentum*>(&v);
898 if (!mv) {
899 string err = "deltaR with scheme RAPIDITY can only be called with FourMomentum objects, not FourVectors";
900 throw std::runtime_error(err);
901 }
902 return deltaR2(*mv, eta2, phi2, scheme);
903 }
904 default:
905 throw std::runtime_error("The specified deltaR scheme is not yet implemented");
906 }
907 }
908
915 inline double deltaR(const FourVector& v,
916 double eta2, double phi2,
917 RapScheme scheme=PSEUDORAPIDITY) {
918 return sqrt(deltaR2(v, eta2, phi2, scheme));
919 }
920
921
928 inline double deltaR2(double eta1, double phi1,
929 const FourVector& v,
930 RapScheme scheme=PSEUDORAPIDITY) {
931 switch (scheme) {
932 case PSEUDORAPIDITY :
933 return deltaR2(eta1, phi1, v.vector3());
934 case RAPIDITY:
935 {
936 const FourMomentum* mv = dynamic_cast<const FourMomentum*>(&v);
937 if (!mv) {
938 string err = "deltaR with scheme RAPIDITY can only be called with FourMomentum objects, not FourVectors";
939 throw std::runtime_error(err);
940 }
941 return deltaR2(eta1, phi1, *mv, scheme);
942 }
943 default:
944 throw std::runtime_error("The specified deltaR scheme is not yet implemented");
945 }
946 }
947
954 inline double deltaR(double eta1, double phi1,
955 const FourVector& v,
956 RapScheme scheme=PSEUDORAPIDITY) {
957 return sqrt(deltaR2(eta1, phi1, v, scheme));
958 }
959
960
967 inline double deltaR2(const FourMomentum& a, const FourMomentum& b,
968 RapScheme scheme=PSEUDORAPIDITY) {
969 switch (scheme) {
970 case PSEUDORAPIDITY:
971 return deltaR2(a.vector3(), b.vector3());
972 case RAPIDITY:
973 return deltaR2(a.rapidity(), a.azimuthalAngle(), b.rapidity(), b.azimuthalAngle());
974 default:
975 throw std::runtime_error("The specified deltaR scheme is not yet implemented");
976 }
977 }
978
985 inline double deltaR(const FourMomentum& a, const FourMomentum& b,
986 RapScheme scheme=PSEUDORAPIDITY) {
987 return sqrt(deltaR2(a, b, scheme));
988 }
989
990
996 inline double deltaR2(const FourMomentum& v,
997 double eta2, double phi2,
998 RapScheme scheme=PSEUDORAPIDITY) {
999 switch (scheme) {
1000 case PSEUDORAPIDITY:
1001 return deltaR2(v.vector3(), eta2, phi2);
1002 case RAPIDITY:
1003 return deltaR2(v.rapidity(), v.azimuthalAngle(), eta2, phi2);
1004 default:
1005 throw std::runtime_error("The specified deltaR scheme is not yet implemented");
1006 }
1007 }
1008
1014 inline double deltaR(const FourMomentum& v,
1015 double eta2, double phi2,
1016 RapScheme scheme=PSEUDORAPIDITY) {
1017 return sqrt(deltaR2(v, eta2, phi2, scheme));
1018 }
1019
1020
1026 inline double deltaR2(double eta1, double phi1,
1027 const FourMomentum& v,
1028 RapScheme scheme=PSEUDORAPIDITY) {
1029 switch (scheme) {
1030 case PSEUDORAPIDITY:
1031 return deltaR2(eta1, phi1, v.vector3());
1032 case RAPIDITY:
1033 return deltaR2(eta1, phi1, v.rapidity(), v.azimuthalAngle());
1034 default:
1035 throw std::runtime_error("The specified deltaR scheme is not yet implemented");
1036 }
1037 }
1038
1044 inline double deltaR(double eta1, double phi1,
1045 const FourMomentum& v,
1046 RapScheme scheme=PSEUDORAPIDITY) {
1047 return sqrt(deltaR2(eta1, phi1, v, scheme));
1048 }
1049
1050
1056 inline double deltaR2(const FourMomentum& a, const FourVector& b,
1057 RapScheme scheme=PSEUDORAPIDITY) {
1058 switch (scheme) {
1059 case PSEUDORAPIDITY:
1060 return deltaR2(a.vector3(), b.vector3());
1061 case RAPIDITY:
1063 default:
1064 throw std::runtime_error("The specified deltaR scheme is not yet implemented");
1065 }
1066 }
1067
1073 inline double deltaR(const FourMomentum& a, const FourVector& b,
1074 RapScheme scheme=PSEUDORAPIDITY) {
1075 return sqrt(deltaR2(a, b, scheme));
1076 }
1077
1078
1084 inline double deltaR2(const FourVector& a, const FourMomentum& b,
1085 RapScheme scheme=PSEUDORAPIDITY) {
1086 return deltaR2(b, a, scheme); //< note reversed args
1087 }
1088
1094 inline double deltaR(const FourVector& a, const FourMomentum& b,
1095 RapScheme scheme=PSEUDORAPIDITY) {
1096 return deltaR(b, a, scheme); //< note reversed args
1097 }
1098
1099
1102 inline double deltaR2(const FourMomentum& a, const Vector3& b) {
1103 return deltaR2(a.vector3(), b);
1104 }
1105
1108 inline double deltaR(const FourMomentum& a, const Vector3& b) {
1109 return deltaR(a.vector3(), b);
1110 }
1111
1114 inline double deltaR2(const Vector3& a, const FourMomentum& b) {
1115 return deltaR2(a, b.vector3());
1116 }
1117
1120 inline double deltaR(const Vector3& a, const FourMomentum& b) {
1121 return deltaR(a, b.vector3());
1122 }
1123
1126 inline double deltaR2(const FourVector& a, const Vector3& b) {
1127 return deltaR2(a.vector3(), b);
1128 }
1129
1132 inline double deltaR(const FourVector& a, const Vector3& b) {
1133 return deltaR(a.vector3(), b);
1134 }
1135
1138 inline double deltaR2(const Vector3& a, const FourVector& b) {
1139 return deltaR2(a, b.vector3());
1140 }
1141
1144 inline double deltaR(const Vector3& a, const FourVector& b) {
1145 return deltaR(a, b.vector3());
1146 }
1147
1149
1150
1152
1153
1156
1158 inline double deltaPhi(const FourMomentum& a, const FourMomentum& b, bool sign=false) {
1159 return deltaPhi(a.vector3(), b.vector3(), sign);
1160 }
1161
1163 inline double deltaPhi(const FourMomentum& v, double phi2, bool sign=false) {
1164 return deltaPhi(v.vector3(), phi2, sign);
1165 }
1166
1168 inline double deltaPhi(double phi1, const FourMomentum& v, bool sign=false) {
1169 return deltaPhi(phi1, v.vector3(), sign);
1170 }
1171
1173 inline double deltaPhi(const FourVector& a, const FourVector& b, bool sign=false) {
1174 return deltaPhi(a.vector3(), b.vector3(), sign);
1175 }
1176
1178 inline double deltaPhi(const FourVector& v, double phi2, bool sign=false) {
1179 return deltaPhi(v.vector3(), phi2, sign);
1180 }
1181
1183 inline double deltaPhi(double phi1, const FourVector& v, bool sign=false) {
1184 return deltaPhi(phi1, v.vector3(), sign);
1185 }
1186
1188 inline double deltaPhi(const FourVector& a, const FourMomentum& b, bool sign=false) {
1189 return deltaPhi(a.vector3(), b.vector3(), sign);
1190 }
1191
1193 inline double deltaPhi(const FourMomentum& a, const FourVector& b, bool sign=false) {
1194 return deltaPhi(a.vector3(), b.vector3(), sign);
1195 }
1196
1198 inline double deltaPhi(const FourVector& a, const Vector3& b, bool sign=false) {
1199 return deltaPhi(a.vector3(), b, sign);
1200 }
1201
1203 inline double deltaPhi(const Vector3& a, const FourVector& b, bool sign=false) {
1204 return deltaPhi(a, b.vector3(), sign);
1205 }
1206
1208 inline double deltaPhi(const FourMomentum& a, const Vector3& b, bool sign=false) {
1209 return deltaPhi(a.vector3(), b, sign);
1210 }
1211
1213 inline double deltaPhi(const Vector3& a, const FourMomentum& b, bool sign=false) {
1214 return deltaPhi(a, b.vector3(), sign);
1215 }
1216
1218
1219
1221
1222
1225
1227 inline double deltaEta(const FourMomentum& a, const FourMomentum& b, bool sign=false) {
1228 return deltaEta(a.vector3(), b.vector3(), sign);
1229 }
1230
1232 inline double deltaEta(const FourMomentum& v, double eta2, bool sign=false) {
1233 return deltaEta(v.vector3(), eta2, sign);
1234 }
1235
1237 inline double deltaEta(double eta1, const FourMomentum& v, bool sign=false) {
1238 return deltaEta(eta1, v.vector3(), sign);
1239 }
1240
1242 inline double deltaEta(const FourVector& a, const FourVector& b, bool sign=false) {
1243 return deltaEta(a.vector3(), b.vector3(), sign);
1244 }
1245
1247 inline double deltaEta(const FourVector& v, double eta2, bool sign=false) {
1248 return deltaEta(v.vector3(), eta2, sign);
1249 }
1250
1252 inline double deltaEta(double eta1, const FourVector& v, bool sign=false) {
1253 return deltaEta(eta1, v.vector3(), sign);
1254 }
1255
1257 inline double deltaEta(const FourVector& a, const FourMomentum& b, bool sign=false) {
1258 return deltaEta(a.vector3(), b.vector3(), sign);
1259 }
1260
1262 inline double deltaEta(const FourMomentum& a, const FourVector& b, bool sign=false) {
1263 return deltaEta(a.vector3(), b.vector3(), sign);
1264 }
1265
1267 inline double deltaEta(const FourVector& a, const Vector3& b, bool sign=false) {
1268 return deltaEta(a.vector3(), b, sign);
1269 }
1270
1272 inline double deltaEta(const Vector3& a, const FourVector& b, bool sign=false) {
1273 return deltaEta(a, b.vector3(), sign);
1274 }
1275
1277 inline double deltaEta(const FourMomentum& a, const Vector3& b, bool sign=false) {
1278 return deltaEta(a.vector3(), b, sign);
1279 }
1280
1282 inline double deltaEta(const Vector3& a, const FourMomentum& b, bool sign=false) {
1283 return deltaEta(a, b.vector3(), sign);
1284 }
1285
1287
1288
1291
1293 inline double deltaRap(const FourMomentum& a, const FourMomentum& b, bool sign=false) {
1294 return deltaRap(a.rapidity(), b.rapidity(), sign);
1295 }
1296
1298 inline double deltaRap(const FourMomentum& v, double y2, bool sign=false) {
1299 return deltaRap(v.rapidity(), y2, sign);
1300 }
1301
1303 inline double deltaRap(double y1, const FourMomentum& v, bool sign=false) {
1304 return deltaRap(y1, v.rapidity(), sign);
1305 }
1306
1308
1309
1311
1312
1315
1318
1320 inline bool cmpMomByPt(const FourMomentum& a, const FourMomentum& b) {
1321 return a.pt() > b.pt();
1322 }
1324 inline bool cmpMomByAscPt(const FourMomentum& a, const FourMomentum& b) {
1325 return a.pt() < b.pt();
1326 }
1327
1329 inline bool cmpMomByP(const FourMomentum& a, const FourMomentum& b) {
1330 return a.vector3().mod() > b.vector3().mod();
1331 }
1333 inline bool cmpMomByAscP(const FourMomentum& a, const FourMomentum& b) {
1334 return a.vector3().mod() < b.vector3().mod();
1335 }
1336
1338 inline bool cmpMomByEt(const FourMomentum& a, const FourMomentum& b) {
1339 return a.Et() > b.Et();
1340 }
1342 inline bool cmpMomByAscEt(const FourMomentum& a, const FourMomentum& b) {
1343 return a.Et() < b.Et();
1344 }
1345
1347 inline bool cmpMomByE(const FourMomentum& a, const FourMomentum& b) {
1348 return a.E() > b.E();
1349 }
1351 inline bool cmpMomByAscE(const FourMomentum& a, const FourMomentum& b) {
1352 return a.E() < b.E();
1353 }
1354
1356 inline bool cmpMomByMass(const FourMomentum& a, const FourMomentum& b) {
1357 return a.mass() > b.mass();
1358 }
1360 inline bool cmpMomByAscMass(const FourMomentum& a, const FourMomentum& b) {
1361 return a.mass() < b.mass();
1362 }
1363
1365 inline bool cmpMomByEta(const FourMomentum& a, const FourMomentum& b) {
1366 return a.eta() < b.eta();
1367 }
1368
1370 inline bool cmpMomByDescEta(const FourMomentum& a, const FourMomentum& b) {
1371 return a.pseudorapidity() > b.pseudorapidity();
1372 }
1373
1375 inline bool cmpMomByAbsEta(const FourMomentum& a, const FourMomentum& b) {
1376 return fabs(a.eta()) < fabs(b.eta());
1377 }
1378
1380 inline bool cmpMomByDescAbsEta(const FourMomentum& a, const FourMomentum& b) {
1381 return fabs(a.eta()) > fabs(b.eta());
1382 }
1383
1385 inline bool cmpMomByRap(const FourMomentum& a, const FourMomentum& b) {
1386 return a.rapidity() < b.rapidity();
1387 }
1388
1390 inline bool cmpMomByDescRap(const FourMomentum& a, const FourMomentum& b) {
1391 return a.rapidity() > b.rapidity();
1392 }
1393
1395 inline bool cmpMomByAbsRap(const FourMomentum& a, const FourMomentum& b) {
1396 return fabs(a.rapidity()) < fabs(b.rapidity());
1397 }
1398
1400 inline bool cmpMomByDescAbsRap(const FourMomentum& a, const FourMomentum& b) {
1401 return fabs(a.rapidity()) > fabs(b.rapidity());
1402 }
1403
1405
1406
1408 template<typename MOMS, typename CMP>
1409 inline MOMS& isortBy(MOMS& pbs, const CMP& cmp) {
1410 std::sort(pbs.begin(), pbs.end(), cmp);
1411 return pbs;
1412 }
1414 template<typename MOMS, typename CMP>
1415 inline MOMS sortBy(const MOMS& pbs, const CMP& cmp) {
1416 MOMS rtn = pbs;
1417 std::sort(rtn.begin(), rtn.end(), cmp);
1418 return rtn;
1419 }
1420
1422 template<typename MOMS>
1423 inline MOMS& isortByPt(MOMS& pbs) {
1424 return isortBy(pbs, cmpMomByPt);
1425 }
1427 template<typename MOMS>
1428 inline MOMS sortByPt(const MOMS& pbs) {
1429 return sortBy(pbs, cmpMomByPt);
1430 }
1431
1433 template<typename MOMS>
1434 inline MOMS& isortByE(MOMS& pbs) {
1435 return isortBy(pbs, cmpMomByE);
1436 }
1438 template<typename MOMS>
1439 inline MOMS sortByE(const MOMS& pbs) {
1440 return sortBy(pbs, cmpMomByE);
1441 }
1442
1444 template<typename MOMS>
1445 inline MOMS& isortByEt(MOMS& pbs) {
1446 return isortBy(pbs, cmpMomByEt);
1447 }
1449 template<typename MOMS>
1450 inline MOMS sortByEt(const MOMS& pbs) {
1451 return sortBy(pbs, cmpMomByEt);
1452 }
1453
1455
1456
1459
1461 inline double mass(const FourMomentum& a, const FourMomentum& b) {
1462 return (a + b).mass();
1463 }
1464
1466 inline double mass2(const FourMomentum& a, const FourMomentum& b) {
1467 return (a + b).mass2();
1468 }
1469
1474 inline double mT(const FourMomentum& vis, const FourMomentum& invis) {
1475 return mT(vis.p3(), invis.p3());
1476 }
1477
1482 inline double mT(const FourMomentum& vis, const Vector3& invis) {
1483 return mT(vis.p3(), invis);
1484 }
1485
1490 inline double mT(const Vector3& vis, const FourMomentum& invis) {
1491 return mT(vis, invis.p3());
1492 }
1493
1495 inline double pT(const FourMomentum& vis, const FourMomentum& invis) {
1496 return pT(vis.p3(), invis.p3());
1497 }
1498
1500 inline double pT(const FourMomentum& vis, const Vector3& invis) {
1501 return pT(vis.p3(), invis);
1502 }
1503
1505 inline double pT(const Vector3& vis, const FourMomentum& invis) {
1506 return pT(vis, invis.p3());
1507 }
1508
1510
1511
1513
1514
1517
1519 inline std::string toString(const FourVector& lv) {
1520 std::ostringstream out;
1521 out << "(" << (fabs(lv.t()) < 1E-30 ? 0.0 : lv.t())
1522 << "; " << (fabs(lv.x()) < 1E-30 ? 0.0 : lv.x())
1523 << ", " << (fabs(lv.y()) < 1E-30 ? 0.0 : lv.y())
1524 << ", " << (fabs(lv.z()) < 1E-30 ? 0.0 : lv.z())
1525 << ")";
1526 return out.str();
1527 }
1528
1530 inline std::ostream& operator<<(std::ostream& out, const FourVector& lv) {
1531 out << toString(lv);
1532 return out;
1533 }
1534
1536
1539 typedef std::vector<FourVector> FourVectors;
1540 typedef std::vector<FourMomentum> FourMomenta;
1542
1544
1545
1546}
1547
1548#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:731
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:777
Vector3 gammaVec() const
Definition Vector4.hh:674
double pt() const
Calculate the transverse momentum .
Definition Vector4.hh:647
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:767
double beta() const
Definition Vector4.hh:680
double Et() const
Calculate the transverse energy .
Definition Vector4.hh:656
double Et2() const
Calculate the transverse energy .
Definition Vector4.hh:652
static FourMomentum mkXYZE(double px, double py, double pz, double E)
Make a vector from (px,py,pz,E) coordinates.
Definition Vector4.hh:747
FourMomentum & setPy(double py)
Set y-component of momentum .
Definition Vector4.hh:364
double rap() const
Alias for rapidity.
Definition Vector4.hh:611
FourMomentum & setPx(double px)
Set x-component of momentum .
Definition Vector4.hh:358
double pt2() const
Calculate the squared transverse momentum .
Definition Vector4.hh:638
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:787
Vector3 pTvec() const
Calculate the transverse momentum vector .
Definition Vector4.hh:625
FourMomentum & setRapPhiME(double y, double phi, double mass, double E)
Definition Vector4.hh:453
Vector3 ptvec() const
Synonym for pTvec.
Definition Vector4.hh:629
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:620
double pT2() const
Calculate the squared transverse momentum .
Definition Vector4.hh:634
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:757
double px2() const
Get x-squared .
Definition Vector4.hh:560
Vector3 betaVec() const
Definition Vector4.hh:686
FourMomentum & operator-=(const FourMomentum &v)
Subtract from this 4-vector. NB time as well as space components are subtracted.
Definition Vector4.hh:718
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:752
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:712
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:762
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:772
double pT() const
Calculate the transverse momentum .
Definition Vector4.hh:643
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:782
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:724
double absrapidity() const
Absolute rapidity.
Definition Vector4.hh:616
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:706
FourMomentum & setE(double E)
Set energy (time component of momentum).
Definition Vector4.hh:352
double gamma() const
Definition Vector4.hh:668
FourMomentum & operator*=(double a)
Multiply by a scalar.
Definition Vector4.hh:700
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:219
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:1409
MOMS & isortByE(MOMS &pbs)
Sort a container of momenta by E (decreasing) and return by reference for non-const inputs.
Definition Vector4.hh:1434
MOMS & isortByEt(MOMS &pbs)
Sort a container of momenta by Et (decreasing) and return by reference for non-const inputs.
Definition Vector4.hh:1445
bool cmpMomByDescEta(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by decreasing eta (pseudorapidity)
Definition Vector4.hh:1370
bool cmpMomByMass(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by decreasing mass.
Definition Vector4.hh:1356
bool cmpMomByP(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by decreasing 3-momentum magnitude |p|.
Definition Vector4.hh:1329
bool cmpMomByAscEt(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by increasing transverse energy.
Definition Vector4.hh:1342
bool cmpMomByRap(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by increasing rapidity.
Definition Vector4.hh:1385
bool cmpMomByE(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by decreasing energy.
Definition Vector4.hh:1347
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:1415
bool cmpMomByAbsRap(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by increasing absolute rapidity.
Definition Vector4.hh:1395
MOMS sortByPt(const MOMS &pbs)
Sort a container of momenta by pT (decreasing) and return by value for const inputs.
Definition Vector4.hh:1428
MOMS sortByEt(const MOMS &pbs)
Sort a container of momenta by Et (decreasing) and return by value for const inputs.
Definition Vector4.hh:1450
bool cmpMomByEt(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by decreasing transverse energy.
Definition Vector4.hh:1338
bool cmpMomByAscPt(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by increasing pT.
Definition Vector4.hh:1324
bool cmpMomByAscP(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by increasing 3-momentum magnitude |p|.
Definition Vector4.hh:1333
bool cmpMomByAscMass(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by increasing mass.
Definition Vector4.hh:1360
MOMS sortByE(const MOMS &pbs)
Sort a container of momenta by E (decreasing) and return by value for const inputs.
Definition Vector4.hh:1439
bool cmpMomByPt(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by decreasing pT.
Definition Vector4.hh:1320
bool cmpMomByAbsEta(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by increasing absolute eta (pseudorapidity)
Definition Vector4.hh:1375
bool cmpMomByDescAbsEta(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by increasing absolute eta (pseudorapidity)
Definition Vector4.hh:1380
bool cmpMomByEta(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by increasing eta (pseudorapidity)
Definition Vector4.hh:1365
bool cmpMomByAscE(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by increasing energy.
Definition Vector4.hh:1351
bool cmpMomByDescAbsRap(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by decreasing absolute rapidity.
Definition Vector4.hh:1400
MOMS & isortByPt(MOMS &pbs)
Sort a container of momenta by pT (decreasing) and return by reference for non-const inputs.
Definition Vector4.hh:1423
bool cmpMomByDescRap(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by decreasing rapidity.
Definition Vector4.hh:1390
double mass(const FourMomentum &a, const FourMomentum &b)
Calculate mass of two 4-vectors.
Definition Vector4.hh:1461
double mass2(const FourMomentum &a, const FourMomentum &b)
Calculate mass^2 of two 4-vectors.
Definition Vector4.hh:1466
double pT(const Vector3 &a, const Vector3 &b)
Calculate transverse momentum of pair of 3-vectors.
Definition Vector3.hh:639
std::vector< FourVector > FourVectors
Definition Vector4.hh:1539
string to_str(const T &x)
Convert any object to a string.
Definition Utils.hh:67
Definition MC_CENT_PPB_Projections.hh:10
double deltaR(double rap1, double phi1, double rap2, double phi2)
Definition MathUtils.hh:698
double deltaPhi(double phi1, double phi2, bool sign=false)
Calculate the difference between two angles in radians.
Definition MathUtils.hh:668
double deltaEta(double eta1, double eta2, bool sign=false)
Definition MathUtils.hh:676
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:691
double mT(double pT1, double pT2, double dphi)
Definition MathUtils.hh:720
std::ostream & operator<<(std::ostream &os, const AnalysisInfo &ai)
Stream an AnalysisInfo as a text description.
Definition AnalysisInfo.hh:362
double invariant(const FourVector &lv)
Definition Vector4.hh:292
constexpr std::enable_if< std::is_arithmetic< NUM >::value, int >::type sign(NUM val)
Find the sign of a number.
Definition MathUtils.hh:265
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:684
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.
std::enable_if< std::is_arithmetic< NUM >::value, NUM >::type sqr(NUM a)
Named number-type squaring operation.
Definition MathUtils.hh:219
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:24
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:703