Rivet API documentation

Rivet 4.1.3
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
40 template<typename V4TYPE, typename std::enable_if<HasXYZT<V4TYPE>::value, int>::type DUMMY=0>
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
48 FourVector(const Vector<4>& other)
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 FourVector(double t, const Vector3& v3) {
59 this->setT(t);
60 this->setX(v3.x());
61 this->setY(v3.y());
62 this->setZ(v3.z());
63 }
64
65 virtual ~FourVector() { }
66
71 operator fastjet::PseudoJet () const;
72
73
74 public:
75
76 double t() const { return get(0); }
77 double t2() const { return sqr(t()); }
78 FourVector& setT(const double t) { set(0, t); return *this; }
79
80 double x() const { return get(1); }
81 double x2() const { return sqr(x()); }
82 FourVector& setX(const double x) { set(1, x); return *this; }
83
84 double y() const { return get(2); }
85 double y2() const { return sqr(y()); }
86 FourVector& setY(const double y) { set(2, y); return *this; }
87
88 double z() const { return get(3); }
89 double z2() const { return sqr(z()); }
90 FourVector& setZ(const double z) { set(3, z); return *this; }
91
92 double invariant() const {
93 // Done this way for numerical precision
94 return (t() + z())*(t() - z()) - x()*x() - y()*y();
95 }
96
97 bool isNull() const {
98 return Rivet::isZero(invariant());
99 }
100
102 double angle(const FourVector& v) const {
103 return vector3().angle( v.vector3() );
104 }
105
106 double angle(const Vector3& v3) const {
107 return vector3().angle(v3);
108 }
109
113 double polarRadius2() const {
114 return vector3().polarRadius2();
115 }
116
117 double perp2() const {
118 return vector3().perp2();
119 }
120
121 double rho2() const {
122 return vector3().rho2();
123 }
124
126 double polarRadius() const {
127 return vector3().polarRadius();
128 }
129
130 double perp() const {
131 return vector3().perp();
132 }
133
134 double rho() const {
135 return vector3().rho();
136 }
137
140 return vector3().polarVec();
141 }
142
143 Vector3 perpVec() const {
144 return vector3().perpVec();
145 }
146
147 Vector3 rhoVec() const {
148 return vector3().rhoVec();
149 }
150
152 double azimuthalAngle(const PhiMapping mapping=ZERO_2PI) const {
153 return vector3().azimuthalAngle(mapping);
154 }
155
156 double phi(const PhiMapping mapping=ZERO_2PI) const {
157 return vector3().phi(mapping);
158 }
159
161 double polarAngle() const {
162 return vector3().polarAngle();
163 }
164
165 double theta() const {
166 return vector3().theta();
167 }
168
170 double pseudorapidity() const {
171 return vector3().pseudorapidity();
172 }
173
174 double eta() const {
175 return vector3().eta();
176 }
177
179 double abspseudorapidity() const { return fabs(eta()); }
181 double abseta() const { return fabs(eta()); }
182
184 Vector3 vector3() const {
185 return Vector3(get(1), get(2), get(3));
186 }
187
188 Vector3 v3() const {
189 return vector3();
190 }
191
193 operator Vector3 () const { return vector3(); }
194
195
196 public:
197
199 double contract(const FourVector& v) const {
200 const double result = t()*v.t() - x()*v.x() - y()*v.y() - z()*v.z();
201 return result;
202 }
203
205 double dot(const FourVector& v) const {
206 return contract(v);
207 }
208
210 double operator * (const FourVector& v) const {
211 return contract(v);
212 }
213
215 FourVector& operator *= (double a) {
216 _vec = multiply(a, *this)._vec;
217 return *this;
218 }
219
221 FourVector& operator /= (double a) {
222 _vec = multiply(1.0/a, *this)._vec;
223 return *this;
224 }
225
227 FourVector& operator += (const FourVector& v) {
228 _vec = add(*this, v)._vec;
229 return *this;
230 }
231
233 FourVector& operator -= (const FourVector& v) {
234 _vec = add(*this, -v)._vec;
235 return *this;
236 }
237
239 FourVector operator - () const {
240 FourVector result;
241 result._vec = -_vec;
242 return result;
243 }
244
246 FourVector reverse() const {
247 FourVector result = -*this;
248 result.setT(-result.t());
249 return result;
250 }
251
252 };
253
254
256 inline double contract(const FourVector& a, const FourVector& b) {
257 return a.contract(b);
258 }
259
261 inline double dot(const FourVector& a, const FourVector& b) {
262 return contract(a, b);
263 }
264
265 inline FourVector multiply(const double a, const FourVector& v) {
266 FourVector result;
267 result._vec = a * v._vec;
268 return result;
269 }
270
271 inline FourVector multiply(const FourVector& v, const double a) {
272 return multiply(a, v);
273 }
274
275 inline FourVector operator * (const double a, const FourVector& v) {
276 return multiply(a, v);
277 }
278
279 inline FourVector operator * (const FourVector& v, const double a) {
280 return multiply(a, v);
281 }
282
283 inline FourVector operator / (const FourVector& v, const double a) {
284 return multiply(1.0/a, v);
285 }
286
287 inline FourVector add(const FourVector& a, const FourVector& b) {
288 FourVector result;
289 result._vec = a._vec + b._vec;
290 return result;
291 }
292
293 inline FourVector operator+(const FourVector& a, const FourVector& b) {
294 return add(a, b);
295 }
296
297 inline FourVector operator-(const FourVector& a, const FourVector& b) {
298 return add(a, -b);
299 }
300
303 inline double invariant(const FourVector& lv) {
304 return lv.invariant();
305 }
306
308 inline double angle(const FourVector& a, const FourVector& b) {
309 return a.angle(b);
310 }
311
313 inline double angle(const Vector3& a, const FourVector& b) {
314 return angle( a, b.vector3() );
315 }
316
318 inline double angle(const FourVector& a, const Vector3& b) {
319 return a.angle(b);
320 }
321
322
324
325
327 class FourMomentum : public FourVector {
328 friend FourMomentum multiply(const double a, const FourMomentum& v);
329 friend FourMomentum multiply(const FourMomentum& v, const double a);
330 friend FourMomentum add(const FourMomentum& a, const FourMomentum& b);
331 friend FourMomentum transform(const LorentzTransform& lt, const FourMomentum& v4);
332
333 public:
334 FourMomentum() { }
335
336 template<typename V4TYPE, typename std::enable_if<HasXYZT<V4TYPE>::value, int>::type DUMMY=0>
337 FourMomentum(const V4TYPE& other) {
338 this->setE(other.t());
339 this->setPx(other.x());
340 this->setPy(other.y());
341 this->setPz(other.z());
342 }
343
344 FourMomentum(const Vector<4>& other)
345 : FourVector(other) { }
346
347 FourMomentum(const double E, const double px, const double py, const double pz) {
348 this->setE(E);
349 this->setPx(px);
350 this->setPy(py);
351 this->setPz(pz);
352 }
353
354 FourMomentum(const Vector3& v3, double m) {
355 const double e = sqrt(v3.mod2() + m*m);
356 this->setE(e);
357 this->setPx(v3.x());
358 this->setPy(v3.y());
359 this->setPz(v3.z());
360 }
361
362 ~FourMomentum() {}
363
364 public:
365
366
369
371 FourMomentum& setE(double E) {
372 setT(E);
373 return *this;
374 }
375
377 FourMomentum& setPx(double px) {
378 setX(px);
379 return *this;
380 }
381
383 FourMomentum& setPy(double py) {
384 setY(py);
385 return *this;
386 }
387
389 FourMomentum& setPz(double pz) {
390 setZ(pz);
391 return *this;
392 }
393
394
396 FourMomentum& setPE(double px, double py, double pz, double E) {
397 if (E < 0)
398 throw std::invalid_argument("Negative energy given as argument: " + to_str(E));
399 setPx(px); setPy(py); setPz(pz); setE(E);
400 return *this;
401 }
402
403 FourMomentum& setXYZE(double px, double py, double pz, double E) {
404 return setPE(px, py, pz, E);
405 }
406 // /// Near-alias with switched arg order
407 // FourMomentum& setEP(double E, double px, double py, double pz) {
408 // return setPE(px, py, pz, E);
409 // }
410 // /// Alias for setEP
411 // FourMomentum& setEXYZ(double E, double px, double py, double pz) {
412 // return setEP(E, px, py, pz);
413 // }
414
415
417 FourMomentum& setPM(double px, double py, double pz, double mass) {
418 if (mass < 0)
419 throw std::invalid_argument("Negative mass given as argument: " + to_str(mass));
420 const double E = sqrt( sqr(mass) + sqr(px) + sqr(py) + sqr(pz) );
421 // setPx(px); setPy(py); setPz(pz); setE(E);
422 return setPE(px, py, pz, E);
423 }
424
425 FourMomentum& setXYZM(double px, double py, double pz, double mass) {
426 return setPM(px, py, pz, mass);
427 }
428
429
434 FourMomentum& setEtaPhiME(double eta, double phi, double mass, double E) {
435 if (mass < 0)
436 throw std::invalid_argument("Negative mass given as argument");
437 if (E < 0)
438 throw std::invalid_argument("Negative energy given as argument");
439 const double theta = 2 * atan(exp(-eta));
440 if (theta < 0 || theta > M_PI)
441 throw std::domain_error("Polar angle outside 0..pi in calculation");
443 return *this;
444 }
445
450 FourMomentum& setEtaPhiMPt(double eta, double phi, double mass, double pt) {
451 if (mass < 0)
452 throw std::invalid_argument("Negative mass given as argument");
453 if (pt < 0)
454 throw std::invalid_argument("Negative transverse momentum given as argument");
455 const double theta = 2 * atan(exp(-eta));
456 if (theta < 0 || theta > M_PI)
457 throw std::domain_error("Polar angle outside 0..pi in calculation");
458 const double p = pt / sin(theta);
459 const double E = sqrt( sqr(p) + sqr(mass) );
461 return *this;
462 }
463
472 FourMomentum& setRapPhiME(double y, double phi, double mass, double E) {
473 if (mass < 0)
474 throw std::invalid_argument("Negative mass given as argument");
475 if (E < 0)
476 throw std::invalid_argument("Negative energy given as argument");
477 const double sqrt_pt2_m2 = E / cosh(y);
478 const double pt = sqrt( sqr(sqrt_pt2_m2) - sqr(mass) );
479 if (pt < 0)
480 throw std::domain_error("Negative transverse momentum in calculation");
481 const double pz = sqrt_pt2_m2 * sinh(y);
482 const double px = pt * cos(phi);
483 const double py = pt * sin(phi);
484 setPE(px, py, pz, E);
485 return *this;
486 }
487
492 FourMomentum& setRapPhiMPt(double y, double phi, double mass, double pt) {
493 if (mass < 0)
494 throw std::invalid_argument("Negative mass given as argument");
495 if (pt < 0)
496 throw std::invalid_argument("Negative transverse mass given as argument");
497 const double E = sqrt( sqr(pt) + sqr(mass) ) * cosh(y);
498 if (E < 0)
499 throw std::domain_error("Negative energy in calculation");
500 setRapPhiME(y, phi, mass, E);
501 return *this;
502 }
503
509 FourMomentum& setThetaPhiME(double theta, double phi, double mass, double E) {
510 if (theta < 0 || theta > M_PI)
511 throw std::invalid_argument("Polar angle outside 0..pi given as argument");
512 if (mass < 0)
513 throw std::invalid_argument("Negative mass given as argument");
514 if (E < 0)
515 throw std::invalid_argument("Negative energy given as argument");
516 const double p = sqrt( sqr(E) - sqr(mass) );
517 const double pz = p * cos(theta);
518 const double pt = p * sin(theta);
519 if (pt < 0)
520 throw std::invalid_argument("Negative transverse momentum in calculation");
521 const double px = pt * cos(phi);
522 const double py = pt * sin(phi);
523 setPE(px, py, pz, E);
524 return *this;
525 }
526
532 FourMomentum& setThetaPhiMPt(double theta, double phi, double mass, double pt) {
533 if (theta < 0 || theta > M_PI)
534 throw std::invalid_argument("Polar angle outside 0..pi given as argument");
535 if (mass < 0)
536 throw std::invalid_argument("Negative mass given as argument");
537 if (pt < 0)
538 throw std::invalid_argument("Negative transverse momentum given as argument");
539 const double p = pt / sin(theta);
540 const double px = pt * cos(phi);
541 const double py = pt * sin(phi);
542 const double pz = p * cos(theta);
543 const double E = sqrt( sqr(p) + sqr(mass) );
544 setPE(px, py, pz, E);
545 return *this;
546 }
547
551 FourMomentum& setPtPhiME(double pt, double phi, double mass, double E) {
552 if (pt < 0)
553 throw std::invalid_argument("Negative transverse momentum given as argument");
554 if (mass < 0)
555 throw std::invalid_argument("Negative mass given as argument");
556 if (E < 0)
557 throw std::invalid_argument("Negative energy given as argument");
558 const double px = pt * cos(phi);
559 const double py = pt * sin(phi);
560 const double pz = sqrt(sqr(E) - sqr(mass) - sqr(pt));
561 setPE(px, py, pz, E);
562 return *this;
563 }
564
566
567
570
572 double E() const { return t(); }
574 double E2() const { return t2(); }
575
577 double px() const { return x(); }
579 double px2() const { return x2(); }
580
582 double py() const { return y(); }
584 double py2() const { return y2(); }
585
587 double pz() const { return z(); }
589 double pz2() const { return z2(); }
590
591
595 double mass() const {
596 // assert(Rivet::isZero(mass2()) || mass2() > 0);
597 // if (Rivet::isZero(mass2())) {
598 // return 0.0;
599 // } else {
600 // return sqrt(mass2());
601 // }
602 return sign(mass2()) * sqrt(fabs(mass2()));
603 }
604
606 double mass2() const {
607 return invariant();
608 }
609
610
612 Vector3 p3() const { return vector3(); }
613
615 double p() const {
616 return p3().mod();
617 }
618
620 double p2() const {
621 return p3().mod2();
622 }
623
624
626 double rapidity() const {
627 if (E() == 0.0) return 0.0;
628 if (E() == fabs(pz())) return std::copysign(INF, pz());
629 return 0.5 * std::log( (E() + pz()) / (E() - pz()) );
630 }
631
632 double rap() const {
633 return rapidity();
634 }
635
637 double absrapidity() const {
638 return fabs(rapidity());
639 }
640
641 double absrap() const {
642 return fabs(rap());
643 }
644
646 Vector3 pTvec() const {
647 return p3().polarVec();
648 }
649
650 Vector3 ptvec() const {
651 return pTvec();
652 }
653
655 double pT2() const {
656 return vector3().polarRadius2();
657 }
658
659 double pt2() const {
660 return vector3().polarRadius2();
661 }
662
664 double pT() const {
665 return sqrt(pT2());
666 }
667
668 double pt() const {
669 return sqrt(pT2());
670 }
671
673 double Et2() const {
674 return Et() * Et();
675 }
676
677 double Et() const {
678 return E() * sin(polarAngle());
679 }
680
682
683
686
689 double gamma() const {
690 return sqrt(E2()/mass2());
691 }
692
696 return gamma() * p3().unit();
697 }
698
701 double beta() const {
702 return p()/E();
703 }
704
707 Vector3 betaVec() const {
708 // return Vector3(px()/E(), py()/E(), pz()/E());
709 return p3()/E();
710 }
711
713
714
716
719
721 FourMomentum& operator*=(double a) {
722 _vec = multiply(a, *this)._vec;
723 return *this;
724 }
725
727 FourMomentum& operator/=(double a) {
728 _vec = multiply(1.0/a, *this)._vec;
729 return *this;
730 }
731
733 FourMomentum& operator+=(const FourMomentum& v) {
734 _vec = add(*this, v)._vec;
735 return *this;
736 }
737
739 FourMomentum& operator-=(const FourMomentum& v) {
740 _vec = add(*this, -v)._vec;
741 return *this;
742 }
743
745 FourMomentum operator-() const {
746 FourMomentum result;
747 result._vec = -_vec;
748 return result;
749 }
750
752 FourMomentum reverse() const {
753 FourMomentum result = -*this;
754 result.setE(-result.E());
755 return result;
756 }
757
759
760
762
763
766
768 static FourMomentum mkXYZE(double px, double py, double pz, double E) {
769 return FourMomentum().setPE(px, py, pz, E);
770 }
771
773 static FourMomentum mkXYZM(double px, double py, double pz, double mass) {
774 return FourMomentum().setPM(px, py, pz, mass);
775 }
776
778 static FourMomentum mkEtaPhiME(double eta, double phi, double mass, double E) {
779 return FourMomentum().setEtaPhiME(eta, phi, mass, E);
780 }
781
783 static FourMomentum mkEtaPhiMPt(double eta, double phi, double mass, double pt) {
784 return FourMomentum().setEtaPhiMPt(eta, phi, mass, pt);
785 }
786
788 static FourMomentum mkRapPhiME(double y, double phi, double mass, double E) {
789 return FourMomentum().setRapPhiME(y, phi, mass, E);
790 }
791
793 static FourMomentum mkRapPhiMPt(double y, double phi, double mass, double pt) {
794 return FourMomentum().setRapPhiMPt(y, phi, mass, pt);
795 }
796
798 static FourMomentum mkThetaPhiME(double theta, double phi, double mass, double E) {
799 return FourMomentum().setThetaPhiME(theta, phi, mass, E);
800 }
801
803 static FourMomentum mkThetaPhiMPt(double theta, double phi, double mass, double pt) {
804 return FourMomentum().setThetaPhiMPt(theta, phi, mass, pt);
805 }
806
808 static FourMomentum mkPtPhiME(double pt, double phi, double mass, double E) {
809 return FourMomentum().setPtPhiME(pt, phi, mass, E);
810 }
811
813
814
815 };
816
817
818 inline FourMomentum multiply(const double a, const FourMomentum& v) {
819 FourMomentum result;
820 result._vec = a * v._vec;
821 return result;
822 }
823
824 inline FourMomentum multiply(const FourMomentum& v, const double a) {
825 return multiply(a, v);
826 }
827
828 inline FourMomentum operator*(const double a, const FourMomentum& v) {
829 return multiply(a, v);
830 }
831
832 inline FourMomentum operator*(const FourMomentum& v, const double a) {
833 return multiply(a, v);
834 }
835
836 inline FourMomentum operator/(const FourMomentum& v, const double a) {
837 return multiply(1.0/a, v);
838 }
839
840 inline FourMomentum add(const FourMomentum& a, const FourMomentum& b) {
841 FourMomentum result;
842 result._vec = a._vec + b._vec;
843 return result;
844 }
845
846 inline FourMomentum operator+(const FourMomentum& a, const FourMomentum& b) {
847 return add(a, b);
848 }
849
850 inline FourMomentum operator-(const FourMomentum& a, const FourMomentum& b) {
851 return add(a, -b);
852 }
853
854
856
857
860
869 inline double deltaR2(const FourVector& a, const FourVector& b,
870 RapScheme scheme=PSEUDORAPIDITY) {
871 switch (scheme) {
872 case PSEUDORAPIDITY :
873 return deltaR2(a.vector3(), b.vector3());
874 case RAPIDITY:
875 {
876 const FourMomentum* ma = dynamic_cast<const FourMomentum*>(&a);
877 const FourMomentum* mb = dynamic_cast<const FourMomentum*>(&b);
878 if (!ma || !mb) {
879 string err = "deltaR with scheme RAPIDITY can only be called with FourMomentum objects, not FourVectors";
880 throw std::runtime_error(err);
881 }
882 return deltaR2(*ma, *mb, scheme);
883 }
884 default:
885 throw std::runtime_error("The specified deltaR scheme is not yet implemented");
886 }
887 }
888
897 inline double deltaR(const FourVector& a, const FourVector& b,
898 RapScheme scheme=PSEUDORAPIDITY) {
899 return sqrt(deltaR2(a, b, scheme));
900 }
901
902
903
910 inline double deltaR2(const FourVector& v,
911 double eta2, double phi2,
912 RapScheme scheme=PSEUDORAPIDITY) {
913 switch (scheme) {
914 case PSEUDORAPIDITY :
915 return deltaR2(v.vector3(), eta2, phi2);
916 case RAPIDITY:
917 {
918 const FourMomentum* mv = dynamic_cast<const FourMomentum*>(&v);
919 if (!mv) {
920 string err = "deltaR with scheme RAPIDITY can only be called with FourMomentum objects, not FourVectors";
921 throw std::runtime_error(err);
922 }
923 return deltaR2(*mv, eta2, phi2, scheme);
924 }
925 default:
926 throw std::runtime_error("The specified deltaR scheme is not yet implemented");
927 }
928 }
929
936 inline double deltaR(const FourVector& v,
937 double eta2, double phi2,
938 RapScheme scheme=PSEUDORAPIDITY) {
939 return sqrt(deltaR2(v, eta2, phi2, scheme));
940 }
941
942
949 inline double deltaR2(double eta1, double phi1,
950 const FourVector& v,
951 RapScheme scheme=PSEUDORAPIDITY) {
952 switch (scheme) {
953 case PSEUDORAPIDITY :
954 return deltaR2(eta1, phi1, v.vector3());
955 case RAPIDITY:
956 {
957 const FourMomentum* mv = dynamic_cast<const FourMomentum*>(&v);
958 if (!mv) {
959 string err = "deltaR with scheme RAPIDITY can only be called with FourMomentum objects, not FourVectors";
960 throw std::runtime_error(err);
961 }
962 return deltaR2(eta1, phi1, *mv, scheme);
963 }
964 default:
965 throw std::runtime_error("The specified deltaR scheme is not yet implemented");
966 }
967 }
968
975 inline double deltaR(double eta1, double phi1,
976 const FourVector& v,
977 RapScheme scheme=PSEUDORAPIDITY) {
978 return sqrt(deltaR2(eta1, phi1, v, scheme));
979 }
980
981
988 inline double deltaR2(const FourMomentum& a, const FourMomentum& b,
989 RapScheme scheme=PSEUDORAPIDITY) {
990 switch (scheme) {
991 case PSEUDORAPIDITY:
992 return deltaR2(a.vector3(), b.vector3());
993 case RAPIDITY:
994 return deltaR2(a.rapidity(), a.azimuthalAngle(), b.rapidity(), b.azimuthalAngle());
995 default:
996 throw std::runtime_error("The specified deltaR scheme is not yet implemented");
997 }
998 }
999
1006 inline double deltaR(const FourMomentum& a, const FourMomentum& b,
1007 RapScheme scheme=PSEUDORAPIDITY) {
1008 return sqrt(deltaR2(a, b, scheme));
1009 }
1010
1011
1017 inline double deltaR2(const FourMomentum& v,
1018 double eta2, double phi2,
1019 RapScheme scheme=PSEUDORAPIDITY) {
1020 switch (scheme) {
1021 case PSEUDORAPIDITY:
1022 return deltaR2(v.vector3(), eta2, phi2);
1023 case RAPIDITY:
1024 return deltaR2(v.rapidity(), v.azimuthalAngle(), eta2, phi2);
1025 default:
1026 throw std::runtime_error("The specified deltaR scheme is not yet implemented");
1027 }
1028 }
1029
1035 inline double deltaR(const FourMomentum& v,
1036 double eta2, double phi2,
1037 RapScheme scheme=PSEUDORAPIDITY) {
1038 return sqrt(deltaR2(v, eta2, phi2, scheme));
1039 }
1040
1041
1047 inline double deltaR2(double eta1, double phi1,
1048 const FourMomentum& v,
1049 RapScheme scheme=PSEUDORAPIDITY) {
1050 switch (scheme) {
1051 case PSEUDORAPIDITY:
1052 return deltaR2(eta1, phi1, v.vector3());
1053 case RAPIDITY:
1054 return deltaR2(eta1, phi1, v.rapidity(), v.azimuthalAngle());
1055 default:
1056 throw std::runtime_error("The specified deltaR scheme is not yet implemented");
1057 }
1058 }
1059
1065 inline double deltaR(double eta1, double phi1,
1066 const FourMomentum& v,
1067 RapScheme scheme=PSEUDORAPIDITY) {
1068 return sqrt(deltaR2(eta1, phi1, v, scheme));
1069 }
1070
1071
1077 inline double deltaR2(const FourMomentum& a, const FourVector& b,
1078 RapScheme scheme=PSEUDORAPIDITY) {
1079 switch (scheme) {
1080 case PSEUDORAPIDITY:
1081 return deltaR2(a.vector3(), b.vector3());
1082 case RAPIDITY:
1084 default:
1085 throw std::runtime_error("The specified deltaR scheme is not yet implemented");
1086 }
1087 }
1088
1094 inline double deltaR(const FourMomentum& a, const FourVector& b,
1095 RapScheme scheme=PSEUDORAPIDITY) {
1096 return sqrt(deltaR2(a, b, scheme));
1097 }
1098
1099
1105 inline double deltaR2(const FourVector& a, const FourMomentum& b,
1106 RapScheme scheme=PSEUDORAPIDITY) {
1107 return deltaR2(b, a, scheme); //< note reversed args
1108 }
1109
1115 inline double deltaR(const FourVector& a, const FourMomentum& b,
1116 RapScheme scheme=PSEUDORAPIDITY) {
1117 return deltaR(b, a, scheme); //< note reversed args
1118 }
1119
1120
1123 inline double deltaR2(const FourMomentum& a, const Vector3& b) {
1124 return deltaR2(a.vector3(), b);
1125 }
1126
1129 inline double deltaR(const FourMomentum& a, const Vector3& b) {
1130 return deltaR(a.vector3(), b);
1131 }
1132
1135 inline double deltaR2(const Vector3& a, const FourMomentum& b) {
1136 return deltaR2(a, b.vector3());
1137 }
1138
1141 inline double deltaR(const Vector3& a, const FourMomentum& b) {
1142 return deltaR(a, b.vector3());
1143 }
1144
1147 inline double deltaR2(const FourVector& a, const Vector3& b) {
1148 return deltaR2(a.vector3(), b);
1149 }
1150
1153 inline double deltaR(const FourVector& a, const Vector3& b) {
1154 return deltaR(a.vector3(), b);
1155 }
1156
1159 inline double deltaR2(const Vector3& a, const FourVector& b) {
1160 return deltaR2(a, b.vector3());
1161 }
1162
1165 inline double deltaR(const Vector3& a, const FourVector& b) {
1166 return deltaR(a, b.vector3());
1167 }
1168
1170
1171
1173
1174
1177
1179 inline double deltaPhi(const FourMomentum& a, const FourMomentum& b, bool sign=false) {
1180 return deltaPhi(a.vector3(), b.vector3(), sign);
1181 }
1182
1184 inline double deltaPhi(const FourMomentum& v, double phi2, bool sign=false) {
1185 return deltaPhi(v.vector3(), phi2, sign);
1186 }
1187
1189 inline double deltaPhi(double phi1, const FourMomentum& v, bool sign=false) {
1190 return deltaPhi(phi1, v.vector3(), sign);
1191 }
1192
1194 inline double deltaPhi(const FourVector& a, const FourVector& b, bool sign=false) {
1195 return deltaPhi(a.vector3(), b.vector3(), sign);
1196 }
1197
1199 inline double deltaPhi(const FourVector& v, double phi2, bool sign=false) {
1200 return deltaPhi(v.vector3(), phi2, sign);
1201 }
1202
1204 inline double deltaPhi(double phi1, const FourVector& v, bool sign=false) {
1205 return deltaPhi(phi1, v.vector3(), sign);
1206 }
1207
1209 inline double deltaPhi(const FourVector& a, const FourMomentum& b, bool sign=false) {
1210 return deltaPhi(a.vector3(), b.vector3(), sign);
1211 }
1212
1214 inline double deltaPhi(const FourMomentum& a, const FourVector& b, bool sign=false) {
1215 return deltaPhi(a.vector3(), b.vector3(), sign);
1216 }
1217
1219 inline double deltaPhi(const FourVector& a, const Vector3& b, bool sign=false) {
1220 return deltaPhi(a.vector3(), b, sign);
1221 }
1222
1224 inline double deltaPhi(const Vector3& a, const FourVector& b, bool sign=false) {
1225 return deltaPhi(a, b.vector3(), sign);
1226 }
1227
1229 inline double deltaPhi(const FourMomentum& a, const Vector3& b, bool sign=false) {
1230 return deltaPhi(a.vector3(), b, sign);
1231 }
1232
1234 inline double deltaPhi(const Vector3& a, const FourMomentum& b, bool sign=false) {
1235 return deltaPhi(a, b.vector3(), sign);
1236 }
1237
1239
1240
1242
1243
1246
1248 inline double deltaEta(const FourMomentum& a, const FourMomentum& b, bool sign=false) {
1249 return deltaEta(a.vector3(), b.vector3(), sign);
1250 }
1251
1253 inline double deltaEta(const FourMomentum& v, double eta2, bool sign=false) {
1254 return deltaEta(v.vector3(), eta2, sign);
1255 }
1256
1258 inline double deltaEta(double eta1, const FourMomentum& v, bool sign=false) {
1259 return deltaEta(eta1, v.vector3(), sign);
1260 }
1261
1263 inline double deltaEta(const FourVector& a, const FourVector& b, bool sign=false) {
1264 return deltaEta(a.vector3(), b.vector3(), sign);
1265 }
1266
1268 inline double deltaEta(const FourVector& v, double eta2, bool sign=false) {
1269 return deltaEta(v.vector3(), eta2, sign);
1270 }
1271
1273 inline double deltaEta(double eta1, const FourVector& v, bool sign=false) {
1274 return deltaEta(eta1, v.vector3(), sign);
1275 }
1276
1278 inline double deltaEta(const FourVector& a, const FourMomentum& b, bool sign=false) {
1279 return deltaEta(a.vector3(), b.vector3(), sign);
1280 }
1281
1283 inline double deltaEta(const FourMomentum& a, const FourVector& b, bool sign=false) {
1284 return deltaEta(a.vector3(), b.vector3(), sign);
1285 }
1286
1288 inline double deltaEta(const FourVector& a, const Vector3& b, bool sign=false) {
1289 return deltaEta(a.vector3(), b, sign);
1290 }
1291
1293 inline double deltaEta(const Vector3& a, const FourVector& b, bool sign=false) {
1294 return deltaEta(a, b.vector3(), sign);
1295 }
1296
1298 inline double deltaEta(const FourMomentum& a, const Vector3& b, bool sign=false) {
1299 return deltaEta(a.vector3(), b, sign);
1300 }
1301
1303 inline double deltaEta(const Vector3& a, const FourMomentum& b, bool sign=false) {
1304 return deltaEta(a, b.vector3(), sign);
1305 }
1306
1308
1309
1312
1314 inline double deltaRap(const FourMomentum& a, const FourMomentum& b, bool sign=false) {
1315 return deltaRap(a.rapidity(), b.rapidity(), sign);
1316 }
1317
1319 inline double deltaRap(const FourMomentum& v, double y2, bool sign=false) {
1320 return deltaRap(v.rapidity(), y2, sign);
1321 }
1322
1324 inline double deltaRap(double y1, const FourMomentum& v, bool sign=false) {
1325 return deltaRap(y1, v.rapidity(), sign);
1326 }
1327
1329
1330
1332
1333
1336
1339
1341 inline bool cmpMomByPt(const FourMomentum& a, const FourMomentum& b) {
1342 return a.pt() > b.pt();
1343 }
1344
1345 inline bool cmpMomByAscPt(const FourMomentum& a, const FourMomentum& b) {
1346 return a.pt() < b.pt();
1347 }
1348
1350 inline bool cmpMomByP(const FourMomentum& a, const FourMomentum& b) {
1351 return a.vector3().mod() > b.vector3().mod();
1352 }
1353
1354 inline bool cmpMomByAscP(const FourMomentum& a, const FourMomentum& b) {
1355 return a.vector3().mod() < b.vector3().mod();
1356 }
1357
1359 inline bool cmpMomByEt(const FourMomentum& a, const FourMomentum& b) {
1360 return a.Et() > b.Et();
1361 }
1362
1363 inline bool cmpMomByAscEt(const FourMomentum& a, const FourMomentum& b) {
1364 return a.Et() < b.Et();
1365 }
1366
1368 inline bool cmpMomByE(const FourMomentum& a, const FourMomentum& b) {
1369 return a.E() > b.E();
1370 }
1371
1372 inline bool cmpMomByAscE(const FourMomentum& a, const FourMomentum& b) {
1373 return a.E() < b.E();
1374 }
1375
1377 inline bool cmpMomByMass(const FourMomentum& a, const FourMomentum& b) {
1378 return a.mass() > b.mass();
1379 }
1380
1381 inline bool cmpMomByAscMass(const FourMomentum& a, const FourMomentum& b) {
1382 return a.mass() < b.mass();
1383 }
1384
1386 inline bool cmpMomByEta(const FourMomentum& a, const FourMomentum& b) {
1387 return a.eta() < b.eta();
1388 }
1389
1391 inline bool cmpMomByDescEta(const FourMomentum& a, const FourMomentum& b) {
1392 return a.pseudorapidity() > b.pseudorapidity();
1393 }
1394
1396 inline bool cmpMomByAbsEta(const FourMomentum& a, const FourMomentum& b) {
1397 return fabs(a.eta()) < fabs(b.eta());
1398 }
1399
1401 inline bool cmpMomByDescAbsEta(const FourMomentum& a, const FourMomentum& b) {
1402 return fabs(a.eta()) > fabs(b.eta());
1403 }
1404
1406 inline bool cmpMomByRap(const FourMomentum& a, const FourMomentum& b) {
1407 return a.rapidity() < b.rapidity();
1408 }
1409
1411 inline bool cmpMomByDescRap(const FourMomentum& a, const FourMomentum& b) {
1412 return a.rapidity() > b.rapidity();
1413 }
1414
1416 inline bool cmpMomByAbsRap(const FourMomentum& a, const FourMomentum& b) {
1417 return fabs(a.rapidity()) < fabs(b.rapidity());
1418 }
1419
1421 inline bool cmpMomByDescAbsRap(const FourMomentum& a, const FourMomentum& b) {
1422 return fabs(a.rapidity()) > fabs(b.rapidity());
1423 }
1424
1426
1427
1429 template<typename MOMS, typename CMP>
1430 inline MOMS& isortBy(MOMS& pbs, const CMP& cmp) {
1431 std::sort(pbs.begin(), pbs.end(), cmp);
1432 return pbs;
1433 }
1434
1435 template<typename MOMS, typename CMP>
1436 inline MOMS sortBy(const MOMS& pbs, const CMP& cmp) {
1437 MOMS rtn = pbs;
1438 std::sort(rtn.begin(), rtn.end(), cmp);
1439 return rtn;
1440 }
1441
1443 template<typename MOMS>
1444 inline MOMS& isortByPt(MOMS& pbs) {
1445 return isortBy(pbs, cmpMomByPt);
1446 }
1447
1448 template<typename MOMS>
1449 inline MOMS sortByPt(const MOMS& pbs) {
1450 return sortBy(pbs, cmpMomByPt);
1451 }
1452
1454 template<typename MOMS>
1455 inline MOMS& isortByE(MOMS& pbs) {
1456 return isortBy(pbs, cmpMomByE);
1457 }
1458
1459 template<typename MOMS>
1460 inline MOMS sortByE(const MOMS& pbs) {
1461 return sortBy(pbs, cmpMomByE);
1462 }
1463
1465 template<typename MOMS>
1466 inline MOMS& isortByEt(MOMS& pbs) {
1467 return isortBy(pbs, cmpMomByEt);
1468 }
1469
1470 template<typename MOMS>
1471 inline MOMS sortByEt(const MOMS& pbs) {
1472 return sortBy(pbs, cmpMomByEt);
1473 }
1474
1476
1477
1480
1482 inline double mass(const FourMomentum& a, const FourMomentum& b) {
1483 return (a + b).mass();
1484 }
1485
1487 inline double mass2(const FourMomentum& a, const FourMomentum& b) {
1488 return (a + b).mass2();
1489 }
1490
1497 inline double mT(const FourMomentum& vis, const FourMomentum& invis) {
1498 return mT(vis.p3(), invis.p3());
1499 }
1500
1507 inline double mT(const FourMomentum& vis, const Vector3& invis) {
1508 return mT(vis.p3(), invis);
1509 }
1510
1517 inline double mT(const Vector3& vis, const FourMomentum& invis) {
1518 return mT(vis, invis.p3());
1519 }
1520
1527 inline double mT(const FourMomentum& vis, const Vector4& invis) {
1528 return mT(vis.p3(), invis.vector3());
1529 }
1530
1537 inline double mT(const Vector4& vis, const FourMomentum& invis) {
1538 return mT(vis.vector3(), invis.p3());
1539 }
1540
1542 inline double pT(const FourMomentum& vis, const FourMomentum& invis) {
1543 return pT(vis.p3(), invis.p3());
1544 }
1545
1547 inline double pT(const FourMomentum& vis, const Vector3& invis) {
1548 return pT(vis.p3(), invis);
1549 }
1550
1552 inline double pT(const Vector3& vis, const FourMomentum& invis) {
1553 return pT(vis, invis.p3());
1554 }
1555
1557
1558
1560
1561
1564
1566 inline std::string toString(const FourVector& lv) {
1567 std::ostringstream out;
1568 out << "(" << (fabs(lv.t()) < 1E-30 ? 0.0 : lv.t())
1569 << "; " << (fabs(lv.x()) < 1E-30 ? 0.0 : lv.x())
1570 << ", " << (fabs(lv.y()) < 1E-30 ? 0.0 : lv.y())
1571 << ", " << (fabs(lv.z()) < 1E-30 ? 0.0 : lv.z())
1572 << ")";
1573 return out.str();
1574 }
1575
1577 inline std::ostream& operator<<(std::ostream& out, const FourVector& lv) {
1578 out << toString(lv);
1579 return out;
1580 }
1581
1583
1586 typedef std::vector<FourVector> FourVectors;
1587 typedef std::vector<FourMomentum> FourMomenta;
1589
1591
1592
1593}
1594
1595#endif
Specialized version of the FourVector with momentum/energy functionality.
Definition Vector4.hh:327
double p2() const
Get the modulus-squared of the 3-momentum.
Definition Vector4.hh:620
double E2() const
Get energy-squared .
Definition Vector4.hh:574
double mass() const
Get the mass (the Lorentz self-invariant).
Definition Vector4.hh:595
double pz2() const
Get z-squared .
Definition Vector4.hh:589
double rapidity() const
Calculate the rapidity.
Definition Vector4.hh:626
FourMomentum reverse() const
Multiply space components only by -1.
Definition Vector4.hh:752
double px() const
Get x-component of momentum .
Definition Vector4.hh:577
FourMomentum & setPtPhiME(double pt, double phi, double mass, double E)
Definition Vector4.hh:551
FourMomentum & setRapPhiMPt(double y, double phi, double mass, double pt)
Definition Vector4.hh:492
Vector3 p3() const
Get 3-momentum part, .
Definition Vector4.hh:612
FourMomentum & setPz(double pz)
Set z-component of momentum .
Definition Vector4.hh:389
double p() const
Get the modulus of the 3-momentum.
Definition Vector4.hh:615
FourMomentum & setXYZE(double px, double py, double pz, double E)
Alias for setPE.
Definition Vector4.hh:403
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:798
Vector3 gammaVec() const
Definition Vector4.hh:695
double pt() const
Calculate the transverse momentum .
Definition Vector4.hh:668
FourMomentum & setPE(double px, double py, double pz, double E)
Set the p coordinates and energy simultaneously.
Definition Vector4.hh:396
double py() const
Get y-component of momentum .
Definition Vector4.hh:582
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:788
double beta() const
Definition Vector4.hh:701
double Et() const
Calculate the transverse energy .
Definition Vector4.hh:677
double Et2() const
Calculate the transverse energy .
Definition Vector4.hh:673
static FourMomentum mkXYZE(double px, double py, double pz, double E)
Make a vector from (px,py,pz,E) coordinates.
Definition Vector4.hh:768
FourMomentum & setPy(double py)
Set y-component of momentum .
Definition Vector4.hh:383
double rap() const
Alias for rapidity.
Definition Vector4.hh:632
FourMomentum & setPx(double px)
Set x-component of momentum .
Definition Vector4.hh:377
double pt2() const
Calculate the squared transverse momentum .
Definition Vector4.hh:659
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:808
Vector3 pTvec() const
Calculate the transverse momentum vector .
Definition Vector4.hh:646
FourMomentum & setRapPhiME(double y, double phi, double mass, double E)
Definition Vector4.hh:472
Vector3 ptvec() const
Synonym for pTvec.
Definition Vector4.hh:650
double mass2() const
Get the squared mass (the Lorentz self-invariant).
Definition Vector4.hh:606
FourMomentum & setEtaPhiMPt(double eta, double phi, double mass, double pt)
Definition Vector4.hh:450
double pz() const
Get z-component of momentum .
Definition Vector4.hh:587
double py2() const
Get y-squared .
Definition Vector4.hh:584
double absrap() const
Absolute rapidity.
Definition Vector4.hh:641
double pT2() const
Calculate the squared transverse momentum .
Definition Vector4.hh:655
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:778
double px2() const
Get x-squared .
Definition Vector4.hh:579
Vector3 betaVec() const
Definition Vector4.hh:707
FourMomentum & operator-=(const FourMomentum &v)
Subtract from this 4-vector. NB time as well as space components are subtracted.
Definition Vector4.hh:739
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:773
FourMomentum & setThetaPhiME(double theta, double phi, double mass, double E)
Definition Vector4.hh:509
double E() const
Get energy (time component of momentum).
Definition Vector4.hh:572
FourMomentum & setXYZM(double px, double py, double pz, double mass)
Alias for setPM.
Definition Vector4.hh:425
FourMomentum & operator+=(const FourMomentum &v)
Add to this 4-vector. NB time as well as space components are added.
Definition Vector4.hh:733
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:783
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:793
double pT() const
Calculate the transverse momentum .
Definition Vector4.hh:664
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:803
FourMomentum & setThetaPhiMPt(double theta, double phi, double mass, double pt)
Definition Vector4.hh:532
FourMomentum operator-() const
Multiply all components (time and space) by -1.
Definition Vector4.hh:745
double absrapidity() const
Absolute rapidity.
Definition Vector4.hh:637
FourMomentum & setEtaPhiME(double eta, double phi, double mass, double E)
Definition Vector4.hh:434
FourMomentum & operator/=(double a)
Divide by a scalar.
Definition Vector4.hh:727
FourMomentum & setE(double E)
Set energy (time component of momentum).
Definition Vector4.hh:371
double gamma() const
Definition Vector4.hh:689
FourMomentum & operator*=(double a)
Multiply by a scalar.
Definition Vector4.hh:721
FourMomentum & setPM(double px, double py, double pz, double mass)
Set the p coordinates and mass simultaneously.
Definition Vector4.hh:417
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:239
double perp2() const
Synonym for polarRadius2.
Definition Vector4.hh:117
FourVector & operator/=(double a)
Divide by a scalar.
Definition Vector4.hh:221
double rho() const
Synonym for polarRadius.
Definition Vector4.hh:134
Vector3 perpVec() const
Synonym for polarVec.
Definition Vector4.hh:143
Vector3 vector3() const
Get the spatial part of the 4-vector as a 3-vector.
Definition Vector4.hh:184
double angle(const FourVector &v) const
Angle between this vector and another.
Definition Vector4.hh:102
double phi(const PhiMapping mapping=ZERO_2PI) const
Synonym for azimuthalAngle.
Definition Vector4.hh:156
double contract(const FourVector &v) const
Contract two 4-vectors, with metric signature (+ - - -).
Definition Vector4.hh:199
double abseta() const
Get the directly (alias).
Definition Vector4.hh:181
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:152
double dot(const FourVector &v) const
Contract two 4-vectors, with metric signature (+ - - -).
Definition Vector4.hh:205
double perp() const
Synonym for polarRadius.
Definition Vector4.hh:130
FourVector & operator*=(double a)
Multiply by a scalar.
Definition Vector4.hh:215
double operator*(const FourVector &v) const
Contract two 4-vectors, with metric signature (+ - - -).
Definition Vector4.hh:210
double eta() const
Synonym for pseudorapidity.
Definition Vector4.hh:174
double theta() const
Synonym for polarAngle.
Definition Vector4.hh:165
double pseudorapidity() const
Pseudorapidity (defined purely by the 3-vector components).
Definition Vector4.hh:170
Vector3 v3() const
Short "p3-like" alias for vector3.
Definition Vector4.hh:188
double polarRadius() const
Magnitude of projection of 3-vector on to the plane.
Definition Vector4.hh:126
Vector3 rhoVec() const
Synonym for polarVec.
Definition Vector4.hh:147
double angle(const Vector3 &v3) const
Angle between this vector and another (3-vector).
Definition Vector4.hh:106
Vector3 polarVec() const
Projection of 3-vector on to the plane.
Definition Vector4.hh:139
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:113
FourVector reverse() const
Multiply space components only by -1.
Definition Vector4.hh:246
double polarAngle() const
Angle subtended by the 3-vector and the z-axis.
Definition Vector4.hh:161
FourVector & operator-=(const FourVector &v)
Subtract from this 4-vector. NB time as well as space components are subtracted.
Definition Vector4.hh:233
double abspseudorapidity() const
Get the directly.
Definition Vector4.hh:179
FourVector & operator+=(const FourVector &v)
Add to this 4-vector.
Definition Vector4.hh:227
double rho2() const
Synonym for polarRadius2.
Definition Vector4.hh:121
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:139
double polarRadius() const
Polar radius.
Definition Vector3.hh:157
double rho() const
Synonym for polarRadius.
Definition Vector3.hh:165
double eta() const
Synonym for pseudorapidity.
Definition Vector3.hh:229
double perp() const
Synonym for polarRadius.
Definition Vector3.hh:161
double perp2() const
Synonym for polarRadius2.
Definition Vector3.hh:148
double pseudorapidity() const
Purely geometric approximation to rapidity.
Definition Vector3.hh:221
Vector3 unit() const
Synonym for unitVec.
Definition Vector3.hh:124
double theta() const
Synonym for polarAngle.
Definition Vector3.hh:209
double phi(const PhiMapping mapping=ZERO_2PI) const
Synonym for azimuthalAngle.
Definition Vector3.hh:182
double polarRadius2() const
Square of the polar radius (.
Definition Vector3.hh:144
Vector3 perpVec() const
Synonym for polarVec.
Definition Vector3.hh:135
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:173
double rho2() const
Synonym for polarRadius2.
Definition Vector3.hh:152
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:202
Vector3 polarVec() const
Polar projection of this vector into the x-y plane.
Definition Vector3.hh:129
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)
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:1430
MOMS & isortByE(MOMS &pbs)
Sort a container of momenta by E (decreasing) and return by reference for non-const inputs.
Definition Vector4.hh:1455
MOMS & isortByEt(MOMS &pbs)
Sort a container of momenta by Et (decreasing) and return by reference for non-const inputs.
Definition Vector4.hh:1466
bool cmpMomByDescEta(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by decreasing eta (pseudorapidity).
Definition Vector4.hh:1391
bool cmpMomByMass(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by decreasing mass.
Definition Vector4.hh:1377
bool cmpMomByP(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by decreasing 3-momentum magnitude |p|.
Definition Vector4.hh:1350
bool cmpMomByAscEt(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by increasing transverse energy.
Definition Vector4.hh:1363
bool cmpMomByRap(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by increasing rapidity.
Definition Vector4.hh:1406
bool cmpMomByE(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by decreasing energy.
Definition Vector4.hh:1368
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:1436
bool cmpMomByAbsRap(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by increasing absolute rapidity.
Definition Vector4.hh:1416
MOMS sortByPt(const MOMS &pbs)
Sort a container of momenta by pT (decreasing) and return by value for const inputs.
Definition Vector4.hh:1449
MOMS sortByEt(const MOMS &pbs)
Sort a container of momenta by Et (decreasing) and return by value for const inputs.
Definition Vector4.hh:1471
bool cmpMomByEt(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by decreasing transverse energy.
Definition Vector4.hh:1359
bool cmpMomByAscPt(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by increasing pT.
Definition Vector4.hh:1345
bool cmpMomByAscP(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by increasing 3-momentum magnitude |p|.
Definition Vector4.hh:1354
bool cmpMomByAscMass(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by increasing mass.
Definition Vector4.hh:1381
MOMS sortByE(const MOMS &pbs)
Sort a container of momenta by E (decreasing) and return by value for const inputs.
Definition Vector4.hh:1460
bool cmpMomByPt(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by decreasing pT.
Definition Vector4.hh:1341
bool cmpMomByAbsEta(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by increasing absolute eta (pseudorapidity).
Definition Vector4.hh:1396
bool cmpMomByDescAbsEta(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by increasing absolute eta (pseudorapidity).
Definition Vector4.hh:1401
bool cmpMomByEta(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by increasing eta (pseudorapidity).
Definition Vector4.hh:1386
bool cmpMomByAscE(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by increasing energy.
Definition Vector4.hh:1372
bool cmpMomByDescAbsRap(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by decreasing absolute rapidity.
Definition Vector4.hh:1421
MOMS & isortByPt(MOMS &pbs)
Sort a container of momenta by pT (decreasing) and return by reference for non-const inputs.
Definition Vector4.hh:1444
bool cmpMomByDescRap(const FourMomentum &a, const FourMomentum &b)
Comparison to give a sorting by decreasing rapidity.
Definition Vector4.hh:1411
double mass(const FourMomentum &a, const FourMomentum &b)
Calculate mass of two 4-vectors.
Definition Vector4.hh:1482
double mass2(const FourMomentum &a, const FourMomentum &b)
Calculate mass^2 of two 4-vectors.
Definition Vector4.hh:1487
double pT(const Vector3 &a, const Vector3 &b)
Calculate transverse momentum of pair of 3-vectors.
Definition Vector3.hh:649
std::vector< FourVector > FourVectors
Definition Vector4.hh:1586
double E(const ParticleBase &p)
Unbound function access to E.
Definition ParticleBaseUtils.hh:659
string to_str(const T &x)
Convert any object to a string.
Definition Utils.hh:78
Definition MC_CENT_PPB_Projections.hh:10
constexpr std::enable_if_t< std::is_arithmetic_v< NUM >, int > sign(NUM val)
Find the sign of a number.
Definition MathUtils.hh:275
double deltaR(double rap1, double phi1, double rap2, double phi2)
Definition MathUtils.hh:708
double deltaPhi(double phi1, double phi2, bool sign=false)
Calculate the difference between two angles in radians.
Definition MathUtils.hh:678
double deltaEta(double eta1, double eta2, bool sign=false)
Definition MathUtils.hh:686
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:701
double mT(double pT1, double pT2, double dphi)
Definition MathUtils.hh:731
std::ostream & operator<<(std::ostream &os, const AnalysisInfo &ai)
Stream an AnalysisInfo as a text description.
Definition AnalysisInfo.hh:362
std::enable_if_t< std::is_floating_point_v< NUM >, bool > isZero(NUM val, double tolerance=1e-8)
Compare a number to zero.
Definition MathUtils.hh:24
std::enable_if_t< std::is_arithmetic_v< NUM >, NUM > sqr(NUM a)
Named number-type squaring operation.
Definition MathUtils.hh:218
double invariant(const FourVector &lv)
Definition Vector4.hh:303
double contract(const FourVector &a, const FourVector &b)
Contract two 4-vectors, with metric signature (+ - - -).
Definition Vector4.hh:256
double add(double a, double b, double tolerance=1e-5)
Add two numbers with FP fuzziness.
Definition MathUtils.hh:229
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:694
Cmp< T > cmp(const T &t1, const T &t2)
Global helper function for easy creation of Cmp objects.
Definition Cmp.hh:255
std::string toString(const AnalysisInfo &ai)
String representation.
double angle(const Vector2 &a, const Vector2 &b)
Angle (in radians) between two 2-vectors.
Definition Vector2.hh:177
double rapidity(double E, double pz)
Calculate a rapidity value from the supplied energy E and longitudinal momentum pz.
Definition MathUtils.hh:713