rivet is hosted by Hepforge, IPPP Durham
Rivet 3.1.6
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
10namespace Rivet {
11
12
13 class FourVector;
14 typedef FourVector Vector4;
15 typedef FourVector V4;
16
17 class FourMomentum;
18 typedef FourMomentum P4;
19
20 class LorentzTransform;
21 FourVector transform(const LorentzTransform& lt, const FourVector& v4);
22
23
27 class FourVector : public Vector<4> {
28 friend FourVector multiply(const double a, const FourVector& v);
29 friend FourVector multiply(const FourVector& v, const double a);
30 friend FourVector add(const FourVector& a, const FourVector& b);
31 friend FourVector transform(const LorentzTransform& lt, const FourVector& v4);
32
33 public:
34
35 FourVector() : Vector<4>() { }
36
37 template<typename V4TYPE, typename std::enable_if<HasXYZT<V4TYPE>::value, int>::type DUMMY=0>
38 FourVector(const V4TYPE& other) {
39 this->setT(other.t());
40 this->setX(other.x());
41 this->setY(other.y());
42 this->setZ(other.z());
43 }
44
45 FourVector(const Vector<4>& other)
46 : Vector<4>(other) { }
47
48 FourVector(const double t, const double x, const double y, const double z) {
49 this->setT(t);
50 this->setX(x);
51 this->setY(y);
52 this->setZ(z);
53 }
54
55 virtual ~FourVector() { }
56
57 public:
58
59 double t() const { return get(0); }
60 double t2() const { return sqr(t()); }
61 FourVector& setT(const double t) { set(0, t); return *this; }
62
63 double x() const { return get(1); }
64 double x2() const { return sqr(x()); }
65 FourVector& setX(const double x) { set(1, x); return *this; }
66
67 double y() const { return get(2); }
68 double y2() const { return sqr(y()); }
69 FourVector& setY(const double y) { set(2, y); return *this; }
70
71 double z() const { return get(3); }
72 double z2() const { return sqr(z()); }
73 FourVector& setZ(const double z) { set(3, z); return *this; }
74
75 double invariant() const {
76 // Done this way for numerical precision
77 return (t() + z())*(t() - z()) - x()*x() - y()*y();
78 }
79
80 bool isNull() const {
81 return Rivet::isZero(invariant());
82 }
83
85 double angle(const FourVector& v) const {
86 return vector3().angle( v.vector3() );
87 }
89 double angle(const Vector3& v3) const {
90 return vector3().angle(v3);
91 }
92
96 double polarRadius2() const {
97 return vector3().polarRadius2();
98 }
100 double perp2() const {
101 return vector3().perp2();
102 }
104 double rho2() const {
105 return vector3().rho2();
106 }
107
109 double polarRadius() const {
110 return vector3().polarRadius();
111 }
113 double perp() const {
114 return vector3().perp();
115 }
117 double rho() const {
118 return vector3().rho();
119 }
120
123 return vector3().polarVec();
124 }
126 Vector3 perpVec() const {
127 return vector3().perpVec();
128 }
130 Vector3 rhoVec() const {
131 return vector3().rhoVec();
132 }
133
135 double azimuthalAngle(const PhiMapping mapping=ZERO_2PI) const {
136 return vector3().azimuthalAngle(mapping);
137 }
139 double phi(const PhiMapping mapping=ZERO_2PI) const {
140 return vector3().phi(mapping);
141 }
142
144 double polarAngle() const {
145 return vector3().polarAngle();
146 }
148 double theta() const {
149 return vector3().theta();
150 }
151
153 double pseudorapidity() const {
154 return vector3().pseudorapidity();
155 }
157 double eta() const {
158 return vector3().eta();
159 }
160
162 double abspseudorapidity() const { return fabs(eta()); }
164 double abseta() const { return fabs(eta()); }
165
167 Vector3 vector3() const {
168 return Vector3(get(1), get(2), get(3));
169 }
170
172 operator Vector3 () const { return vector3(); }
173
174
175 public:
176
178 double contract(const FourVector& v) const {
179 const double result = t()*v.t() - x()*v.x() - y()*v.y() - z()*v.z();
180 return result;
181 }
182
184 double dot(const FourVector& v) const {
185 return contract(v);
186 }
187
189 double operator * (const FourVector& v) const {
190 return contract(v);
191 }
192
195 _vec = multiply(a, *this)._vec;
196 return *this;
197 }
198
201 _vec = multiply(1.0/a, *this)._vec;
202 return *this;
203 }
204
207 _vec = add(*this, v)._vec;
208 return *this;
209 }
210
213 _vec = add(*this, -v)._vec;
214 return *this;
215 }
216
219 FourVector result;
220 result._vec = -_vec;
221 return result;
222 }
223
226 FourVector result = -*this;
227 result.setT(-result.t());
228 return result;
229 }
230
231 };
232
233
235 inline double contract(const FourVector& a, const FourVector& b) {
236 return a.contract(b);
237 }
238
240 inline double dot(const FourVector& a, const FourVector& b) {
241 return contract(a, b);
242 }
243
244 inline FourVector multiply(const double a, const FourVector& v) {
245 FourVector result;
246 result._vec = a * v._vec;
247 return result;
248 }
249
250 inline FourVector multiply(const FourVector& v, const double a) {
251 return multiply(a, v);
252 }
253
254 inline FourVector operator * (const double a, const FourVector& v) {
255 return multiply(a, v);
256 }
257
258 inline FourVector operator * (const FourVector& v, const double a) {
259 return multiply(a, v);
260 }
261
262 inline FourVector operator / (const FourVector& v, const double a) {
263 return multiply(1.0/a, v);
264 }
265
266 inline FourVector add(const FourVector& a, const FourVector& b) {
267 FourVector result;
268 result._vec = a._vec + b._vec;
269 return result;
270 }
271
272 inline FourVector operator+(const FourVector& a, const FourVector& b) {
273 return add(a, b);
274 }
275
276 inline FourVector operator-(const FourVector& a, const FourVector& b) {
277 return add(a, -b);
278 }
279
282 inline double invariant(const FourVector& lv) {
283 return lv.invariant();
284 }
285
287 inline double angle(const FourVector& a, const FourVector& b) {
288 return a.angle(b);
289 }
290
292 inline double angle(const Vector3& a, const FourVector& b) {
293 return angle( a, b.vector3() );
294 }
295
297 inline double angle(const FourVector& a, const Vector3& b) {
298 return a.angle(b);
299 }
300
301
303
304
306 class FourMomentum : public FourVector {
307 friend FourMomentum multiply(const double a, const FourMomentum& v);
308 friend FourMomentum multiply(const FourMomentum& v, const double a);
309 friend FourMomentum add(const FourMomentum& a, const FourMomentum& b);
310 friend FourMomentum transform(const LorentzTransform& lt, const FourMomentum& v4);
311
312 public:
313 FourMomentum() { }
314
315 template<typename V4TYPE, typename std::enable_if<HasXYZT<V4TYPE>::value, int>::type DUMMY=0>
316 FourMomentum(const V4TYPE& other) {
317 this->setE(other.t());
318 this->setPx(other.x());
319 this->setPy(other.y());
320 this->setPz(other.z());
321 }
322
323 FourMomentum(const Vector<4>& other)
324 : FourVector(other) { }
325
326 FourMomentum(const double E, const double px, const double py, const double pz) {
327 this->setE(E);
328 this->setPx(px);
329 this->setPy(py);
330 this->setPz(pz);
331 }
332
333 ~FourMomentum() {}
334
335 public:
336
337
339
340
342 FourMomentum& setE(double E) {
343 setT(E);
344 return *this;
345 }
346
349 setX(px);
350 return *this;
351 }
352
355 setY(py);
356 return *this;
357 }
358
361 setZ(pz);
362 return *this;
363 }
364
365
367 FourMomentum& setPE(double px, double py, double pz, double E) {
368 if (E < 0)
369 throw std::invalid_argument("Negative energy given as argument: " + to_str(E));
370 setPx(px); setPy(py); setPz(pz); setE(E);
371 return *this;
372 }
374 FourMomentum& setXYZE(double px, double py, double pz, double E) {
375 return setPE(px, py, pz, E);
376 }
377 // /// Near-alias with switched arg order
378 // FourMomentum& setEP(double E, double px, double py, double pz) {
379 // return setPE(px, py, pz, E);
380 // }
381 // /// Alias for setEP
382 // FourMomentum& setEXYZ(double E, double px, double py, double pz) {
383 // return setEP(E, px, py, pz);
384 // }
385
386
388 FourMomentum& setPM(double px, double py, double pz, double mass) {
389 if (mass < 0)
390 throw std::invalid_argument("Negative mass given as argument: " + to_str(mass));
391 const double E = sqrt( sqr(mass) + sqr(px) + sqr(py) + sqr(pz) );
392 // setPx(px); setPy(py); setPz(pz); setE(E);
393 return setPE(px, py, pz, E);
394 }
396 FourMomentum& setXYZM(double px, double py, double pz, double mass) {
397 return setPM(px, py, pz, mass);
398 }
399
400
405 FourMomentum& setEtaPhiME(double eta, double phi, double mass, double E) {
406 if (mass < 0)
407 throw std::invalid_argument("Negative mass given as argument");
408 if (E < 0)
409 throw std::invalid_argument("Negative energy given as argument");
410 const double theta = 2 * atan(exp(-eta));
411 if (theta < 0 || theta > M_PI)
412 throw std::domain_error("Polar angle outside 0..pi in calculation");
414 return *this;
415 }
416
421 FourMomentum& setEtaPhiMPt(double eta, double phi, double mass, double pt) {
422 if (mass < 0)
423 throw std::invalid_argument("Negative mass given as argument");
424 if (pt < 0)
425 throw std::invalid_argument("Negative transverse momentum given as argument");
426 const double theta = 2 * atan(exp(-eta));
427 if (theta < 0 || theta > M_PI)
428 throw std::domain_error("Polar angle outside 0..pi in calculation");
429 const double p = pt / sin(theta);
430 const double E = sqrt( sqr(p) + sqr(mass) );
432 return *this;
433 }
434
443 FourMomentum& setRapPhiME(double y, double phi, double mass, double E) {
444 if (mass < 0)
445 throw std::invalid_argument("Negative mass given as argument");
446 if (E < 0)
447 throw std::invalid_argument("Negative energy given as argument");
448 const double sqrt_pt2_m2 = E / cosh(y);
449 const double pt = sqrt( sqr(sqrt_pt2_m2) - sqr(mass) );
450 if (pt < 0)
451 throw std::domain_error("Negative transverse momentum in calculation");
452 const double pz = sqrt_pt2_m2 * sinh(y);
453 const double px = pt * cos(phi);
454 const double py = pt * sin(phi);
455 setPE(px, py, pz, E);
456 return *this;
457 }
458
463 FourMomentum& setRapPhiMPt(double y, double phi, double mass, double pt) {
464 if (mass < 0)
465 throw std::invalid_argument("Negative mass given as argument");
466 if (pt < 0)
467 throw std::invalid_argument("Negative transverse mass given as argument");
468 const double E = sqrt( sqr(pt) + sqr(mass) ) * cosh(y);
469 if (E < 0)
470 throw std::domain_error("Negative energy in calculation");
471 setRapPhiME(y, phi, mass, E);
472 return *this;
473 }
474
480 FourMomentum& setThetaPhiME(double theta, double phi, double mass, double E) {
481 if (theta < 0 || theta > M_PI)
482 throw std::invalid_argument("Polar angle outside 0..pi given as argument");
483 if (mass < 0)
484 throw std::invalid_argument("Negative mass given as argument");
485 if (E < 0)
486 throw std::invalid_argument("Negative energy given as argument");
487 const double p = sqrt( sqr(E) - sqr(mass) );
488 const double pz = p * cos(theta);
489 const double pt = p * sin(theta);
490 if (pt < 0)
491 throw std::invalid_argument("Negative transverse momentum in calculation");
492 const double px = pt * cos(phi);
493 const double py = pt * sin(phi);
494 setPE(px, py, pz, E);
495 return *this;
496 }
497
503 FourMomentum& setThetaPhiMPt(double theta, double phi, double mass, double pt) {
504 if (theta < 0 || theta > M_PI)
505 throw std::invalid_argument("Polar angle outside 0..pi given as argument");
506 if (mass < 0)
507 throw std::invalid_argument("Negative mass given as argument");
508 if (pt < 0)
509 throw std::invalid_argument("Negative transverse momentum given as argument");
510 const double p = pt / sin(theta);
511 const double px = pt * cos(phi);
512 const double py = pt * sin(phi);
513 const double pz = p * cos(theta);
514 const double E = sqrt( sqr(p) + sqr(mass) );
515 setPE(px, py, pz, E);
516 return *this;
517 }
518
522 FourMomentum& setPtPhiME(double pt, double phi, double mass, double E) {
523 if (pt < 0)
524 throw std::invalid_argument("Negative transverse momentum given as argument");
525 if (mass < 0)
526 throw std::invalid_argument("Negative mass given as argument");
527 if (E < 0)
528 throw std::invalid_argument("Negative energy given as argument");
529 const double px = pt * cos(phi);
530 const double py = pt * sin(phi);
531 const double pz = sqrt(sqr(E) - sqr(mass) - sqr(pt));
532 setPE(px, py, pz, E);
533 return *this;
534 }
535
537
538
540
541
543 double E() const { return t(); }
545 double E2() const { return t2(); }
546
548 double px() const { return x(); }
550 double px2() const { return x2(); }
551
553 double py() const { return y(); }
555 double py2() const { return y2(); }
556
558 double pz() const { return z(); }
560 double pz2() const { return z2(); }
561
562
566 double mass() const {
567 // assert(Rivet::isZero(mass2()) || mass2() > 0);
568 // if (Rivet::isZero(mass2())) {
569 // return 0.0;
570 // } else {
571 // return sqrt(mass2());
572 // }
573 return sign(mass2()) * sqrt(fabs(mass2()));
574 }
575
577 double mass2() const {
578 return invariant();
579 }
580
581
583 Vector3 p3() const { return vector3(); }
584
586 double p() const {
587 return p3().mod();
588 }
589
591 double p2() const {
592 return p3().mod2();
593 }
594
595
597 double rapidity() const {
598 return 0.5 * std::log( (E() + pz()) / (E() - pz()) );
599 }
601 double rap() const {
602 return rapidity();
603 }
604
606 double absrapidity() const {
607 return fabs(rapidity());
608 }
610 double absrap() const {
611 return fabs(rap());
612 }
613
615 Vector3 pTvec() const {
616 return p3().polarVec();
617 }
619 Vector3 ptvec() const {
620 return pTvec();
621 }
622
624 double pT2() const {
625 return vector3().polarRadius2();
626 }
628 double pt2() const {
629 return vector3().polarRadius2();
630 }
631
633 double pT() const {
634 return sqrt(pT2());
635 }
637 double pt() const {
638 return sqrt(pT2());
639 }
640
642 double Et2() const {
643 return Et() * Et();
644 }
646 double Et() const {
647 return E() * sin(polarAngle());
648 }
649
651
652
654
655
658 double gamma() const {
659 return sqrt(E2()/mass2());
660 }
661
665 return gamma() * p3().unit();
666 }
667
670 double beta() const {
671 return p()/E();
672 }
673
676 Vector3 betaVec() const {
677 // return Vector3(px()/E(), py()/E(), pz()/E());
678 return p3()/E();
679 }
680
682
683
685
686
688
691 struct byEAscending {
692 bool operator()(const FourMomentum& left, const FourMomentum& right) const{
693 const double pt2left = left.E();
694 const double pt2right = right.E();
695 return pt2left < pt2right;
696 }
697
698 bool operator()(const FourMomentum* left, const FourMomentum* right) const{
699 return (*this)(*left, *right);
700 }
701 };
702
703
706 struct byEDescending {
707 bool operator()(const FourMomentum& left, const FourMomentum& right) const{
708 return byEAscending()(right, left);
709 }
710
711 bool operator()(const FourMomentum* left, const FourVector* right) const{
712 return (*this)(*left, *right);
713 }
714 };
715
717
718
720
721
723
724
727 _vec = multiply(a, *this)._vec;
728 return *this;
729 }
730
733 _vec = multiply(1.0/a, *this)._vec;
734 return *this;
735 }
736
739 _vec = add(*this, v)._vec;
740 return *this;
741 }
742
745 _vec = add(*this, -v)._vec;
746 return *this;
747 }
748
751 FourMomentum result;
752 result._vec = -_vec;
753 return result;
754 }
755
758 FourMomentum result = -*this;
759 result.setE(-result.E());
760 return result;
761 }
762
764
765
767
768
770
771
773 static FourMomentum mkXYZE(double px, double py, double pz, double E) {
774 return FourMomentum().setPE(px, py, pz, E);
775 }
776
778 static FourMomentum mkXYZM(double px, double py, double pz, double mass) {
779 return FourMomentum().setPM(px, py, pz, mass);
780 }
781
783 static FourMomentum mkEtaPhiME(double eta, double phi, double mass, double E) {
784 return FourMomentum().setEtaPhiME(eta, phi, mass, E);
785 }
786
788 static FourMomentum mkEtaPhiMPt(double eta, double phi, double mass, double pt) {
789 return FourMomentum().setEtaPhiMPt(eta, phi, mass, pt);
790 }
791
793 static FourMomentum mkRapPhiME(double y, double phi, double mass, double E) {
794 return FourMomentum().setRapPhiME(y, phi, mass, E);
795 }
796
798 static FourMomentum mkRapPhiMPt(double y, double phi, double mass, double pt) {
799 return FourMomentum().setRapPhiMPt(y, phi, mass, pt);
800 }
801
803 static FourMomentum mkThetaPhiME(double theta, double phi, double mass, double E) {
805 }
806
808 static FourMomentum mkThetaPhiMPt(double theta, double phi, double mass, double pt) {
810 }
811
813 static FourMomentum mkPtPhiME(double pt, double phi, double mass, double E) {
814 return FourMomentum().setPtPhiME(pt, phi, mass, E);
815 }
816
818
819
820 };
821
822
823 inline FourMomentum multiply(const double a, const FourMomentum& v) {
824 FourMomentum result;
825 result._vec = a * v._vec;
826 return result;
827 }
828
829 inline FourMomentum multiply(const FourMomentum& v, const double a) {
830 return multiply(a, v);
831 }
832
833 inline FourMomentum operator*(const double a, const FourMomentum& v) {
834 return multiply(a, v);
835 }
836
837 inline FourMomentum operator*(const FourMomentum& v, const double a) {
838 return multiply(a, v);
839 }
840
841 inline FourMomentum operator/(const FourMomentum& v, const double a) {
842 return multiply(1.0/a, v);
843 }
844
845 inline FourMomentum add(const FourMomentum& a, const FourMomentum& b) {
846 FourMomentum result;
847 result._vec = a._vec + b._vec;
848 return result;
849 }
850
851 inline FourMomentum operator+(const FourMomentum& a, const FourMomentum& b) {
852 return add(a, b);
853 }
854
855 inline FourMomentum operator-(const FourMomentum& a, const FourMomentum& b) {
856 return add(a, -b);
857 }
858
859
861
862
864
865
874 inline double deltaR2(const FourVector& a, const FourVector& b,
875 RapScheme scheme=PSEUDORAPIDITY) {
876 switch (scheme) {
877 case PSEUDORAPIDITY :
878 return deltaR2(a.vector3(), b.vector3());
879 case RAPIDITY:
880 {
881 const FourMomentum* ma = dynamic_cast<const FourMomentum*>(&a);
882 const FourMomentum* mb = dynamic_cast<const FourMomentum*>(&b);
883 if (!ma || !mb) {
884 string err = "deltaR with scheme RAPIDITY can only be called with FourMomentum objects, not FourVectors";
885 throw std::runtime_error(err);
886 }
887 return deltaR2(*ma, *mb, scheme);
888 }
889 default:
890 throw std::runtime_error("The specified deltaR scheme is not yet implemented");
891 }
892 }
893
902 inline double deltaR(const FourVector& a, const FourVector& b,
903 RapScheme scheme=PSEUDORAPIDITY) {
904 return sqrt(deltaR2(a, b, scheme));
905 }
906
907
908
915 inline double deltaR2(const FourVector& v,
916 double eta2, double phi2,
917 RapScheme scheme=PSEUDORAPIDITY) {
918 switch (scheme) {
919 case PSEUDORAPIDITY :
920 return deltaR2(v.vector3(), eta2, phi2);
921 case RAPIDITY:
922 {
923 const FourMomentum* mv = dynamic_cast<const FourMomentum*>(&v);
924 if (!mv) {
925 string err = "deltaR with scheme RAPIDITY can only be called with FourMomentum objects, not FourVectors";
926 throw std::runtime_error(err);
927 }
928 return deltaR2(*mv, eta2, phi2, scheme);
929 }
930 default:
931 throw std::runtime_error("The specified deltaR scheme is not yet implemented");
932 }
933 }
934
941 inline double deltaR(const FourVector& v,
942 double eta2, double phi2,
943 RapScheme scheme=PSEUDORAPIDITY) {
944 return sqrt(deltaR2(v, eta2, phi2, scheme));
945 }
946
947
954 inline double deltaR2(double eta1, double phi1,
955 const FourVector& v,
956 RapScheme scheme=PSEUDORAPIDITY) {
957 switch (scheme) {
958 case PSEUDORAPIDITY :
959 return deltaR2(eta1, phi1, v.vector3());
960 case RAPIDITY:
961 {
962 const FourMomentum* mv = dynamic_cast<const FourMomentum*>(&v);
963 if (!mv) {
964 string err = "deltaR with scheme RAPIDITY can only be called with FourMomentum objects, not FourVectors";
965 throw std::runtime_error(err);
966 }
967 return deltaR2(eta1, phi1, *mv, scheme);
968 }
969 default:
970 throw std::runtime_error("The specified deltaR scheme is not yet implemented");
971 }
972 }
973
980 inline double deltaR(double eta1, double phi1,
981 const FourVector& v,
982 RapScheme scheme=PSEUDORAPIDITY) {
983 return sqrt(deltaR2(eta1, phi1, v, scheme));
984 }
985
986
993 inline double deltaR2(const FourMomentum& a, const FourMomentum& b,
994 RapScheme scheme=PSEUDORAPIDITY) {
995 switch (scheme) {
996 case PSEUDORAPIDITY:
997 return deltaR2(a.vector3(), b.vector3());
998 case RAPIDITY:
999 return deltaR2(a.rapidity(), a.azimuthalAngle(), b.rapidity(), b.azimuthalAngle());
1000 default:
1001 throw std::runtime_error("The specified deltaR scheme is not yet implemented");
1002 }
1003 }
1004
1011 inline double deltaR(const FourMomentum& a, const FourMomentum& b,
1012 RapScheme scheme=PSEUDORAPIDITY) {
1013 return sqrt(deltaR2(a, b, scheme));
1014 }
1015
1016
1022 inline double deltaR2(const FourMomentum& v,
1023 double eta2, double phi2,
1024 RapScheme scheme=PSEUDORAPIDITY) {
1025 switch (scheme) {
1026 case PSEUDORAPIDITY:
1027 return deltaR2(v.vector3(), eta2, phi2);
1028 case RAPIDITY:
1029 return deltaR2(v.rapidity(), v.azimuthalAngle(), eta2, phi2);
1030 default:
1031 throw std::runtime_error("The specified deltaR scheme is not yet implemented");
1032 }
1033 }
1034
1040 inline double deltaR(const FourMomentum& v,
1041 double eta2, double phi2,
1042 RapScheme scheme=PSEUDORAPIDITY) {
1043 return sqrt(deltaR2(v, eta2, phi2, scheme));
1044 }
1045
1046
1052 inline double deltaR2(double eta1, double phi1,
1053 const FourMomentum& v,
1054 RapScheme scheme=PSEUDORAPIDITY) {
1055 switch (scheme) {
1056 case PSEUDORAPIDITY:
1057 return deltaR2(eta1, phi1, v.vector3());
1058 case RAPIDITY:
1059 return deltaR2(eta1, phi1, v.rapidity(), v.azimuthalAngle());
1060 default:
1061 throw std::runtime_error("The specified deltaR scheme is not yet implemented");
1062 }
1063 }
1064
1070 inline double deltaR(double eta1, double phi1,
1071 const FourMomentum& v,
1072 RapScheme scheme=PSEUDORAPIDITY) {
1073 return sqrt(deltaR2(eta1, phi1, v, scheme));
1074 }
1075
1076
1082 inline double deltaR2(const FourMomentum& a, const FourVector& b,
1083 RapScheme scheme=PSEUDORAPIDITY) {
1084 switch (scheme) {
1085 case PSEUDORAPIDITY:
1086 return deltaR2(a.vector3(), b.vector3());
1087 case RAPIDITY:
1089 default:
1090 throw std::runtime_error("The specified deltaR scheme is not yet implemented");
1091 }
1092 }
1093
1099 inline double deltaR(const FourMomentum& a, const FourVector& b,
1100 RapScheme scheme=PSEUDORAPIDITY) {
1101 return sqrt(deltaR2(a, b, scheme));
1102 }
1103
1104
1110 inline double deltaR2(const FourVector& a, const FourMomentum& b,
1111 RapScheme scheme=PSEUDORAPIDITY) {
1112 return deltaR2(b, a, scheme); //< note reversed args
1113 }
1114
1120 inline double deltaR(const FourVector& a, const FourMomentum& b,
1121 RapScheme scheme=PSEUDORAPIDITY) {
1122 return deltaR(b, a, scheme); //< note reversed args
1123 }
1124
1125
1128 inline double deltaR2(const FourMomentum& a, const Vector3& b) {
1129 return deltaR2(a.vector3(), b);
1130 }
1131
1134 inline double deltaR(const FourMomentum& a, const Vector3& b) {
1135 return deltaR(a.vector3(), b);
1136 }
1137
1140 inline double deltaR2(const Vector3& a, const FourMomentum& b) {
1141 return deltaR2(a, b.vector3());
1142 }
1143
1146 inline double deltaR(const Vector3& a, const FourMomentum& b) {
1147 return deltaR(a, b.vector3());
1148 }
1149
1152 inline double deltaR2(const FourVector& a, const Vector3& b) {
1153 return deltaR2(a.vector3(), b);
1154 }
1155
1158 inline double deltaR(const FourVector& a, const Vector3& b) {
1159 return deltaR(a.vector3(), b);
1160 }
1161
1164 inline double deltaR2(const Vector3& a, const FourVector& b) {
1165 return deltaR2(a, b.vector3());
1166 }
1167
1170 inline double deltaR(const Vector3& a, const FourVector& b) {
1171 return deltaR(a, b.vector3());
1172 }
1173
1175
1176
1178
1179
1181
1182
1184 inline double deltaPhi(const FourMomentum& a, const FourMomentum& b, bool sign=false) {
1185 return deltaPhi(a.vector3(), b.vector3(), sign);
1186 }
1187
1189 inline double deltaPhi(const FourMomentum& v, double phi2, bool sign=false) {
1190 return deltaPhi(v.vector3(), phi2, sign);
1191 }
1192
1194 inline double deltaPhi(double phi1, const FourMomentum& v, bool sign=false) {
1195 return deltaPhi(phi1, v.vector3(), sign);
1196 }
1197
1199 inline double deltaPhi(const FourVector& a, const FourVector& b, bool sign=false) {
1200 return deltaPhi(a.vector3(), b.vector3(), sign);
1201 }
1202
1204 inline double deltaPhi(const FourVector& v, double phi2, bool sign=false) {
1205 return deltaPhi(v.vector3(), phi2, sign);
1206 }
1207
1209 inline double deltaPhi(double phi1, const FourVector& v, bool sign=false) {
1210 return deltaPhi(phi1, v.vector3(), sign);
1211 }
1212
1214 inline double deltaPhi(const FourVector& a, const FourMomentum& b, bool sign=false) {
1215 return deltaPhi(a.vector3(), b.vector3(), sign);
1216 }
1217
1219 inline double deltaPhi(const FourMomentum& a, const FourVector& b, bool sign=false) {
1220 return deltaPhi(a.vector3(), b.vector3(), sign);
1221 }
1222
1224 inline double deltaPhi(const FourVector& a, const Vector3& b, bool sign=false) {
1225 return deltaPhi(a.vector3(), b, sign);
1226 }
1227
1229 inline double deltaPhi(const Vector3& a, const FourVector& b, bool sign=false) {
1230 return deltaPhi(a, b.vector3(), sign);
1231 }
1232
1234 inline double deltaPhi(const FourMomentum& a, const Vector3& b, bool sign=false) {
1235 return deltaPhi(a.vector3(), b, sign);
1236 }
1237
1239 inline double deltaPhi(const Vector3& a, const FourMomentum& b, bool sign=false) {
1240 return deltaPhi(a, b.vector3(), sign);
1241 }
1242
1244
1245
1247
1248
1250
1251
1253 inline double deltaEta(const FourMomentum& a, const FourMomentum& b, bool sign=false) {
1254 return deltaEta(a.vector3(), b.vector3(), sign);
1255 }
1256
1258 inline double deltaEta(const FourMomentum& v, double eta2, bool sign=false) {
1259 return deltaEta(v.vector3(), eta2, sign);
1260 }
1261
1263 inline double deltaEta(double eta1, const FourMomentum& v, bool sign=false) {
1264 return deltaEta(eta1, v.vector3(), sign);
1265 }
1266
1268 inline double deltaEta(const FourVector& a, const FourVector& b, bool sign=false) {
1269 return deltaEta(a.vector3(), b.vector3(), sign);
1270 }
1271
1273 inline double deltaEta(const FourVector& v, double eta2, bool sign=false) {
1274 return deltaEta(v.vector3(), eta2, sign);
1275 }
1276
1278 inline double deltaEta(double eta1, const FourVector& v, bool sign=false) {
1279 return deltaEta(eta1, v.vector3(), sign);
1280 }
1281
1283 inline double deltaEta(const FourVector& a, const FourMomentum& b, bool sign=false) {
1284 return deltaEta(a.vector3(), b.vector3(), sign);
1285 }
1286
1288 inline double deltaEta(const FourMomentum& a, const FourVector& b, bool sign=false) {
1289 return deltaEta(a.vector3(), b.vector3(), sign);
1290 }
1291
1293 inline double deltaEta(const FourVector& a, const Vector3& b, bool sign=false) {
1294 return deltaEta(a.vector3(), b, sign);
1295 }
1296
1298 inline double deltaEta(const Vector3& a, const FourVector& b, bool sign=false) {
1299 return deltaEta(a, b.vector3(), sign);
1300 }
1301
1303 inline double deltaEta(const FourMomentum& a, const Vector3& b, bool sign=false) {
1304 return deltaEta(a.vector3(), b, sign);
1305 }
1306
1308 inline double deltaEta(const Vector3& a, const FourMomentum& b, bool sign=false) {
1309 return deltaEta(a, b.vector3(), sign);
1310 }
1311
1313
1314
1316
1317
1319 inline double deltaRap(const FourMomentum& a, const FourMomentum& b, bool sign=false) {
1320 return deltaRap(a.rapidity(), b.rapidity(), sign);
1321 }
1322
1324 inline double deltaRap(const FourMomentum& v, double y2, bool sign=false) {
1325 return deltaRap(v.rapidity(), y2, sign);
1326 }
1327
1329 inline double deltaRap(double y1, const FourMomentum& v, bool sign=false) {
1330 return deltaRap(y1, v.rapidity(), sign);
1331 }
1332
1334
1335
1337
1338
1341
1344
1346 inline bool cmpMomByPt(const FourMomentum& a, const FourMomentum& b) {
1347 return a.pt() > b.pt();
1348 }
1350 inline bool cmpMomByAscPt(const FourMomentum& a, const FourMomentum& b) {
1351 return a.pt() < b.pt();
1352 }
1353
1355 inline bool cmpMomByP(const FourMomentum& a, const FourMomentum& b) {
1356 return a.vector3().mod() > b.vector3().mod();
1357 }
1359 inline bool cmpMomByAscP(const FourMomentum& a, const FourMomentum& b) {
1360 return a.vector3().mod() < b.vector3().mod();
1361 }
1362
1364 inline bool cmpMomByEt(const FourMomentum& a, const FourMomentum& b) {
1365 return a.Et() > b.Et();
1366 }
1368 inline bool cmpMomByAscEt(const FourMomentum& a, const FourMomentum& b) {
1369 return a.Et() < b.Et();
1370 }
1371
1373 inline bool cmpMomByE(const FourMomentum& a, const FourMomentum& b) {
1374 return a.E() > b.E();
1375 }
1377 inline bool cmpMomByAscE(const FourMomentum& a, const FourMomentum& b) {
1378 return a.E() < b.E();
1379 }
1380
1382 inline bool cmpMomByMass(const FourMomentum& a, const FourMomentum& b) {
1383 return a.mass() > b.mass();
1384 }
1386 inline bool cmpMomByAscMass(const FourMomentum& a, const FourMomentum& b) {
1387 return a.mass() < b.mass();
1388 }
1389
1391 inline bool cmpMomByEta(const FourMomentum& a, const FourMomentum& b) {
1392 return a.eta() < b.eta();
1393 }
1394
1396 inline bool cmpMomByDescEta(const FourMomentum& a, const FourMomentum& b) {
1397 return a.pseudorapidity() > b.pseudorapidity();
1398 }
1399
1401 inline bool cmpMomByAbsEta(const FourMomentum& a, const FourMomentum& b) {
1402 return fabs(a.eta()) < fabs(b.eta());
1403 }
1404
1406 inline bool cmpMomByDescAbsEta(const FourMomentum& a, const FourMomentum& b) {
1407 return fabs(a.eta()) > fabs(b.eta());
1408 }
1409
1411 inline bool cmpMomByRap(const FourMomentum& a, const FourMomentum& b) {
1412 return a.rapidity() < b.rapidity();
1413 }
1414
1416 inline bool cmpMomByDescRap(const FourMomentum& a, const FourMomentum& b) {
1417 return a.rapidity() > b.rapidity();
1418 }
1419
1421 inline bool cmpMomByAbsRap(const FourMomentum& a, const FourMomentum& b) {
1422 return fabs(a.rapidity()) < fabs(b.rapidity());
1423 }
1424
1426 inline bool cmpMomByDescAbsRap(const FourMomentum& a, const FourMomentum& b) {
1427 return fabs(a.rapidity()) > fabs(b.rapidity());
1428 }
1429
1431
1432
1434 template<typename MOMS, typename CMP>
1435 inline MOMS& isortBy(MOMS& pbs, const CMP& cmp) {
1436 std::sort(pbs.begin(), pbs.end(), cmp);
1437 return pbs;
1438 }
1440 template<typename MOMS, typename CMP>
1441 inline MOMS sortBy(const MOMS& pbs, const CMP& cmp) {
1442 MOMS rtn = pbs;
1443 std::sort(rtn.begin(), rtn.end(), cmp);
1444 return rtn;
1445 }
1446
1448 template<typename MOMS>
1449 inline MOMS& isortByPt(MOMS& pbs) {
1450 return isortBy(pbs, cmpMomByPt);
1451 }
1453 template<typename MOMS>
1454 inline MOMS sortByPt(const MOMS& pbs) {
1455 return sortBy(pbs, cmpMomByPt);
1456 }
1457
1459 template<typename MOMS>
1460 inline MOMS& isortByE(MOMS& pbs) {
1461 return isortBy(pbs, cmpMomByE);
1462 }
1464 template<typename MOMS>
1465 inline MOMS sortByE(const MOMS& pbs) {
1466 return sortBy(pbs, cmpMomByE);
1467 }
1468
1470 template<typename MOMS>
1471 inline MOMS& isortByEt(MOMS& pbs) {
1472 return isortBy(pbs, cmpMomByEt);
1473 }
1475 template<typename MOMS>
1476 inline MOMS sortByEt(const MOMS& pbs) {
1477 return sortBy(pbs, cmpMomByEt);
1478 }
1479
1481
1482
1485
1487 inline double mT(const FourMomentum& vis, const FourMomentum& invis) {
1488 return mT(vis.p3(), invis.p3());
1489 }
1490
1492 inline double mT(const FourMomentum& vis, const Vector3& invis) {
1493 return mT(vis.p3(), invis);
1494 }
1495
1497 inline double mT(const Vector3& vis, const FourMomentum& invis) {
1498 return mT(vis, invis.p3());
1499 }
1500
1502
1503
1505
1506
1509
1511 inline std::string toString(const FourVector& lv) {
1512 std::ostringstream out;
1513 out << "(" << (fabs(lv.t()) < 1E-30 ? 0.0 : lv.t())
1514 << "; " << (fabs(lv.x()) < 1E-30 ? 0.0 : lv.x())
1515 << ", " << (fabs(lv.y()) < 1E-30 ? 0.0 : lv.y())
1516 << ", " << (fabs(lv.z()) < 1E-30 ? 0.0 : lv.z())
1517 << ")";
1518 return out.str();
1519 }
1520
1522 inline std::ostream& operator<<(std::ostream& out, const FourVector& lv) {
1523 out << toString(lv);
1524 return out;
1525 }
1526
1528
1531 typedef std::vector<FourVector> FourVectors;
1532 typedef std::vector<FourMomentum> FourMomenta;
1534
1536
1537
1538}
1539
1540#endif
Specialized version of the FourVector with momentum/energy functionality.
Definition: Vector4.hh:306
Specialisation of VectorN to a general (non-momentum) Lorentz 4-vector.
Definition: Vector4.hh:27
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
const CONTAINER2 & transform(const CONTAINER1 &in, CONTAINER2 &out, const FN &f)
A single-container-arg version of std::transform, aka map.
Definition: Utils.hh:408
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:1435
MOMS & isortByE(MOMS &pbs)
Sort a container of momenta by E (decreasing) and return by reference for non-const inputs.
Definition: Vector4.hh:1460
MOMS & isortByEt(MOMS &pbs)
Sort a container of momenta by Et (decreasing) and return by reference for non-const inputs.
Definition: Vector4.hh:1471
bool cmpMomByDescEta(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by decreasing eta (pseudorapidity)
Definition: Vector4.hh:1396
bool cmpMomByMass(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by decreasing mass.
Definition: Vector4.hh:1382
bool cmpMomByP(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by decreasing 3-momentum magnitude |p|.
Definition: Vector4.hh:1355
bool cmpMomByAscEt(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by increasing transverse energy.
Definition: Vector4.hh:1368
bool cmpMomByRap(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by increasing rapidity.
Definition: Vector4.hh:1411
bool cmpMomByE(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by decreasing energy.
Definition: Vector4.hh:1373
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:1441
bool cmpMomByAbsRap(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by increasing absolute rapidity.
Definition: Vector4.hh:1421
MOMS sortByPt(const MOMS &pbs)
Sort a container of momenta by pT (decreasing) and return by value for const inputs.
Definition: Vector4.hh:1454
MOMS sortByEt(const MOMS &pbs)
Sort a container of momenta by Et (decreasing) and return by value for const inputs.
Definition: Vector4.hh:1476
bool cmpMomByEt(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by decreasing transverse energy.
Definition: Vector4.hh:1364
bool cmpMomByAscPt(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by increasing pT.
Definition: Vector4.hh:1350
bool cmpMomByAscP(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by increasing 3-momentum magnitude |p|.
Definition: Vector4.hh:1359
bool cmpMomByAscMass(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by increasing mass.
Definition: Vector4.hh:1386
MOMS sortByE(const MOMS &pbs)
Sort a container of momenta by E (decreasing) and return by value for const inputs.
Definition: Vector4.hh:1465
bool cmpMomByPt(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by decreasing pT.
Definition: Vector4.hh:1346
bool cmpMomByAbsEta(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by increasing absolute eta (pseudorapidity)
Definition: Vector4.hh:1401
bool cmpMomByDescAbsEta(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by increasing absolute eta (pseudorapidity)
Definition: Vector4.hh:1406
bool cmpMomByEta(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by increasing eta (pseudorapidity)
Definition: Vector4.hh:1391
bool cmpMomByAscE(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by increasing energy.
Definition: Vector4.hh:1377
bool cmpMomByDescAbsRap(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by decreasing absolute rapidity.
Definition: Vector4.hh:1426
MOMS & isortByPt(MOMS &pbs)
Sort a container of momenta by pT (decreasing) and return by reference for non-const inputs.
Definition: Vector4.hh:1449
std::vector< FourVector > FourVectors
Definition: Vector4.hh:1531
bool cmpMomByDescRap(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by decreasing rapidity.
Definition: Vector4.hh:1416
FourVector operator-() const
Multiply all components (space and time) by -1.
Definition: Vector4.hh:218
double p2() const
Get the modulus-squared of the 3-momentum.
Definition: Vector4.hh:591
double E2() const
Get energy-squared .
Definition: Vector4.hh:545
double mass() const
Get the mass (the Lorentz self-invariant).
Definition: Vector4.hh:566
double pz2() const
Get z-squared .
Definition: Vector4.hh:560
double rapidity() const
Calculate the rapidity.
Definition: Vector4.hh:597
FourMomentum reverse() const
Multiply space components only by -1.
Definition: Vector4.hh:757
double perp2() const
Synonym for polarRadius2.
Definition: Vector4.hh:100
double px() const
Get x-component of momentum .
Definition: Vector4.hh:548
FourMomentum & setPtPhiME(double pt, double phi, double mass, double E)
Definition: Vector4.hh:522
FourMomentum & setRapPhiMPt(double y, double phi, double mass, double pt)
Definition: Vector4.hh:463
Vector3 p3() const
Get 3-momentum part, .
Definition: Vector4.hh:583
FourVector & operator/=(double a)
Divide by a scalar.
Definition: Vector4.hh:200
double rho() const
Synonym for polarRadius.
Definition: Vector4.hh:117
FourMomentum & setPz(double pz)
Set z-component of momentum .
Definition: Vector4.hh:360
double p() const
Get the modulus of the 3-momentum.
Definition: Vector4.hh:586
FourMomentum & setXYZE(double px, double py, double pz, double E)
Alias for setPE.
Definition: Vector4.hh:374
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:803
Vector3 perpVec() const
Synonym for polarVec.
Definition: Vector4.hh:126
Vector3 vector3() const
Get the spatial part of the 4-vector as a 3-vector.
Definition: Vector4.hh:167
double angle(const FourVector &v) const
Angle between this vector and another.
Definition: Vector4.hh:85
Vector3 gammaVec() const
Definition: Vector4.hh:664
double pt() const
Calculate the transverse momentum .
Definition: Vector4.hh:637
double phi(const PhiMapping mapping=ZERO_2PI) const
Synonym for azimuthalAngle.
Definition: Vector4.hh:139
double contract(const FourVector &v) const
Contract two 4-vectors, with metric signature (+ - - -).
Definition: Vector4.hh:178
FourMomentum & setPE(double px, double py, double pz, double E)
Set the p coordinates and energy simultaneously.
Definition: Vector4.hh:367
double py() const
Get y-component of momentum .
Definition: Vector4.hh:553
double abseta() const
Get the directly (alias).
Definition: Vector4.hh:164
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:793
double beta() const
Definition: Vector4.hh:670
double Et() const
Calculate the transverse energy .
Definition: Vector4.hh:646
double Et2() const
Calculate the transverse energy .
Definition: Vector4.hh:642
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:135
static FourMomentum mkXYZE(double px, double py, double pz, double E)
Make a vector from (px,py,pz,E) coordinates.
Definition: Vector4.hh:773
FourMomentum & setPy(double py)
Set y-component of momentum .
Definition: Vector4.hh:354
double rap() const
Alias for rapidity.
Definition: Vector4.hh:601
FourMomentum & setPx(double px)
Set x-component of momentum .
Definition: Vector4.hh:348
double pt2() const
Calculate the squared transverse momentum .
Definition: Vector4.hh:628
double dot(const FourVector &v) const
Contract two 4-vectors, with metric signature (+ - - -).
Definition: Vector4.hh:184
double perp() const
Synonym for polarRadius.
Definition: Vector4.hh:113
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:813
Vector3 pTvec() const
Calculate the transverse momentum vector .
Definition: Vector4.hh:615
FourVector & operator*=(double a)
Multiply by a scalar.
Definition: Vector4.hh:194
double operator*(const FourVector &v) const
Contract two 4-vectors, with metric signature (+ - - -).
Definition: Vector4.hh:189
double eta() const
Synonym for pseudorapidity.
Definition: Vector4.hh:157
FourMomentum & setRapPhiME(double y, double phi, double mass, double E)
Definition: Vector4.hh:443
double theta() const
Synonym for polarAngle.
Definition: Vector4.hh:148
Vector3 ptvec() const
Synonym for pTvec.
Definition: Vector4.hh:619
double mass2() const
Get the squared mass (the Lorentz self-invariant).
Definition: Vector4.hh:577
FourMomentum & setEtaPhiMPt(double eta, double phi, double mass, double pt)
Definition: Vector4.hh:421
double pz() const
Get z-component of momentum .
Definition: Vector4.hh:558
double pseudorapidity() const
Pseudorapidity (defined purely by the 3-vector components)
Definition: Vector4.hh:153
double py2() const
Get y-squared .
Definition: Vector4.hh:555
double absrap() const
Absolute rapidity.
Definition: Vector4.hh:610
double pT2() const
Calculate the squared transverse momentum .
Definition: Vector4.hh:624
double polarRadius() const
Magnitude of projection of 3-vector on to the plane.
Definition: Vector4.hh:109
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:783
double px2() const
Get x-squared .
Definition: Vector4.hh:550
Vector3 betaVec() const
Definition: Vector4.hh:676
Vector3 rhoVec() const
Synonym for polarVec.
Definition: Vector4.hh:130
double angle(const Vector3 &v3) const
Angle between this vector and another (3-vector)
Definition: Vector4.hh:89
Vector3 polarVec() const
Projection of 3-vector on to the plane.
Definition: Vector4.hh:122
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:96
FourVector reverse() const
Multiply space components only by -1.
Definition: Vector4.hh:225
FourMomentum & operator-=(const FourMomentum &v)
Subtract from this 4-vector. NB time as well as space components are subtracted.
Definition: Vector4.hh:744
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:778
FourMomentum & setThetaPhiME(double theta, double phi, double mass, double E)
Definition: Vector4.hh:480
double E() const
Get energy (time component of momentum).
Definition: Vector4.hh:543
FourMomentum & setXYZM(double px, double py, double pz, double mass)
Alias for setPM.
Definition: Vector4.hh:396
FourMomentum & operator+=(const FourMomentum &v)
Add to this 4-vector. NB time as well as space components are added.
Definition: Vector4.hh:738
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:788
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:798
double pT() const
Calculate the transverse momentum .
Definition: Vector4.hh:633
double polarAngle() const
Angle subtended by the 3-vector and the z-axis.
Definition: Vector4.hh:144
FourVector & operator-=(const FourVector &v)
Subtract from this 4-vector. NB time as well as space components are subtracted.
Definition: Vector4.hh:212
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:808
double abspseudorapidity() const
Get the directly.
Definition: Vector4.hh:162
FourVector & operator+=(const FourVector &v)
Add to this 4-vector.
Definition: Vector4.hh:206
FourMomentum & setThetaPhiMPt(double theta, double phi, double mass, double pt)
Definition: Vector4.hh:503
FourMomentum operator-() const
Multiply all components (time and space) by -1.
Definition: Vector4.hh:750
double absrapidity() const
Absolute rapidity.
Definition: Vector4.hh:606
FourMomentum & setEtaPhiME(double eta, double phi, double mass, double E)
Definition: Vector4.hh:405
FourMomentum & operator/=(double a)
Divide by a scalar.
Definition: Vector4.hh:732
double rho2() const
Synonym for polarRadius2.
Definition: Vector4.hh:104
FourMomentum & setE(double E)
Set energy (time component of momentum).
Definition: Vector4.hh:342
double gamma() const
Definition: Vector4.hh:658
FourMomentum & operator*=(double a)
Multiply by a scalar.
Definition: Vector4.hh:726
FourMomentum & setPM(double px, double py, double pz, double mass)
Set the p coordinates and mass simultaneously.
Definition: Vector4.hh:388
string to_str(const T &x)
Convert any object to a string.
Definition: Utils.hh:75
Definition: MC_Cent_pPb.hh:10
double deltaR(double rap1, double phi1, double rap2, double phi2)
Definition: MathUtils.hh:659
double deltaPhi(double phi1, double phi2, bool sign=false)
Calculate the difference between two angles in radians.
Definition: MathUtils.hh:629
double deltaEta(double eta1, double eta2, bool sign=false)
Definition: MathUtils.hh:637
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:652
double mT(double pT1, double pT2, double dphi)
Definition: MathUtils.hh:681
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:282
double contract(const FourVector &a, const FourVector &b)
Contract two 4-vectors, with metric signature (+ - - -).
Definition: Vector4.hh:235
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:645
std::enable_if< std::is_arithmetic< NUM >::value, int >::type sign(NUM val)
Find the sign of a number.
Definition: MathUtils.hh:266
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:664