|
|
Rivet 4.0.0
|
1#ifndef RIVET_PARTICLEBASEUTILS_HH
2#define RIVET_PARTICLEBASEUTILS_HH
4#include "Rivet/Tools/Utils.hh"
5#include "Rivet/ParticleBase.hh"
27 virtual bool operator()( const ParticleBase& p) const = 0;
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) { }
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) { }
87 bool operator()( const ParticleBase& p) const { return p.eta() >= etacut.first && p.eta() < etacut.second; }
88 pair<double,double> etacut;
95 AbsEtaGtr( double abseta) : absetacut(abseta) { }
97 bool operator()( const ParticleBase& p) const { return p.abseta() > absetacut; }
105 AbsEtaLess( double abseta) : absetacut(abseta) { }
107 bool operator()( const ParticleBase& p) const { return p.abseta() < absetacut; }
115 AbsEtaInRange( const pair<double, double>& absetacuts) : absetacut(absetacuts) { }
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) { }
148 bool operator()( const ParticleBase& p) const { return p.rap() >= rapcut.first && p.rap() < rapcut.second; }
149 pair<double,double> rapcut;
156 AbsRapGtr( double absrap) : absrapcut(absrap) { }
158 bool operator()( const ParticleBase& p) const { return p.absrap() > absrapcut; }
166 AbsRapLess( double absrap) : absrapcut(absrap) { }
168 bool operator()( const ParticleBase& p) const { return p.absrap() < absrapcut; }
176 AbsRapInRange( const pair<double, double>& absrapcuts) : absrapcut(absrapcuts) { }
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) { }
196 : drcut(dr), rapscheme(PSEUDORAPIDITY) { refvec. setPx(vec.x()); refvec. setPy(vec.y()); refvec. setPz(vec.z()); }
197 bool operator()( const ParticleBase& p) const { return deltaR(p, refvec, rapscheme) > drcut; }
207 : refvec(vec. mom()), drcut(dr), rapscheme(scheme) { }
209 : refvec(vec), drcut(dr), rapscheme(scheme) { }
211 : drcut(dr), rapscheme(PSEUDORAPIDITY) { refvec. setPx(vec.x()); refvec. setPy(vec.y()); refvec. setPz(vec.z()); }
212 bool operator()( const ParticleBase& p) const { return deltaR(p, refvec, rapscheme) < drcut; }
222 : refvec(vec. mom()), drcut(dr), rapscheme(scheme) { }
226 : refvec(vec), drcut(dr), rapscheme(scheme) { }
230 : drcut(dr), rapscheme(PSEUDORAPIDITY) { refvec. setPx(vec.x()); refvec. setPy(vec.y()); refvec. setPz(vec.z()); }
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) { }
251 : refvec(vec), dphicut(dphi) { }
261 : refvec(vec. p3()), dphicut(dphi) { }
263 : refvec(vec. p3()), dphicut(dphi) { }
265 : refvec(vec), dphicut(dphi) { }
275 : refvec(vec. mom()), dphicut(dphi) { }
279 : refvec(vec), dphicut(dphi) { }
283 : refvec(vec), dphicut(dphi) { }
287 const double dphi = deltaPhi(p, refvec);
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) { }
303 : refvec(vec), detacut(deta) { }
304 bool operator()( const ParticleBase& p) const { return std::abs( deltaEta(p, refvec)) > detacut; }
313 : refvec(vec. p3()), detacut(deta) { }
315 : refvec(vec. p3()), detacut(deta) { }
317 : refvec(vec), detacut(deta) { }
318 bool operator()( const ParticleBase& p) const { return std::abs( deltaEta(p, refvec)) < detacut; }
327 : refvec(vec. mom()), detacut(deta) { }
331 : refvec(vec), detacut(deta) { }
335 : refvec(vec), detacut(deta) { }
339 const double deta = deltaEta(p, refvec);
340 return deta >= detacut.first && deta < detacut.second;
343 pair<double,double> detacut;
351 : refvec(vec. mom()), drapcut(drap) { }
353 : refvec(vec), drapcut(drap) { }
354 bool operator()( const ParticleBase& p) const { return std::abs( deltaRap(p, refvec)) > drapcut; }
363 : refvec(vec. mom()), drapcut(drap) { }
365 : refvec(vec), drapcut(drap) { }
366 bool operator()( const ParticleBase& p) const { return std::abs( deltaRap(p, refvec)) < drapcut; }
375 : refvec(vec. mom()), drapcut(drap) { }
379 : refvec(vec), drapcut(drap) { }
383 const double drap = deltaRap(p, refvec);
384 return drap >= drapcut.first && drap < drapcut.second;
387 pair<double,double> drapcut;
402 virtual double operator()( const ParticleBase& p) const = 0;
410 DeltaRWRT( const Vector3& p3) : p(p3. mod(), p3.x(), p3.y(), p3.z()), rapscheme(PSEUDORAPIDITY) {}
413 double operator()( const Vector3& p3) const { return deltaR(p, p3); }
450 double operator()( const Vector3& p3) const { return fabs( deltaEta(p, p3)); }
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);
642 inline FourMomentum mom( const ParticleBase& p) { return p.mom(); }
644 inline FourMomentum p4( const ParticleBase& p) { return p.mom(); }
647 inline Vector3 p3( const ParticleBase& p) { return p.p3(); }
650 inline Vector3 pTvec( const ParticleBase& p) { return p.pTvec(); }
653 inline double p( const ParticleBase& p) { return p.p(); }
656 inline double pT( const ParticleBase& p) { return p.pT(); }
659 inline double Et( const ParticleBase& p) { return p.Et(); }
662 inline double eta( const ParticleBase& p) { return p.eta(); }
665 inline double abseta( const ParticleBase& p) { return p.abseta(); }
668 inline double rap( const ParticleBase& p) { return p.rap(); }
671 inline double absrap( const ParticleBase& p) { return p.absrap(); }
674 inline double mass( const ParticleBase& p) { return p.mass(); }
685 inline double mass( const ParticleBase& p, const FourMomentum& p4) {
686 return mass(p.mom(), p4);
690 inline double mass( const FourMomentum& p4, const ParticleBase& p) {
691 return mass(p4, p.mom());
695 inline double mass( const ParticleBase& p1, const ParticleBase& p2) {
696 return mass(p1.mom(), p2.mom());
700 inline double mass2( const ParticleBase& p, const P4& p4) {
701 return mass2(p.mom(), p4);
705 inline double mass2( const P4& p4, const ParticleBase& p) {
706 return mass2(p4, p.mom());
710 inline double mass2( const ParticleBase& p1, const ParticleBase& p2) {
711 return mass2(p1.mom(), p2.mom());
718 inline double mT( const ParticleBase& p, const P4& p4) {
719 return mT(p.mom(), p4);
726 inline double mT( const P4& p4, const ParticleBase& p) {
727 return mT(p4, p.mom());
734 inline double mT( const ParticleBase& p1, const ParticleBase& p2) {
735 return mT(p1.mom(), p2.mom());
739 inline double pT( const ParticleBase& p, const P4& p4) {
740 return pT(p.mom(), p4);
744 inline double pT( const P4& p4, const ParticleBase& p) {
745 return pT(p4, p.mom());
749 inline double pT( const ParticleBase& p1, const ParticleBase& p2) {
750 return pT(p1.mom(), p2.mom());
769 template < typename CONTAINER, typename = isCIterable<CONTAINER>>
770 inline int closestMassIndex(CONTAINER&& c, double mtarget, double mmin=-DBL_MAX, double mmax=DBL_MAX) {
771 return closestMatchIndex(std::forward<CONTAINER>(c), Kin::mass, mtarget, mmin, mmax);
778 template < typename CONTAINER1, typename CONTAINER2, typename = isCIterable<CONTAINER1, CONTAINER2>>
780 double mtarget, double mmin=-DBL_MAX, double mmax=DBL_MAX) {
782 std::forward<CONTAINER2>(c2), Kin::mass, mtarget, mmin, mmax);
789 template < typename CONTAINER, typename T, typename = isCIterable<CONTAINER>>
791 double mtarget, double mmin=-DBL_MAX, double mmax=DBL_MAX) {
793 return closestMatchIndex<CONTAINER, T, FN>(std::forward<CONTAINER>(c), x, Kin::mass, mtarget, mmin, mmax);
801 template < typename CONTAINER, typename T, typename = isCIterable<CONTAINER>>
803 double mtarget, double mmin=-DBL_MAX, double mmax=DBL_MAX) {
804 return closestMatchIndex(std::forward<T>(x), std::forward<CONTAINER>(c), Kin::mass, mtarget, mmin, mmax);
Specialized version of the FourVector with momentum/energy functionality. Definition Vector4.hh:316
FourMomentum & setPz(double pz) Set z-component of momentum . Definition Vector4.hh:370
FourMomentum & setPy(double py) Set y-component of momentum . Definition Vector4.hh:364
double rap() const Alias for rapidity. Definition Vector4.hh:611
FourMomentum & setPx(double px) Set x-component of momentum . Definition Vector4.hh:358
double absrap() const Absolute rapidity. Definition Vector4.hh:620
double pT() const Calculate the transverse momentum . Definition Vector4.hh:643
Vector3 vector3() const Get the spatial part of the 4-vector as a 3-vector. Definition Vector4.hh:177
double abseta() const Get the directly (alias). Definition Vector4.hh:174
double eta() const Synonym for pseudorapidity. Definition Vector4.hh:167
Base class for particle-like things like Particle and Jet. Definition ParticleBase.hh:13
const FourMomentum & mom() const Get equivalent single momentum four-vector (const) (alias). Definition ParticleBase.hh:39
ThreeMomentum p3() const Get the 3-momentum directly. Definition ParticleBase.hh:108
Three-dimensional specialisation of Vector. Definition Vector3.hh:40
double mod() const Calculate the modulus of a vector. . Definition VectorN.hh:95
pair< int, int > closestMatchIndices(CONTAINER1 &&c1, 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:788
int closestMatchIndex(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:759
bool any(const CONTAINER &c) Return true if x is true for any x in container c, otherwise false. Definition Utils.hh:330
bool all(const CONTAINER &c) Return true if x is true for all x in container c, otherwise false. Definition Utils.hh:355
Jets & idiscard(Jets &jets, const Cut &c) Filter a jet collection in-place to the subset that fails the supplied Cut.
double mass2(const FourMomentum &a, const FourMomentum &b) Calculate mass^2 of two 4-vectors. Definition Vector4.hh:1466
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
pair< int, int > closestMassIndices(CONTAINER1 &&c1, CONTAINER2 &&c2, double mtarget, double mmin=-DBL_MAX, double mmax=DBL_MAX) Return the indices from two vectors which best match fn(c1[i], c2[j]) to the target value. Definition ParticleBaseUtils.hh:779
int closestMassIndex(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:770
Definition MC_CENT_PPB_Projections.hh:10
double deltaR(double rap1, double phi1, double rap2, double phi2) Definition MathUtils.hh:698
double deltaPhi(double phi1, double phi2, bool sign=false) Calculate the difference between two angles in radians. Definition MathUtils.hh:668
double deltaEta(double eta1, double eta2, bool sign=false) Definition MathUtils.hh:676
double mT(double pT1, double pT2, double dphi) Definition MathUtils.hh:720
RapScheme Enum for rapidity variable to be used in calculating , applying rapidity cuts, etc. Definition MathConstants.hh:46
double deltaRap(double y1, double y2, bool sign=false) Definition MathUtils.hh:684
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
|