1#ifndef RIVET_PARTICLEBASEUTILS_HH
2#define RIVET_PARTICLEBASEUTILS_HH
4#include "Rivet/Tools/Utils.hh"
5#include "Rivet/ParticleBase.hh"
34 PtGtr(
double pt) : ptcut(pt) { }
36 bool operator()(
const ParticleBase&
p)
const {
return p.pT() > ptcut; }
45 PtLess(
double pt) : ptcut(pt) { }
46 bool operator()(
const ParticleBase&
p)
const {
return p.pT() < ptcut; }
54 PtInRange(pair<double, double> ptcuts) : ptcut(ptcuts) { }
55 PtInRange(
double ptlow,
double pthigh) : PtInRange(make_pair(ptlow, pthigh)) { }
57 bool operator()(
const ParticleBase&
p)
const {
return p.pT() >= ptcut.first &&
p.pT() < ptcut.second; }
58 pair<double,double> ptcut;
66 EtaGtr(
double eta) : etacut(
eta) { }
68 bool operator()(
const ParticleBase&
p)
const {
return p.eta() > etacut; }
75 EtaLess(
double eta) : etacut(
eta) { }
77 bool operator()(
const ParticleBase&
p)
const {
return p.eta() < etacut; }
84 EtaInRange(pair<double, double> etacuts) : etacut(etacuts) { }
85 EtaInRange(
double etalow,
double etahigh) : EtaInRange(make_pair(etalow, etahigh)) { }
87 bool operator()(
const ParticleBase&
p)
const {
return p.eta() >= etacut.first &&
p.eta() < etacut.second; }
88 pair<double,double> etacut;
97 bool operator()(
const ParticleBase&
p)
const {
return p.abseta() > absetacut; }
107 bool operator()(
const ParticleBase&
p)
const {
return p.abseta() < absetacut; }
115 AbsEtaInRange(
const pair<double, double>& absetacuts) : absetacut(absetacuts) { }
116 AbsEtaInRange(
double absetalow,
double absetahigh) : AbsEtaInRange(make_pair(absetalow, absetahigh)) { }
118 bool operator()(
const ParticleBase&
p)
const {
return p.abseta() >= absetacut.first &&
p.abseta() < absetacut.second; }
119 pair<double,double> absetacut;
127 RapGtr(
double rap) : rapcut(
rap) { }
129 bool operator()(
const ParticleBase&
p)
const {
return p.rap() > rapcut; }
136 RapLess(
double rap) : rapcut(
rap) { }
138 bool operator()(
const ParticleBase&
p)
const {
return p.rap() < rapcut; }
145 RapInRange(
const pair<double, double>& rapcuts) : rapcut(rapcuts) { }
146 RapInRange(
double raplow,
double raphigh) : RapInRange(make_pair(raplow, raphigh)) { }
148 bool operator()(
const ParticleBase&
p)
const {
return p.rap() >= rapcut.first &&
p.rap() < rapcut.second; }
149 pair<double,double> rapcut;
158 bool operator()(
const ParticleBase&
p)
const {
return p.absrap() > absrapcut; }
168 bool operator()(
const ParticleBase&
p)
const {
return p.absrap() < absrapcut; }
176 AbsRapInRange(
const pair<double, double>& absrapcuts) : absrapcut(absrapcuts) { }
177 AbsRapInRange(
double absraplow,
double absraphigh) : AbsRapInRange(make_pair(absraplow, absraphigh)) { }
179 bool operator()(
const ParticleBase&
p)
const {
return p.absrap() >= absrapcut.first &&
p.absrap() < absrapcut.second; }
180 pair<double,double> absrapcut;
192 : refvec(vec.
mom()), drcut(dr), rapscheme(scheme) { }
194 : refvec(vec), drcut(dr), rapscheme(scheme) { }
195 DeltaRGtr(
const Vector3& vec,
double dr)
196 : drcut(dr), rapscheme(PSEUDORAPIDITY) { refvec.setPx(vec.x()); refvec.setPy(vec.y()); refvec.setPz(vec.z()); }
207 : refvec(vec.
mom()), drcut(dr), rapscheme(scheme) { }
209 : refvec(vec), drcut(dr), rapscheme(scheme) { }
210 DeltaRLess(
const Vector3& vec,
double dr)
211 : drcut(dr), rapscheme(PSEUDORAPIDITY) { refvec.setPx(vec.x()); refvec.setPy(vec.y()); refvec.setPz(vec.z()); }
221 DeltaRInRange(
const ParticleBase& vec,
const pair<double,double>& dr,
RapScheme scheme=PSEUDORAPIDITY)
222 : refvec(vec.
mom()), drcut(dr), rapscheme(scheme) { }
224 : DeltaRInRange(vec, make_pair(drmin, drmax), scheme) { }
225 DeltaRInRange(
const FourMomentum& vec,
const pair<double,double>& dr,
RapScheme scheme=PSEUDORAPIDITY)
226 : refvec(vec), drcut(dr), rapscheme(scheme) { }
228 : DeltaRInRange(vec, make_pair(drmin, drmax), scheme) { }
229 DeltaRInRange(
const Vector3& vec,
const pair<double,double>& dr)
230 : drcut(dr), rapscheme(PSEUDORAPIDITY) { refvec.setPx(vec.x()); refvec.setPy(vec.y()); refvec.setPz(vec.z()); }
231 DeltaRInRange(
const Vector3& vec,
double drmin,
double drmax)
232 : DeltaRInRange(vec, make_pair(drmin, drmax)) { }
234 const double dR =
deltaR(
p, refvec, rapscheme);
235 return dR >= drcut.first && dR < drcut.second;
238 pair<double,double> drcut;
247 : refvec(vec.
p3()), dphicut(dphi) { }
249 : refvec(vec.
p3()), dphicut(dphi) { }
250 DeltaPhiGtr(
const Vector3& vec,
double dphi)
251 : refvec(vec), dphicut(dphi) { }
261 : refvec(vec.
p3()), dphicut(dphi) { }
263 : refvec(vec.
p3()), dphicut(dphi) { }
264 DeltaPhiLess(
const Vector3& vec,
double dphi)
265 : refvec(vec), dphicut(dphi) { }
274 DeltaPhiInRange(
const ParticleBase& vec,
const pair<double,double>& dphi)
275 : refvec(vec.
mom()), dphicut(dphi) { }
276 DeltaPhiInRange(
const ParticleBase& vec,
double dphimin,
double dphimax)
277 : DeltaPhiInRange(vec, make_pair(dphimin, dphimax)) { }
278 DeltaPhiInRange(
const FourMomentum& vec,
const pair<double,double>& dphi)
279 : refvec(vec), dphicut(dphi) { }
280 DeltaPhiInRange(
const FourMomentum& vec,
double dphimin,
double dphimax)
281 : DeltaPhiInRange(vec, make_pair(dphimin, dphimax)) { }
282 DeltaPhiInRange(
const Vector3& vec,
const pair<double,double>& dphi)
283 : refvec(vec), dphicut(dphi) { }
284 DeltaPhiInRange(
const Vector3& vec,
double dphimin,
double dphimax)
285 : DeltaPhiInRange(vec, make_pair(dphimin, dphimax)) { }
288 return dphi >= dphicut.first && dphi < dphicut.second;
291 pair<double,double> dphicut;
299 : refvec(vec.
p3()), detacut(deta) { }
301 : refvec(vec.
p3()), detacut(deta) { }
302 DeltaEtaGtr(
const Vector3& vec,
double deta)
303 : refvec(vec), detacut(deta) { }
313 : refvec(vec.
p3()), detacut(deta) { }
315 : refvec(vec.
p3()), detacut(deta) { }
316 DeltaEtaLess(
const Vector3& vec,
double deta)
317 : refvec(vec), detacut(deta) { }
326 DeltaEtaInRange(
const ParticleBase& vec,
const pair<double,double>& deta)
327 : refvec(vec.
mom()), detacut(deta) { }
328 DeltaEtaInRange(
const ParticleBase& vec,
double detamin,
double detamax)
329 : DeltaEtaInRange(vec, make_pair(detamin, detamax)) { }
330 DeltaEtaInRange(
const FourMomentum& vec,
const pair<double,double>& deta)
331 : refvec(vec), detacut(deta) { }
332 DeltaEtaInRange(
const FourMomentum& vec,
double detamin,
double detamax)
333 : DeltaEtaInRange(vec, make_pair(detamin, detamax)) { }
334 DeltaEtaInRange(
const Vector3& vec,
const pair<double,double>& deta)
335 : refvec(vec), detacut(deta) { }
336 DeltaEtaInRange(
const Vector3& vec,
double detamin,
double detamax)
337 : DeltaEtaInRange(vec, make_pair(detamin, detamax)) { }
340 return deta >= detacut.first && deta < detacut.second;
343 pair<double,double> detacut;
351 : refvec(vec.
mom()), drapcut(drap) { }
353 : refvec(vec), drapcut(drap) { }
363 : refvec(vec.
mom()), drapcut(drap) { }
365 : refvec(vec), drapcut(drap) { }
374 DeltaRapInRange(
const ParticleBase& vec,
const pair<double,double>& drap)
375 : refvec(vec.
mom()), drapcut(drap) { }
376 DeltaRapInRange(
const ParticleBase& vec,
double drapmin,
double drapmax)
377 : DeltaRapInRange(vec, make_pair(drapmin, drapmax)) { }
378 DeltaRapInRange(
const FourMomentum& vec,
const pair<double,double>& drap)
379 : refvec(vec), drapcut(drap) { }
380 DeltaRapInRange(
const FourMomentum& vec,
double drapmin,
double drapmax)
381 : DeltaRapInRange(vec, make_pair(drapmin, drapmax)) { }
384 return drap >= drapcut.first && drap < drapcut.second;
387 pair<double,double> drapcut;
410 DeltaRWRT(
const Vector3&
p3) : p(
p3.mod(),
p3.x(),
p3.y(),
p3.z()), rapscheme(PSEUDORAPIDITY) {}
481 template<
typename PBCONTAINER1,
typename PBCONTAINER2>
482 inline void idiscardIfAny(PBCONTAINER1& tofilter,
const PBCONTAINER2& tocompare,
483 typename std::function<
bool(
const typename PBCONTAINER1::value_type&,
484 const typename PBCONTAINER2::value_type&)> fn) {
485 for (
const auto& pbcmp : tocompare) {
486 idiscard(tofilter, [&](
const typename PBCONTAINER1::value_type& pbfilt){
return fn(pbfilt, pbcmp); });
490 template<
typename PBCONTAINER1,
typename PBCONTAINER2>
491 inline PBCONTAINER1 discardIfAny(
const PBCONTAINER1& tofilter,
const PBCONTAINER2& tocompare,
492 typename std::function<
bool(
const typename PBCONTAINER1::value_type&,
493 const typename PBCONTAINER2::value_type&)> fn) {
494 PBCONTAINER1 tmp{tofilter};
495 idiscardIfAny(tmp, tocompare, fn);
500 template<
typename PBCONTAINER1,
typename PBCONTAINER2>
501 inline PBCONTAINER1 selectIfAny(
const PBCONTAINER1& tofilter,
const PBCONTAINER2& tocompare,
502 typename std::function<
bool(
const typename PBCONTAINER1::value_type&,
503 const typename PBCONTAINER2::value_type&)> fn) {
504 PBCONTAINER1 selected;
505 for (
const auto& pbfilt : tofilter) {
506 if (
any(tocompare, [&](
const typename PBCONTAINER2::value_type& pbcmp){
return fn(pbfilt, pbcmp); })) {
513 template<
typename PBCONTAINER1,
typename PBCONTAINER2>
514 inline void iselectIfAny(PBCONTAINER1& tofilter,
const PBCONTAINER2& tocompare,
515 typename std::function<
bool(
const typename PBCONTAINER1::value_type&,
516 const typename PBCONTAINER2::value_type&)> fn) {
517 tofilter = selectIfAny(tofilter, tocompare, fn);
522 template<
typename PBCONTAINER1,
typename PBCONTAINER2>
523 inline PBCONTAINER1 discardIfAll(
const PBCONTAINER1& tofilter,
const PBCONTAINER2& tocompare,
524 typename std::function<
bool(
const typename PBCONTAINER1::value_type&,
525 const typename PBCONTAINER2::value_type&)> fn) {
526 PBCONTAINER1 selected;
527 for (
const auto& pbfilt : tofilter) {
528 if (!
all(tocompare, [&](
const typename PBCONTAINER2::value_type& pbcmp){
return fn(pbfilt, pbcmp); })) {
535 template<
typename PBCONTAINER1,
typename PBCONTAINER2>
536 inline void idiscardIfAll(PBCONTAINER1& tofilter,
const PBCONTAINER2& tocompare,
537 typename std::function<
bool(
const typename PBCONTAINER1::value_type&,
538 const typename PBCONTAINER2::value_type&)> fn) {
539 tofilter = discardIfAll(tofilter, tocompare, fn);
543 template<
typename PBCONTAINER1,
typename PBCONTAINER2>
544 inline PBCONTAINER1 selectIfAll(
const PBCONTAINER1& tofilter,
const PBCONTAINER2& tocompare,
545 typename std::function<
bool(
const typename PBCONTAINER1::value_type&,
546 const typename PBCONTAINER2::value_type&)> fn) {
547 PBCONTAINER1 selected;
548 for (
const auto& pbfilt : tofilter) {
549 if (
all(tocompare, [&](
const typename PBCONTAINER2::value_type& pbcmp){
return fn(pbfilt, pbcmp); })) {
556 template<
typename PBCONTAINER1,
typename PBCONTAINER2>
557 inline void iselectIfAll(PBCONTAINER1& tofilter,
const PBCONTAINER2& tocompare,
558 typename std::function<
bool(
const typename PBCONTAINER1::value_type&,
559 const typename PBCONTAINER2::value_type&)> fn) {
560 tofilter = selectIfAll(tofilter, tocompare, fn);
569 template<
typename PBCONTAINER1,
typename PBCONTAINER2>
570 inline void idiscardIfAnyDeltaRLess(PBCONTAINER1& tofilter,
const PBCONTAINER2& tocompare,
double dR) {
571 for (
const typename PBCONTAINER2::value_type& pb : tocompare) {
572 idiscard(tofilter, deltaRLess(pb, dR));
576 template<
typename PBCONTAINER1,
typename PBCONTAINER2>
577 inline PBCONTAINER1 discardIfAnyDeltaRLess(
const PBCONTAINER1& tofilter,
const PBCONTAINER2& tocompare,
double dR) {
578 PBCONTAINER1 tmp{tofilter};
579 idiscardIfAnyDeltaRLess(tmp, tocompare, dR);
583 template<
typename PBCONTAINER1,
typename PBCONTAINER2>
584 inline void idiscardIfAnyDeltaPhiLess(PBCONTAINER1& tofilter,
const PBCONTAINER2& tocompare,
double dphi) {
585 for (
const typename PBCONTAINER2::value_type& pb : tocompare) {
586 idiscard(tofilter, deltaPhiLess(pb, dphi));
590 template<
typename PBCONTAINER1,
typename PBCONTAINER2>
591 inline PBCONTAINER1 discardIfAnyDeltaPhiLess(
const PBCONTAINER1& tofilter,
const PBCONTAINER2& tocompare,
double dphi) {
592 PBCONTAINER1 tmp{tofilter};
593 idiscardIfAnyDeltaPhiLess(tmp, tocompare, dphi);
599 template<
typename PBCONTAINER1,
typename PBCONTAINER2>
600 inline PBCONTAINER1 selectIfAnyDeltaRLess(
const PBCONTAINER1& tofilter,
const PBCONTAINER2& tocompare,
double dR) {
601 PBCONTAINER1 selected;
602 for (
const typename PBCONTAINER1::value_type& f : tofilter) {
603 if (
any(tocompare, deltaRLess(f, dR))) selected.push_back(f);
608 template<
typename PBCONTAINER1,
typename PBCONTAINER2>
609 inline void iselectIfAnyDeltaRLess(PBCONTAINER1& tofilter,
const PBCONTAINER2& tocompare,
double dR) {
610 tofilter = selectIfAnyDeltaRLess(tofilter, tocompare, dR);
614 template<
typename PBCONTAINER1,
typename PBCONTAINER2>
615 inline PBCONTAINER1 selectIfAnyDeltaPhiLess(
const PBCONTAINER1& tofilter,
const PBCONTAINER2& tocompare,
double dphi) {
616 PBCONTAINER1 selected;
617 for (
const typename PBCONTAINER1::value_type& f : tofilter) {
618 if (
any(tocompare, deltaPhiLess(f, dphi))) selected.push_back(f);
623 template<
typename PBCONTAINER1,
typename PBCONTAINER2>
624 inline void iselectIfAnyDeltaPhiLess(PBCONTAINER1& tofilter,
const PBCONTAINER2& tocompare,
double dphi) {
625 tofilter = selectIfAnyDeltaPhiLess(tofilter, tocompare, dphi);
722 return mT(
p.mom(),
p4);
730 return mT(
p4,
p.mom());
746 return mT(p1.
mom(), p2);
754 return mT(p1, p2.
mom());
762 return mT(p1.
mom(), p2);
770 return mT(p1, p2.
mom());
775 return pT(
p.mom(),
p4);
780 return pT(
p4,
p.mom());
804 template <
typename CONTAINER,
typename = isCIterable<CONTAINER>>
805 inline int closestMassIndex(
const CONTAINER& c,
double mtarget,
double mmin=-DBL_MAX,
double mmax=DBL_MAX) {
813 template <
typename CONTAINER1,
typename CONTAINER2,
typename = isCIterable<CONTAINER1, CONTAINER2>>
815 double mtarget,
double mmin=-DBL_MAX,
double mmax=DBL_MAX) {
823 template <
typename CONTAINER,
typename T,
typename = isCIterable<CONTAINER>>
825 double mtarget,
double mmin=-DBL_MAX,
double mmax=DBL_MAX) {
835 template <
typename CONTAINER,
typename T,
typename = isCIterable<CONTAINER>>
837 double mtarget,
double mmin=-DBL_MAX,
double mmax=DBL_MAX) {
848 template <
typename CONTAINER,
typename = isCIterable<CONTAINER>>
849 inline double sumPt(
const CONTAINER& c) {
854 template <
typename CONTAINER,
typename = isCIterable<CONTAINER>>
855 inline double sumE(
const CONTAINER& c) {
860 template <
typename CONTAINER,
typename = isCIterable<CONTAINER>>
861 inline double sumEt(
const CONTAINER& c) {
866 template <
typename CONTAINER,
typename = isCIterable<CONTAINER>>
872 template <
typename CONTAINER,
typename = isCIterable<CONTAINER>>
878 template <
typename CONTAINER,
typename = isCIterable<CONTAINER>>
884 template <
typename CONTAINER,
typename = isCIterable<CONTAINER>>
895 template <
typename CONTAINER,
typename = isCIterable<CONTAINER>>
897 double rtn = DBL_MAX;
898 for (
size_t i = 0; i < c.size(); ++i) {
899 for (
size_t j = i; j < c.size(); ++j) {
900 if (i == j)
continue;
901 const double dphi = fabs(
deltaPhi(c[i], c[j]));
902 if (dphi < rtn) rtn = dphi;
915 template <
typename CONTAINER1,
typename CONTAINER2,
typename = isCIterable<CONTAINER1, CONTAINER2>>
916 inline double deltaPhiMin(
const CONTAINER1& c1,
const CONTAINER2& c2) {
917 double rtn = DBL_MAX;
918 for (
size_t i = 0; i < c1.size(); ++i) {
919 for (
size_t j = 0; j < c2.size(); ++j) {
920 const double dphi =
deltaPhi(c1[i], c2[j]);
921 if (dphi < rtn) rtn = dphi;
934 template <
typename CONTAINER,
typename = isCIterable<CONTAINER>>
936 double rtn = DBL_MAX;
937 for (
size_t i = 0; i < c.size(); ++i) {
938 for (
size_t j = i; j < c.size(); ++j) {
939 if (i == j)
continue;
940 const double deta = fabs(
deltaEta(c[i], c[j]));
941 if (deta < rtn) rtn = deta;
954 template <
typename CONTAINER1,
typename CONTAINER2,
typename = isCIterable<CONTAINER1, CONTAINER2>>
955 inline double deltaEtaMin(
const CONTAINER1& c1,
const CONTAINER2& c2) {
956 double rtn = DBL_MAX;
957 for (
size_t i = 0; i < c1.size(); ++i) {
958 for (
size_t j = 0; j < c2.size(); ++j) {
959 const double deta =
deltaEta(c1[i], c2[j]);
960 if (deta < rtn) rtn = deta;
972 template <
typename CONTAINER,
typename = isCIterable<CONTAINER>>
974 double rtn = DBL_MAX;
975 for (
size_t i = 0; i < c.size(); ++i) {
976 for (
size_t j = i; j < c.size(); ++j) {
977 if (i == j)
continue;
978 const double dy = fabs(
deltaRap(c[i], c[j]));
979 if (dy < rtn) rtn = dy;
992 template <
typename CONTAINER1,
typename CONTAINER2,
typename = isCIterable<CONTAINER1, CONTAINER2>>
993 inline double deltaRapMin(
const CONTAINER1& c1,
const CONTAINER2& c2) {
994 double rtn = DBL_MAX;
995 for (
size_t i = 0; i < c1.size(); ++i) {
996 for (
size_t j = 0; j < c2.size(); ++j) {
997 const double dy =
deltaRap(c1[i], c2[j]);
998 if (dy < rtn) rtn = dy;
1011 template <
typename CONTAINER,
typename = isCIterable<CONTAINER>>
1013 double rtn = DBL_MAX;
1014 for (
size_t i = 0; i < c.size(); ++i) {
1015 for (
size_t j = i; j < c.size(); ++j) {
1016 if (i == j)
continue;
1017 const double dr = fabs(
deltaR(c[i], c[j]));
1018 if (dr < rtn) rtn = dr;
1031 template <
typename CONTAINER1,
typename CONTAINER2,
typename = isCIterable<CONTAINER1, CONTAINER2>>
1032 inline double deltaRMin(
const CONTAINER1& c1,
const CONTAINER2& c2) {
1033 double rtn = DBL_MAX;
1034 for (
size_t i = 0; i < c1.size(); ++i) {
1035 for (
size_t j = 0; j < c2.size(); ++j) {
1036 const double dr =
deltaR(c1[i], c2[j]);
1037 if (dr < rtn) rtn = dr;
1045 template <
typename CONTAINER,
typename PVEC,
typename = isCIterable<CONTAINER>>
1046 inline double mTmin(
const CONTAINER& psvis,
const PVEC& pinvis) {
1047 return min(transform(psvis, [&](
const typename CONTAINER::value_type& pvis) {
return mT(pvis, pinvis); }));
Specialized version of the FourVector with momentum/energy functionality.
Definition Vector4.hh:327
double rap() const
Alias for rapidity.
Definition Vector4.hh:632
double absrap() const
Absolute rapidity.
Definition Vector4.hh:641
double pT() const
Calculate the transverse momentum .
Definition Vector4.hh:664
Vector3 vector3() const
Get the spatial part of the 4-vector as a 3-vector.
Definition Vector4.hh:184
double abseta() const
Get the directly (alias).
Definition Vector4.hh:181
double eta() const
Synonym for pseudorapidity.
Definition Vector4.hh:174
Base class for particle-like things like Particle and Jet.
Definition ParticleBase.hh:13
const FourMomentum & mom() const
Get the equivalent momentum four-vector (const) (alias).
Definition ParticleBase.hh:39
ThreeMomentum p3() const
Get the 3-momentum directly.
Definition ParticleBase.hh:110
Three-dimensional specialisation of Vector.
Definition Vector3.hh:40
int closestMatchIndex(const CONTAINER &c, FN &&fn, double target, double minval=-DBL_MAX, double maxval=DBL_MAX)
Return the index from a vector which best matches fn(c[i]) to the target value.
Definition Utils.hh:786
pair< int, int > closestMatchIndices(const CONTAINER1 &c1, const CONTAINER2 &c2, FN &&fn, double target, double minval=-DBL_MAX, double maxval=DBL_MAX)
Return the indices from two vectors which best match fn(c1[i], c2[j]) to the target value.
Definition Utils.hh:815
bool any(const CONTAINER &c)
Return true if x is true for any x in container c, otherwise false.
Definition Utils.hh:358
bool all(const CONTAINER &c)
Return true if x is true for all x in container c, otherwise false.
Definition Utils.hh:383
Jets & idiscard(Jets &jets, const Cut &c)
Filter a jet collection in-place to the subset that fails the supplied Cut.
function< bool(const ParticleBase &)> ParticleBaseSelector
std::function instantiation for functors taking a ParticleBase and returning a bool
Definition ParticleBaseUtils.hh:20
function< bool(const ParticleBase &, const ParticleBase &)> ParticleBaseSorter
std::function instantiation for functors taking two ParticleBase and returning a bool
Definition ParticleBaseUtils.hh:22
double pT(const ParticleBase &p)
Unbound function access to pT.
Definition ParticleBaseUtils.hh:656
FourMomentum p4(const ParticleBase &p)
Unbound function access to momentum.
Definition ParticleBaseUtils.hh:644
Vector3 pTvec(const ParticleBase &p)
Unbound function access to pTvec.
Definition ParticleBaseUtils.hh:650
Vector3 p3(const ParticleBase &p)
Unbound function access to p3.
Definition ParticleBaseUtils.hh:647
FourMomentum mom(const ParticleBase &p)
Unbound function access to momentum.
Definition ParticleBaseUtils.hh:642
double absrap(const ParticleBase &p)
Unbound function access to abs rapidity.
Definition ParticleBaseUtils.hh:674
double E(const ParticleBase &p)
Unbound function access to E.
Definition ParticleBaseUtils.hh:659
double eta(const ParticleBase &p)
Unbound function access to eta.
Definition ParticleBaseUtils.hh:665
double mass(const ParticleBase &p)
Unbound function access to mass.
Definition ParticleBaseUtils.hh:677
double mass2(const ParticleBase &p, const P4 &p4)
Get the mass^2 of a ParticleBase and a P4.
Definition ParticleBaseUtils.hh:703
double rap(const ParticleBase &p)
Unbound function access to rapidity.
Definition ParticleBaseUtils.hh:671
double Et(const ParticleBase &p)
Unbound function access to ET.
Definition ParticleBaseUtils.hh:662
double p(const ParticleBase &p)
Unbound function access to p.
Definition ParticleBaseUtils.hh:653
double abseta(const ParticleBase &p)
Unbound function access to abseta.
Definition ParticleBaseUtils.hh:668
double mT(const ParticleBase &p, const P4 &p4)
Get the transverse mass of a ParticleBase and a P4.
Definition ParticleBaseUtils.hh:721
int closestMassIndex(const CONTAINER &c, double mtarget, double mmin=-DBL_MAX, double mmax=DBL_MAX)
Return the index from a vector which best matches mass(c[i]) to the target value.
Definition ParticleBaseUtils.hh:805
pair< int, int > closestMassIndices(const CONTAINER1 &c1, const CONTAINER2 &c2, double mtarget, double mmin=-DBL_MAX, double mmax=DBL_MAX)
Return the indices from two vectors which best match mass(c1[i], c2[j]) to the target value.
Definition ParticleBaseUtils.hh:814
double mTmin(const CONTAINER &psvis, const PVEC &pinvis)
Get the minimum transverse mass of a pair of vectors, one of them taken from a list.
Definition ParticleBaseUtils.hh:1046
double sumEt(const CONTAINER &c)
Return the ET sum of a collection of particle-like objects.
Definition ParticleBaseUtils.hh:861
double deltaEtaMin(const CONTAINER &c)
Calculate the minimum pseudorapidity separation in a collection of particle-like objects.
Definition ParticleBaseUtils.hh:935
double deltaRapMin(const CONTAINER &c)
Calculate the minimum rapidity separation in a collection of particle-like objects.
Definition ParticleBaseUtils.hh:973
double sumE(const CONTAINER &c)
Return the energy sum of a collection of particle-like objects.
Definition ParticleBaseUtils.hh:855
double deltaRMin(const CONTAINER &c)
Calculate the minimum rapidity separation in a collection of particle-like objects.
Definition ParticleBaseUtils.hh:1012
double sumPt(const CONTAINER &c)
Return the pT sum of a collection of particle-like objects.
Definition ParticleBaseUtils.hh:849
Vector3 sumPtVec(const CONTAINER &c)
Return the four-momentum sum of a collection of particle-like objects.
Definition ParticleBaseUtils.hh:885
FourMomentum sumP4(const CONTAINER &c)
Return the four-momentum sum of a collection of particle-like objects.
Definition ParticleBaseUtils.hh:873
double sumMass(const CONTAINER &c)
Return the mass sum of a collection of particle-like objects.
Definition ParticleBaseUtils.hh:867
double deltaPhiMin(const CONTAINER &c)
Calculate the minimum phi separation in a collection of particle-like objects.
Definition ParticleBaseUtils.hh:896
Vector3 sumP3(const CONTAINER &c)
Return the four-momentum sum of a collection of particle-like objects.
Definition ParticleBaseUtils.hh:879
Definition MC_CENT_PPB_Projections.hh:10
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
double mT(double pT1, double pT2, double dphi)
Definition MathUtils.hh:731
T sum(const DressedLeptons &c, FN &&fn, const T &start=T())
Generic sum function, adding fn(x) for all x in container c, starting with start.
Definition DressedLepton.hh:60
std::enable_if_t< std::is_arithmetic_v< N1 > &&std::is_arithmetic_v< N2 >, signed_if_mixed_t< N1, N2 > > min(N1 a, N2 b)
Get the minimum of two numbers.
Definition MathUtils.hh:104
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
Calculator of with respect to a given momentum.
Definition ParticleBaseUtils.hh:444
Calculator of with respect to a given momentum.
Definition ParticleBaseUtils.hh:466
Abs pseudorapidity greater-than functor.
Definition ParticleBaseUtils.hh:94
Abs pseudorapidity in-range functor.
Definition ParticleBaseUtils.hh:114
Abs pseudorapidity momentum less-than functor.
Definition ParticleBaseUtils.hh:104
Abs rapidity greater-than functor.
Definition ParticleBaseUtils.hh:155
Abs rapidity in-range functor.
Definition ParticleBaseUtils.hh:175
Abs rapidity momentum less-than functor.
Definition ParticleBaseUtils.hh:165
Base type for Particle -> bool functors.
Definition ParticleBaseUtils.hh:26
(with respect to another momentum, vec) greater-than functor
Definition ParticleBaseUtils.hh:297
(with respect to another 4-momentum, vec) in-range functor
Definition ParticleBaseUtils.hh:325
(with respect to another momentum, vec) less-than functor
Definition ParticleBaseUtils.hh:311
Calculator of with respect to a given momentum.
Definition ParticleBaseUtils.hh:432
(with respect to another momentum, vec) greater-than functor
Definition ParticleBaseUtils.hh:245
(with respect to another 4-momentum, vec) in-range functor
Definition ParticleBaseUtils.hh:273
(with respect to another momentum, vec) less-than functor
Definition ParticleBaseUtils.hh:259
Calculator of with respect to a given momentum.
Definition ParticleBaseUtils.hh:420
(with respect to another 4-momentum, vec) greater-than functor
Definition ParticleBaseUtils.hh:190
(with respect to another 4-momentum, vec) in-range functor
Definition ParticleBaseUtils.hh:220
(with respect to another 4-momentum, vec) less-than functor
Definition ParticleBaseUtils.hh:205
Calculator of with respect to a given momentum.
Definition ParticleBaseUtils.hh:407
(with respect to another momentum, vec) greater-than functor
Definition ParticleBaseUtils.hh:349
(with respect to another 4-momentum, vec) in-range functor
Definition ParticleBaseUtils.hh:373
(with respect to another momentum, vec) less-than functor
Definition ParticleBaseUtils.hh:361
Calculator of with respect to a given momentum.
Definition ParticleBaseUtils.hh:456
Base type for Particle -> double functors.
Definition ParticleBaseUtils.hh:401
Pseudorapidity greater-than functor.
Definition ParticleBaseUtils.hh:65
Pseudorapidity in-range functor.
Definition ParticleBaseUtils.hh:83
Pseudorapidity less-than functor.
Definition ParticleBaseUtils.hh:74
Transverse momentum greater-than functor.
Definition ParticleBaseUtils.hh:33
Transverse momentum in-range functor.
Definition ParticleBaseUtils.hh:53
Transverse momentum less-than functor.
Definition ParticleBaseUtils.hh:43
Rapidity greater-than functor.
Definition ParticleBaseUtils.hh:126
Rapidity in-range functor.
Definition ParticleBaseUtils.hh:144
Rapidity momentum less-than functor.
Definition ParticleBaseUtils.hh:135