rivet is hosted by Hepforge, IPPP Durham
Rivet 4.0.2
Particle.hh
1// -*- C++ -*-
2#ifndef RIVET_Particle_HH
3#define RIVET_Particle_HH
4
5#include "Rivet/Particle.fhh"
6#include "Rivet/ParticleBase.hh"
7#include "Rivet/Config/RivetCommon.hh"
8#include "Rivet/Tools/Cuts.hh"
9#include "Rivet/Tools/Utils.hh"
10#include "Rivet/Tools/RivetFastJet.hh"
11#include "Rivet/Math/LorentzTrans.hh"
12// NOTE: Rivet/Tools/ParticleUtils.hh included at the end
13
14namespace Rivet {
15
16
21 class Particles : public std::vector<Particle> {
22 public:
23 using base = std::vector<Particle>; //< using-declarations don't like template syntax
24 using base::base; //< import base-class constructors
25 Particles();
26 Particles(const std::vector<Particle>& vps);
27 FourMomenta moms() const;
28 PseudoJets pseudojets() const;
29 operator FourMomenta () const { return moms(); }
30 operator PseudoJets () const { return pseudojets(); }
31 Particles& operator += (const Particle& p);
32 Particles& operator += (const Particles& ps);
33 };
34
35 Particles operator + (const Particles& a, const Particles& b);
36
38 typedef std::pair<Particle, Particle> ParticlePair;
39
40
42
43
45 class Particle : public ParticleBase {
46 public:
47
50
54 : ParticleBase(),
55 _original(nullptr), _id(PID::ANY), _isDirect(4, std::make_pair(false,false))
56 { }
57
59 Particle(PdgId pid, const FourMomentum& mom, const FourVector& pos=FourVector(), ConstGenParticlePtr gp=nullptr)
60 : ParticleBase(),
61 _original(gp), _id(pid),
62 _momentum(mom), _origin(pos),
63 _isDirect(4, std::make_pair(false,false))
64 { }
65
67 Particle(PdgId pid, const FourMomentum& mom, ConstGenParticlePtr gp, const FourVector& pos=FourVector())
68 : Particle(pid, mom, pos, gp)
69 { }
70
72 Particle(ConstGenParticlePtr gp)
73 : ParticleBase(),
74 _original(gp), _id(gp->pdg_id()),
75 _momentum(gp->momentum()),
76 _isDirect(4, std::make_pair(false,false))
77 {
78 ConstGenVertexPtr vprod = gp->production_vertex();
79 if (vprod != nullptr) {
80 setOrigin(vprod->position().t(), vprod->position().x(), vprod->position().y(), vprod->position().z());
81 }
82 }
83
85 Particle(const RivetHepMC::GenParticle& gp)
86 : Particle(HepMCUtils::getParticlePtr(gp))
87 { }
88
90
91
94
96 const FourMomentum& momentum() const {
97 return _momentum;
98 }
99
102 _momentum = momentum;
103 return *this;
104 }
105
107 Particle& setMomentum(double E, double px, double py, double pz) {
108 _momentum = FourMomentum(E, px, py, pz);
109 return *this;
110 }
111
114
116
117
120
122 const FourVector& origin() const {
123 return _origin;
124 }
127 _origin = position;
128 return *this;
129 }
131 Particle& setOrigin(double t, double x, double y, double z) {
132 _origin = FourMomentum(t, x, y, z);
133 return *this;
134 }
135
137
138
141
144 const FourVector& v0 = origin();
146 const double rho0 = origin().perp() / sin(this->phi() - origin().phi());
147 const double phi0 = M_PI/2 - this->phi();
148 const double x0 = rho0 * cos(phi0);
149 const double y0 = rho0 * sin(phi0);
150 const double z0 = origin().z() - v0.perp()/tan(this->theta());
151 return Vector3(x0, y0, z0);
152 }
153
155
156
159
161 virtual fastjet::PseudoJet pseudojet() const {
162 return fastjet::PseudoJet(mom().px(), mom().py(), mom().pz(), mom().E());
163 }
164
166 operator PseudoJet () const { return pseudojet(); }
167
168
170 Particle& setGenParticle(ConstGenParticlePtr gp) {
171 _original = gp;
172 return *this;
173 }
174
176 ConstGenParticlePtr genParticle() const {
177 return _original;
178 }
179
182 explicit operator ConstGenParticlePtr () const { return genParticle(); }
183
185
186
189
191 PdgId pid() const { return _id; }
193 PdgId abspid() const { return std::abs(_id); }
194
196
197
200
202 double charge() const { return PID::charge(pid()); }
203
205 double abscharge() const { return PID::abscharge(pid()); }
206
208 int charge3() const { return PID::charge3(pid()); }
209
211 int abscharge3() const { return PID::abscharge3(pid()); }
212
214 bool isCharged() const { return charge3() != 0; }
215
217
218
221
223 bool isHadron() const { return PID::isHadron(pid()); }
224
226 bool isMeson() const { return PID::isMeson(pid()); }
227
229 bool isBaryon() const { return PID::isBaryon(pid()); }
230
232 bool isLepton() const { return PID::isLepton(pid()); }
233
235 bool isChargedLepton() const { return PID::isChargedLepton(pid()); }
236
238 bool isNeutrino() const { return PID::isNeutrino(pid()); }
239
241 bool hasBottom() const { return PID::hasBottom(pid()); }
242
244 bool hasCharm() const { return PID::hasCharm(pid()); }
245
246 // /// Does this (hadron) contain an s quark?
247 // bool hasStrange() const { return PID::hasStrange(pid()); }
248
250 bool isVisible() const;
251
253 bool isParton() const { return PID::isParton(pid()); }
254
256
257
260
262 virtual void setConstituents(const Particles& cs, bool setmom=false);
263
265 virtual void addConstituent(const Particle& c, bool addmom=false);
266
268 virtual void addConstituents(const Particles& cs, bool addmom=false);
269
270
272 bool isComposite() const { return !constituents().empty(); }
273
274
279 const Particles& constituents() const { return _constituents; }
280
284 return sortBy(constituents(), sorter);
285 }
286
289 const Particles constituents(const Cut& c) const {
290 return select(constituents(), c);
291 }
292
295 const Particles constituents(const Cut& c, const ParticleSorter& sorter) const {
296 return sortBy(constituents(c), sorter);
297 }
298
301 const Particles constituents(const ParticleSelector& selector) const {
302 return select(constituents(), selector);
303 }
304
307 const Particles constituents(const ParticleSelector& selector, const ParticleSorter& sorter) const {
308 return sortBy(constituents(selector), sorter);
309 }
310
311
315
319 return sortBy(rawConstituents(), sorter);
320 }
321
324 const Particles rawConstituents(const Cut& c) const {
325 return select(rawConstituents(), c);
326 }
327
330 const Particles rawConstituents(const Cut& c, const ParticleSorter& sorter) const {
331 return sortBy(rawConstituents(c), sorter);
332 }
333
336 const Particles rawConstituents(const ParticleSelector& selector) const {
337 return select(rawConstituents(), selector);
338 }
339
342 const Particles rawConstituents(const ParticleSelector& selector, const ParticleSorter& sorter) const {
343 return sortBy(rawConstituents(selector), sorter);
344 }
345
347
348
351
357 Particles parents(const Cut& c=Cuts::OPEN) const;
358
365 return select(parents(), f);
366 }
367
373 bool hasParentWith(const ParticleSelector& f) const {
374 return !parents(f).empty();
375 }
381 bool hasParentWith(const Cut& c) const;
382
388 bool hasParentWithout(const ParticleSelector& f) const {
389 return hasParentWith([&](const Particle& p){ return !f(p); });
390 }
396 bool hasParentWithout(const Cut& c) const;
397
398
405 Particles ancestors(const Cut& c=Cuts::OPEN, bool only_physical=true) const;
406
414 return select(ancestors(Cuts::OPEN, only_physical), f);
415 }
416
422 bool hasAncestorWith(const ParticleSelector& f, bool only_physical=true) const {
423 return !ancestors(f, only_physical).empty();
424 }
430 bool hasAncestorWith(const Cut& c, bool only_physical=true) const;
431
437 bool hasAncestorWithout(const ParticleSelector& f, bool only_physical=true) const {
438 return hasAncestorWith([&](const Particle& p){ return !f(p); }, only_physical);
439 }
445 bool hasAncestorWithout(const Cut& c, bool only_physical=true) const;
446
452 bool fromBottom() const;
453
462 bool fromCharm() const;
463
464 // /// @brief Determine whether the particle is from a s-hadron decay
465 // ///
466 // /// @note If a hadron contains b or c quarks as well as strange it is
467 // /// considered a b or c hadron, but NOT a strange hadron.
468 // ///
469 // /// @note This question is valid in MC, but may not be perfectly answerable
470 // /// experimentally -- use this function with care when replicating
471 // /// experimental analyses!
472 // bool fromStrange() const;
473
479 bool fromHadron() const;
480
486 bool fromTau(bool prompt_taus_only=false) const;
487
493 bool fromPromptTau() const { return fromTau(true); }
494
500 bool fromHadronicTau(bool prompt_taus_only=false) const;
501
512 bool isDirect(bool allow_from_direct_tau=false, bool allow_from_direct_mu=false) const;
513
518
520
521
524
526 bool isStable() const;
527
529
530
532 Particles children(const Cut& c=Cuts::OPEN) const;
533
536 return select(children(), f);
537 }
538
544 bool hasChildWith(const ParticleSelector& f) const {
545 return !children(f).empty();
546 }
552 bool hasChildWith(const Cut& c) const;
553
559 bool hasChildWithout(const ParticleSelector& f) const {
560 return hasChildWith([&](const Particle& p){ return !f(p); });
561 }
567 bool hasChildWithout(const Cut& c) const;
568
569
571 Particles allDescendants(const Cut& c=Cuts::OPEN, bool remove_duplicates=true) const;
572
575 return select(allDescendants(Cuts::OPEN, remove_duplicates), f);
576 }
577
583 bool hasDescendantWith(const ParticleSelector& f, bool remove_duplicates=true) const {
584 return !allDescendants(f, remove_duplicates).empty();
585 }
591 bool hasDescendantWith(const Cut& c, bool remove_duplicates=true) const;
592
599 return hasDescendantWith([&](const Particle& p){ return !f(p); }, remove_duplicates);
600 }
606 bool hasDescendantWithout(const Cut& c, bool remove_duplicates=true) const;
607
608
613 Particles stableDescendants(const Cut& c=Cuts::OPEN) const;
614
617 return select(stableDescendants(), f);
618 }
619
626 return !stableDescendants(f).empty();
627 }
633 bool hasStableDescendantWith(const Cut& c) const;
634
641 return hasStableDescendantWith([&](const Particle& p){ return !f(p); });
642 }
648 bool hasStableDescendantWithout(const Cut& c) const;
649
650
654 double flightLength() const;
655
657
658
661
663 inline bool isFirstWith(const ParticleSelector& f) const {
664 if (!f(*this)) return false; //< This doesn't even meet f, let alone being the last to do so
665 if (any(parents(), f)) return false; //< If a direct parent has this property, this isn't the first
666 return true;
667 }
668
670 inline bool isFirstWithout(const ParticleSelector& f) const {
671 return isFirstWith([&](const Particle& p){ return !f(p); });
672 }
673
675 inline bool isLastWith(const ParticleSelector& f) const {
676 if (!f(*this)) return false; //< This doesn't even meet f, let alone being the last to do so
677 if (any(children(), f)) return false; //< If a child has this property, this isn't the last
678 return true;
679 }
680
682 inline bool isLastWithout(const ParticleSelector& f) const {
683 return isLastWith([&](const Particle& p){ return !f(p); });
684 }
685
687
688
691
695 bool isSame(const Particle& other) const {
696 if (pid() != other.pid()) return false;
697 if (!isZero((mom() - other.mom()).mod())) return false;
698 if (!isZero((origin() - other.origin()).mod())) return false;
699 return true;
700 }
701
703
704
705 protected:
706
708 ConstGenParticlePtr _original;
709
711 Particles _constituents;
712
714 PdgId _id;
715
717 FourMomentum _momentum;
718
720 FourVector _origin;
721
724 mutable std::vector<std::pair<bool,bool> > _isDirect;
725
726 };
727
728
731
733 std::ostream& operator << (std::ostream& os, const Particle& p);
734
736 std::ostream& operator << (std::ostream& os, const ParticlePair& pp);
737
739
740
741}
742
743
744#include "Rivet/Tools/ParticleUtils.hh"
745
746#endif
Specialized version of the FourVector with momentum/energy functionality.
Definition Vector4.hh:316
Specialisation of VectorN to a general (non-momentum) Lorentz 4-vector.
Definition Vector4.hh:30
double perp() const
Synonym for polarRadius.
Definition Vector4.hh:123
Object implementing Lorentz transform calculations and boosts.
Definition LorentzTrans.hh:21
Base class for particle-like things like Particle and Jet.
Definition ParticleBase.hh:13
double py() const
y component of momentum.
Definition ParticleBase.hh:122
const FourMomentum & mom() const
Get equivalent single momentum four-vector (const) (alias).
Definition ParticleBase.hh:39
double phi(const PhiMapping mapping=ZERO_2PI) const
Get the directly.
Definition ParticleBase.hh:105
double p() const
Get the 3-momentum magnitude directly.
Definition ParticleBase.hh:110
double theta() const
Synonym for polarAngle.
Definition ParticleBase.hh:136
double pz() const
z component of momentum.
Definition ParticleBase.hh:124
double px() const
x component of momentum.
Definition ParticleBase.hh:120
double eta() const
Get the directly (alias).
Definition ParticleBase.hh:87
double E() const
Get the energy directly.
Definition ParticleBase.hh:51
Particle representation, either from a HepMC::GenEvent or reconstructed.
Definition Particle.hh:45
Particle & setMomentum(double E, double px, double py, double pz)
Set the momentum via components.
Definition Particle.hh:107
bool isHadron() const
Is this a hadron?
Definition Particle.hh:223
Particles allDescendants(const ParticleSelector &f, bool remove_duplicates=true) const
Get a list of all the descendants from the current particle (with selector function)
Definition Particle.hh:574
const Particles constituents(const ParticleSelector &selector, const ParticleSorter &sorter) const
Direct constituents of this particle, filtered and sorted by functors.
Definition Particle.hh:307
const Particles rawConstituents(const ParticleSorter &sorter) const
Fundamental constituents of this particle, sorted by a functor.
Definition Particle.hh:318
bool hasStableDescendantWithout(const ParticleSelector &f) const
Definition Particle.hh:640
virtual fastjet::PseudoJet pseudojet() const
Converter to FastJet3 PseudoJet.
Definition Particle.hh:161
bool hasChildWith(const ParticleSelector &f) const
Definition Particle.hh:544
Vector3 closestApproach() const
Find the point of closest approach to the primary vertex.
Definition Particle.hh:143
bool isDirect(bool allow_from_direct_tau=false, bool allow_from_direct_mu=false) const
Shorthand definition of 'promptness' based on set definition flags.
bool hasDescendantWith(const Cut &c, bool remove_duplicates=true) const
virtual void addConstituent(const Particle &c, bool addmom=false)
Add a single direct constituent to this particle.
bool isPrompt(bool allow_from_prompt_tau=false, bool allow_from_prompt_mu=false) const
Alias for isDirect.
Definition Particle.hh:515
bool hasAncestorWith(const ParticleSelector &f, bool only_physical=true) const
Definition Particle.hh:422
Particles rawConstituents() const
Fundamental constituents of this particle.
bool hasAncestorWith(const Cut &c, bool only_physical=true) const
const Particles constituents(const Cut &c, const ParticleSorter &sorter) const
Direct constituents of this particle, sorted by a functor.
Definition Particle.hh:295
virtual void addConstituents(const Particles &cs, bool addmom=false)
Add direct constituents to this particle.
bool fromBottom() const
Determine whether the particle is from a b-hadron decay.
bool hasParentWithout(const Cut &c) const
bool isNeutrino() const
Is this a neutrino?
Definition Particle.hh:238
bool isSame(const Particle &other) const
Definition Particle.hh:695
const Particles rawConstituents(const ParticleSelector &selector, const ParticleSorter &sorter) const
Fundamental constituents of this particle, filtered and sorted by functors.
Definition Particle.hh:342
Particle & setOrigin(const FourVector &position)
Set the origin position.
Definition Particle.hh:126
Particle()
Definition Particle.hh:53
const Particles & constituents() const
Direct constituents of this particle, returned by reference.
Definition Particle.hh:279
bool hasChildWithout(const Cut &c) const
bool hasParentWith(const ParticleSelector &f) const
Definition Particle.hh:373
bool isLepton() const
Is this a lepton?
Definition Particle.hh:232
bool isChargedLepton() const
Is this a charged lepton?
Definition Particle.hh:235
bool hasParentWith(const Cut &c) const
bool hasBottom() const
Does this (hadron) contain a b quark?
Definition Particle.hh:241
bool hasChildWithout(const ParticleSelector &f) const
Definition Particle.hh:559
bool isCharged() const
Is this Particle charged?
Definition Particle.hh:214
bool isMeson() const
Is this a meson?
Definition Particle.hh:226
bool isLastWith(const ParticleSelector &f) const
Determine whether a particle is the last in a decay chain to meet the function requirement.
Definition Particle.hh:675
Particles ancestors(const Cut &c=Cuts::OPEN, bool only_physical=true) const
const Particles constituents(const ParticleSelector &selector) const
Direct constituents of this particle, filtered by a selection functor.
Definition Particle.hh:301
double abscharge() const
The absolute charge of this Particle.
Definition Particle.hh:205
bool isFirstWith(const ParticleSelector &f) const
Determine whether a particle is the first in a decay chain to meet the function requirement.
Definition Particle.hh:663
bool hasDescendantWith(const ParticleSelector &f, bool remove_duplicates=true) const
Definition Particle.hh:583
bool hasStableDescendantWith(const Cut &c) const
PdgId pid() const
This Particle's PDG ID code.
Definition Particle.hh:191
bool isLastWithout(const ParticleSelector &f) const
Determine whether a particle is the last in a decay chain not to meet the function requirement.
Definition Particle.hh:682
Particles ancestors(const ParticleSelector &f, bool only_physical=true) const
Definition Particle.hh:413
Particle(const RivetHepMC::GenParticle &gp)
Constructor from a HepMC GenParticle reference.
Definition Particle.hh:85
const Particles constituents(const Cut &c) const
Direct constituents of this particle, filtered by a Cut.
Definition Particle.hh:289
double flightLength() const
bool fromHadron() const
Determine whether the particle is from a hadron decay.
bool hasAncestorWithout(const Cut &c, bool only_physical=true) const
bool hasChildWith(const Cut &c) const
const Particles rawConstituents(const ParticleSelector &selector) const
Fundamental constituents of this particle, filtered by a selection functor.
Definition Particle.hh:336
double charge() const
The charge of this Particle.
Definition Particle.hh:202
Particles allDescendants(const Cut &c=Cuts::OPEN, bool remove_duplicates=true) const
Get a list of all the descendants from the current particle (with optional selection Cut)
bool fromTau(bool prompt_taus_only=false) const
Determine whether the particle is from a tau decay.
Particle & setGenParticle(ConstGenParticlePtr gp)
Set a const pointer to the original GenParticle.
Definition Particle.hh:170
bool hasParentWithout(const ParticleSelector &f) const
Definition Particle.hh:388
Particle & transformBy(const LorentzTransform &lt)
Apply an active Lorentz transform to this particle.
Particles parents(const Cut &c=Cuts::OPEN) const
int abscharge3() const
Three times the absolute charge of this Particle (i.e. integer multiple of smallest quark charge).
Definition Particle.hh:211
Particles stableDescendants(const ParticleSelector &f) const
Get a list of all the stable descendants from the current particle (with selector function)
Definition Particle.hh:616
bool isComposite() const
Determine if this Particle is a composite of other Rivet Particles.
Definition Particle.hh:272
Particle & setOrigin(double t, double x, double y, double z)
Set the origin position via components.
Definition Particle.hh:131
bool fromPromptTau() const
Determine whether the particle is from a prompt tau decay.
Definition Particle.hh:493
bool hasCharm() const
Does this (hadron) contain a c quark?
Definition Particle.hh:244
bool fromHadronicTau(bool prompt_taus_only=false) const
Determine whether the particle is from a tau which decayed hadronically.
bool isFirstWithout(const ParticleSelector &f) const
Determine whether a particle is the first in a decay chain not to meet the function requirement.
Definition Particle.hh:670
bool isStable() const
Whether this particle is stable according to the generator.
const FourMomentum & momentum() const
The momentum.
Definition Particle.hh:96
Particle(ConstGenParticlePtr gp)
Constructor from a HepMC GenParticle pointer.
Definition Particle.hh:72
const Particles rawConstituents(const Cut &c, const ParticleSorter &sorter) const
Fundamental constituents of this particle, sorted by a functor.
Definition Particle.hh:330
bool hasDescendantWithout(const Cut &c, bool remove_duplicates=true) const
Particles children(const ParticleSelector &f) const
Get a list of the direct descendants from the current particle (with selector function)
Definition Particle.hh:535
Particle(PdgId pid, const FourMomentum &mom, ConstGenParticlePtr gp, const FourVector &pos=FourVector())
Constructor from PID, momentum, and a GenParticle for relational links.
Definition Particle.hh:67
PdgId abspid() const
Absolute value of the PDG ID code.
Definition Particle.hh:193
const Particles rawConstituents(const Cut &c) const
Fundamental constituents of this particle, filtered by a Cut.
Definition Particle.hh:324
Particles parents(const ParticleSelector &f) const
Definition Particle.hh:364
Particles children(const Cut &c=Cuts::OPEN) const
Get a list of the direct descendants from the current particle (with optional selection Cut)
bool isVisible() const
Is this particle potentially visible in a detector?
bool hasDescendantWithout(const ParticleSelector &f, bool remove_duplicates=true) const
Definition Particle.hh:598
const Particles constituents(const ParticleSorter &sorter) const
Direct constituents of this particle, sorted by a functor.
Definition Particle.hh:283
bool fromCharm() const
Determine whether the particle is from a c-hadron decay.
bool hasStableDescendantWithout(const Cut &c) const
bool isBaryon() const
Is this a baryon?
Definition Particle.hh:229
Particles stableDescendants(const Cut &c=Cuts::OPEN) const
ConstGenParticlePtr genParticle() const
Get a const pointer to the original GenParticle.
Definition Particle.hh:176
bool hasStableDescendantWith(const ParticleSelector &f) const
Definition Particle.hh:625
bool isParton() const
Is this a parton? (Hopefully not very often... fiducial FTW)
Definition Particle.hh:253
int charge3() const
Three times the charge of this Particle (i.e. integer multiple of smallest quark charge).
Definition Particle.hh:208
bool hasAncestorWithout(const ParticleSelector &f, bool only_physical=true) const
Definition Particle.hh:437
Particle(PdgId pid, const FourMomentum &mom, const FourVector &pos=FourVector(), ConstGenParticlePtr gp=nullptr)
Constructor from PID and momentum.
Definition Particle.hh:59
Particle & setMomentum(const FourMomentum &momentum)
Set the momentum.
Definition Particle.hh:101
virtual void setConstituents(const Particles &cs, bool setmom=false)
Set direct constituents of this particle.
const FourVector & origin() const
The origin position (and time).
Definition Particle.hh:122
Specialised vector of Particle objects.
Definition Particle.hh:21
Three-dimensional specialisation of Vector.
Definition Vector3.hh:40
bool any(const CONTAINER &c)
Return true if x is true for any x in container c, otherwise false.
Definition Utils.hh:330
Jets select(const Jets &jets, const Cut &c)
Filter a jet collection in-place to the subset that passes the supplied Cut.
Definition JetUtils.hh:152
int abscharge3(int pid)
Return the absolute value of 3 times the EM charge.
Definition ParticleIdUtils.hh:928
double charge(int pid)
Return the EM charge (as floating point)
Definition ParticleIdUtils.hh:933
double abscharge(int pid)
Return the EM charge (as floating point)
Definition ParticleIdUtils.hh:938
int charge3(int pid)
Three times the EM charge (as integer)
Definition ParticleIdUtils.hh:858
bool isParton(int pid)
Determine if the PID is that of a parton (quark or gluon)
Definition ParticleIdUtils.hh:164
bool isNeutrino(int pid)
Determine if the PID is that of a neutrino.
Definition ParticleIdUtils.hh:196
bool isChargedLepton(int pid)
Determine if the PID is that of a charged lepton.
Definition ParticleIdUtils.hh:190
bool isLepton(int pid)
Definition ParticleIdUtils.hh:375
bool hasBottom(int pid)
Does this particle contain a bottom quark?
Definition ParticleIdUtils.hh:629
bool hasCharm(int pid)
Does this particle contain a charm quark?
Definition ParticleIdUtils.hh:625
bool isBaryon(int pid)
Check to see if this is a valid baryon.
Definition ParticleIdUtils.hh:297
bool isMeson(int pid)
Check to see if this is a valid meson.
Definition ParticleIdUtils.hh:274
bool isHadron(int pid)
Definition ParticleIdUtils.hh:357
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:1417
Definition MC_CENT_PPB_Projections.hh:10
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::pair< Particle, Particle > ParticlePair
Typedef for a pair of Particle objects.
Definition Particle.hh:38
std::vector< PseudoJet > PseudoJets
Definition RivetFastJet.hh:30
STL namespace.