Rivet API documentation

Rivet 4.1.3
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
123 const FourVector& origin() const {
124 return _origin;
125 }
126
127 Particle& setOrigin(const FourVector& position) {
128 _origin = position;
129 return *this;
130 }
131
132 Particle& setOrigin(double t, double x, double y, double z) {
133 _origin = FourMomentum(t, x, y, z);
134 return *this;
135 }
136
138
140
141
144
152 const Vector3 v0 = origin().vector3();
153 const Vector3 phat = p3().unit();
154 const double a = phat.perpVec().dot(v0.perpVec()) / phat.perp2();
155 const Vector3 rtn = v0 - a * phat;
156 return rtn;
157 }
158
159
160 double impactParam(bool signedip=false) const {
161 const Vector3 closestApproachVector = closestApproach();
162 double ip_trans = closestApproachVector.perp();
163
164 if (!signedip) return ip_trans;
165 const Vector3 productionVertexVector = origin().vector3();
166 const Vector3 particleMomentumVector = momentum().vector3();
167 const Vector3 displacementVector = productionVertexVector - closestApproachVector;
168 const double sign = displacementVector.dot(particleMomentumVector) > 0 ? 1.0 : -1.0;
169 return sign * ip_trans;
170 }
171
172
174
176
177
180
182 virtual fastjet::PseudoJet pseudojet() const {
183 return fastjet::PseudoJet(mom().px(), mom().py(), mom().pz(), mom().E());
184 }
185
187 operator PseudoJet () const { return pseudojet(); }
188
189
191 Particle& setGenParticle(ConstGenParticlePtr gp) {
192 _original = gp;
193 return *this;
194 }
195
197 ConstGenParticlePtr genParticle() const {
198 return _original;
199 }
200
203 explicit operator ConstGenParticlePtr () const { return genParticle(); }
204
206
207
210
212 PdgId pid() const { return _id; }
214 PdgId abspid() const { return std::abs(_id); }
215
217
218
221
223 double charge() const { return PID::charge(pid()); }
224
226 double abscharge() const { return PID::abscharge(pid()); }
227
229 int charge3() const { return PID::charge3(pid()); }
230
232 int abscharge3() const { return PID::abscharge3(pid()); }
233
235 bool isCharged() const { return charge3() != 0; }
236
238
239
242
244 bool isHadron() const { return PID::isHadron(pid()); }
245
247 bool isMeson() const { return PID::isMeson(pid()); }
248
250 bool isBaryon() const { return PID::isBaryon(pid()); }
251
253 bool isLepton() const { return PID::isLepton(pid()); }
254
256 bool isChargedLepton() const { return PID::isChargedLepton(pid()); }
257
259 bool isNeutrino() const { return PID::isNeutrino(pid()); }
260
262 bool hasBottom() const { return PID::hasBottom(pid()); }
263
265 bool hasCharm() const { return PID::hasCharm(pid()); }
266
267 // /// Does this (hadron) contain an s quark?
268 bool hasStrange() const { return PID::hasStrange(pid()); }
269
271 bool isVisible() const;
272
274 bool isParton() const { return PID::isParton(pid()); }
275
277
278
281
283 virtual void setConstituents(const Particles& cs, bool setmom=false);
284
286 virtual void addConstituent(const Particle& c, bool addmom=false);
287
289 virtual void addConstituents(const Particles& cs, bool addmom=false);
290
291
293 bool isComposite() const { return !constituents().empty(); }
294
295
300 const Particles& constituents() const { return _constituents; }
301
305 const Particles constituents(const ParticleSorter& sorter) const {
306 return sortBy(constituents(), sorter);
307 }
308
312 const Particles constituents(const Cut& c) const {
313 return select(constituents(), c);
314 }
315
319 const Particles constituents(const Cut& c, const ParticleSorter& sorter) const {
320 return sortBy(constituents(c), sorter);
321 }
322
326 const Particles constituents(const ParticleSelector& selector) const {
327 return select(constituents(), selector);
328 }
329
333 const Particles constituents(const ParticleSelector& selector, const ParticleSorter& sorter) const {
334 return sortBy(constituents(selector), sorter);
335 }
336
337
342
346 const Particles rawConstituents(const ParticleSorter& sorter) const {
347 return sortBy(rawConstituents(), sorter);
348 }
349
353 const Particles rawConstituents(const Cut& c) const {
354 return select(rawConstituents(), c);
355 }
356
360 const Particles rawConstituents(const Cut& c, const ParticleSorter& sorter) const {
361 return sortBy(rawConstituents(c), sorter);
362 }
363
367 const Particles rawConstituents(const ParticleSelector& selector) const {
368 return select(rawConstituents(), selector);
369 }
370
374 const Particles rawConstituents(const ParticleSelector& selector, const ParticleSorter& sorter) const {
375 return sortBy(rawConstituents(selector), sorter);
376 }
377
379
380
383
389 bool isPhysical(bool include_beams=false) const;
390
395 bool hasPhysicalAncestor(const Cut& c=Cuts::OPEN) const;
396
402 Particles parents(const Cut& c=Cuts::OPEN) const;
403
409 Particles parents(const ParticleSelector& f) const {
410 return select(parents(), f);
411 }
412
418 bool hasParentWith(const ParticleSelector& f) const {
419 return !parents(f).empty();
420 }
421
426 bool hasParentWith(const Cut& c) const;
427
433 bool hasParentWithout(const ParticleSelector& f) const {
434 return hasParentWith([&](const Particle& p){ return !f(p); });
435 }
436
441 bool hasParentWithout(const Cut& c) const;
442
443
450 Particles ancestors(const Cut& c=Cuts::OPEN, bool only_physical=true) const;
451
458 Particles ancestors(const ParticleSelector& f, bool only_physical=true) const {
459 return select(ancestors(Cuts::OPEN, only_physical), f);
460 }
461
467 bool hasAncestorWith(const ParticleSelector& f, bool only_physical=true) const {
468 return !ancestors(f, only_physical).empty();
469 }
470
475 bool hasAncestorWith(const Cut& c, bool only_physical=true) const;
476
482 bool hasAncestorWithout(const ParticleSelector& f, bool only_physical=true) const {
483 return hasAncestorWith([&](const Particle& p){ return !f(p); }, only_physical);
484 }
485
490 bool hasAncestorWithout(const Cut& c, bool only_physical=true) const;
491
497 bool fromBottom() const;
498
507 bool fromCharm() const;
508
509 // /// @brief Determine whether the particle is from a s-hadron decay
510 // ///
511 // /// @note If a hadron contains b or c quarks as well as strange it is
512 // /// considered a b or c hadron, but NOT a strange hadron.
513 // ///
514 // /// @note This question is valid in MC, but may not be perfectly answerable
515 // /// experimentally -- use this function with care when replicating
516 // /// experimental analyses!
517 // bool fromStrange() const;
518
524 bool fromHadron() const;
525
531 bool fromTau(bool prompt_taus_only=false) const;
532
538 bool fromPromptTau() const { return fromTau(true); }
539
545 bool fromHadronicTau(bool prompt_taus_only=false) const;
546
557 bool isDirect(bool allow_from_direct_tau=false, bool allow_from_direct_mu=false) const;
558
560 bool isPrompt(bool allow_from_prompt_tau=false, bool allow_from_prompt_mu=false) const {
561 return isDirect(allow_from_prompt_tau, allow_from_prompt_mu);
562 }
563
565
566
569
571 bool isStable() const;
572
574
575
577 Particles children(const Cut& c=Cuts::OPEN) const;
578
580 Particles children(const ParticleSelector& f) const {
581 return select(children(), f);
582 }
583
589 bool hasChildWith(const ParticleSelector& f) const {
590 return !children(f).empty();
591 }
592
597 bool hasChildWith(const Cut& c) const;
598
604 bool hasChildWithout(const ParticleSelector& f) const {
605 return hasChildWith([&](const Particle& p){ return !f(p); });
606 }
607
612 bool hasChildWithout(const Cut& c) const;
613
614
616 Particles allDescendants(const Cut& c=Cuts::OPEN, bool remove_duplicates=true) const;
617
619 Particles allDescendants(const ParticleSelector& f, bool remove_duplicates=true) const {
620 return select(allDescendants(Cuts::OPEN, remove_duplicates), f);
621 }
622
628 bool hasDescendantWith(const ParticleSelector& f, bool remove_duplicates=true) const {
629 return !allDescendants(f, remove_duplicates).empty();
630 }
631
636 bool hasDescendantWith(const Cut& c, bool remove_duplicates=true) const;
637
643 bool hasDescendantWithout(const ParticleSelector& f, bool remove_duplicates=true) const {
644 return hasDescendantWith([&](const Particle& p){ return !f(p); }, remove_duplicates);
645 }
646
651 bool hasDescendantWithout(const Cut& c, bool remove_duplicates=true) const;
652
653
658 Particles stableDescendants(const Cut& c=Cuts::OPEN) const;
659
661 Particles stableDescendants(const ParticleSelector& f) const {
662 return select(stableDescendants(), f);
663 }
664
670 bool hasStableDescendantWith(const ParticleSelector& f) const {
671 return !stableDescendants(f).empty();
672 }
673
678 bool hasStableDescendantWith(const Cut& c) const;
679
685 bool hasStableDescendantWithout(const ParticleSelector& f) const {
686 return hasStableDescendantWith([&](const Particle& p){ return !f(p); });
687 }
688
693 bool hasStableDescendantWithout(const Cut& c) const;
694
695
699 double flightLength() const;
700
702
703
706
708 inline bool isFirstWith(const ParticleSelector& f) const {
709 if (!f(*this)) return false; //< This doesn't even meet f, let alone being the last to do so
710 if (any(parents(), f)) return false; //< If a direct parent has this property, this isn't the first
711 return true;
712 }
713
715 inline bool isFirstWithout(const ParticleSelector& f) const {
716 return isFirstWith([&](const Particle& p){ return !f(p); });
717 }
718
720 inline bool isLastWith(const ParticleSelector& f) const {
721 if (!f(*this)) return false; //< This doesn't even meet f, let alone being the last to do so
722 if (any(children(), f)) return false; //< If a child has this property, this isn't the last
723 return true;
724 }
725
727 inline bool isLastWithout(const ParticleSelector& f) const {
728 return isLastWith([&](const Particle& p){ return !f(p); });
729 }
730
732
733
736
740 bool isSame(const Particle& other) const {
741 if (pid() != other.pid()) return false;
742 if (!isZero((mom() - other.mom()).mod())) return false;
743 if (!isZero((origin() - other.origin()).mod())) return false;
744 return true;
745 }
746
748
749
750 protected:
751
753 ConstGenParticlePtr _original;
754
756 Particles _constituents;
757
759 PdgId _id;
760
762 FourMomentum _momentum;
763
765 FourVector _origin;
766
769 mutable std::vector<std::pair<bool,bool> > _isDirect;
770
771 };
772
773
776
778 std::ostream& operator << (std::ostream& os, const Particle& p);
779
781 std::ostream& operator << (std::ostream& os, const ParticlePair& pp);
782
784
785
786}
787
788
789#include "Rivet/Tools/ParticleUtils.hh"
790
791#endif
Specialized version of the FourVector with momentum/energy functionality.
Definition Vector4.hh:327
Specialisation of VectorN to a general (non-momentum) Lorentz 4-vector.
Definition Vector4.hh:30
Vector3 vector3() const
Get the spatial part of the 4-vector as a 3-vector.
Definition Vector4.hh:184
Object implementing Lorentz transform calculations and boosts.
Definition LorentzTrans.hh:21
double py() const
y component of momentum.
Definition ParticleBase.hh:124
const FourMomentum & mom() const
Get the equivalent momentum four-vector (const) (alias).
Definition ParticleBase.hh:39
double p() const
Get the 3-momentum magnitude directly.
Definition ParticleBase.hh:112
ParticleBase()
Default constructor.
Definition ParticleBase.hh:17
ThreeMomentum p3() const
Get the 3-momentum directly.
Definition ParticleBase.hh:110
double pz() const
z component of momentum.
Definition ParticleBase.hh:126
double px() const
x component of momentum.
Definition ParticleBase.hh:122
double E() const
Get the energy directly.
Definition ParticleBase.hh:53
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:244
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:619
const Particles constituents(const ParticleSelector &selector, const ParticleSorter &sorter) const
Direct constituents of this particle, filtered and sorted by functors.
Definition Particle.hh:333
const Particles rawConstituents(const ParticleSorter &sorter) const
Fundamental constituents of this particle, sorted by a functor.
Definition Particle.hh:346
bool hasStableDescendantWithout(const ParticleSelector &f) const
Definition Particle.hh:685
virtual fastjet::PseudoJet pseudojet() const
Converter to FastJet3 PseudoJet.
Definition Particle.hh:182
bool hasChildWith(const ParticleSelector &f) const
Definition Particle.hh:589
Vector3 closestApproach() const
Find the point of closest approach (and hence the impact parameter) to the primary vertex.
Definition Particle.hh:151
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:560
bool hasAncestorWith(const ParticleSelector &f, bool only_physical=true) const
Definition Particle.hh:467
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:319
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:259
bool isSame(const Particle &other) const
Definition Particle.hh:740
const Particles rawConstituents(const ParticleSelector &selector, const ParticleSorter &sorter) const
Fundamental constituents of this particle, filtered and sorted by functors.
Definition Particle.hh:374
Particle & setOrigin(const FourVector &position)
Set the origin position.
Definition Particle.hh:127
Particle()
Definition Particle.hh:53
const Particles & constituents() const
Direct constituents of this particle, returned by reference.
Definition Particle.hh:300
bool hasChildWithout(const Cut &c) const
bool hasParentWith(const ParticleSelector &f) const
Definition Particle.hh:418
bool isLepton() const
Is this a lepton?
Definition Particle.hh:253
bool isChargedLepton() const
Is this a charged lepton?
Definition Particle.hh:256
bool hasParentWith(const Cut &c) const
bool hasBottom() const
Does this (hadron) contain a b quark?
Definition Particle.hh:262
bool hasChildWithout(const ParticleSelector &f) const
Definition Particle.hh:604
bool isCharged() const
Is this Particle charged?
Definition Particle.hh:235
bool isMeson() const
Is this a meson?
Definition Particle.hh:247
bool isPhysical(bool include_beams=false) const
Is this particle recorded as physically meaningful in the event record?
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:720
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:326
double abscharge() const
The absolute charge of this Particle.
Definition Particle.hh:226
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:708
bool hasDescendantWith(const ParticleSelector &f, bool remove_duplicates=true) const
Definition Particle.hh:628
bool hasStableDescendantWith(const Cut &c) const
PdgId pid() const
This Particle's PDG ID code.
Definition Particle.hh:212
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:727
Particles ancestors(const ParticleSelector &f, bool only_physical=true) const
Definition Particle.hh:458
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:312
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:367
double charge() const
The charge of this Particle.
Definition Particle.hh:223
bool hasPhysicalAncestor(const Cut &c=Cuts::OPEN) const
Does this (outgoing) particle have physically meaningful (outgoing) parents in the event record?
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:191
bool hasParentWithout(const ParticleSelector &f) const
Definition Particle.hh:433
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:232
Particles stableDescendants(const ParticleSelector &f) const
Get a list of all the stable descendants from the current particle (with selector function).
Definition Particle.hh:661
bool isComposite() const
Determine if this Particle is a composite of other Rivet Particles.
Definition Particle.hh:293
Particle & setOrigin(double t, double x, double y, double z)
Set the origin position via components.
Definition Particle.hh:132
bool fromPromptTau() const
Determine whether the particle is from a prompt tau decay.
Definition Particle.hh:538
bool hasCharm() const
Does this (hadron) contain a c quark?
Definition Particle.hh:265
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:715
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:360
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:580
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:214
const Particles rawConstituents(const Cut &c) const
Fundamental constituents of this particle, filtered by a Cut.
Definition Particle.hh:353
Particles parents(const ParticleSelector &f) const
Definition Particle.hh:409
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:643
const Particles constituents(const ParticleSorter &sorter) const
Direct constituents of this particle, sorted by a functor.
Definition Particle.hh:305
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:250
Particles stableDescendants(const Cut &c=Cuts::OPEN) const
ConstGenParticlePtr genParticle() const
Get a const pointer to the original GenParticle.
Definition Particle.hh:197
bool hasStableDescendantWith(const ParticleSelector &f) const
Definition Particle.hh:670
bool isParton() const
Is this a parton? (Hopefully not very often... fiducial FTW).
Definition Particle.hh:274
int charge3() const
Three times the charge of this Particle (i.e. integer multiple of smallest quark charge).
Definition Particle.hh:229
bool hasAncestorWithout(const ParticleSelector &f, bool only_physical=true) const
Definition Particle.hh:482
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
Definition Particle.hh:123
Specialised vector of Particle objects.
Definition Particle.hh:21
Three-dimensional specialisation of Vector.
Definition Vector3.hh:40
double perp() const
Synonym for polarRadius.
Definition Vector3.hh:161
double perp2() const
Synonym for polarRadius2.
Definition Vector3.hh:148
Vector3 unit() const
Synonym for unitVec.
Definition Vector3.hh:124
double dot(const Vector3 &v) const
Dot-product with another vector.
Definition Vector3.hh:96
Vector3 perpVec() const
Synonym for polarVec.
Definition Vector3.hh:135
bool any(const CONTAINER &c)
Return true if x is true for any x in container c, otherwise false.
Definition Utils.hh:358
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:157
int abscharge3(int pid)
Return the absolute value of 3 times the EM charge.
Definition ParticleIdUtils.hh:972
double charge(int pid)
Return the EM charge (as floating point).
Definition ParticleIdUtils.hh:977
double abscharge(int pid)
Return the EM charge (as floating point).
Definition ParticleIdUtils.hh:982
int charge3(int pid)
Three times the EM charge (as integer).
Definition ParticleIdUtils.hh:895
bool isParton(int pid)
Determine if the PID is that of a parton (quark or gluon).
Definition ParticleIdUtils.hh:180
bool isNeutrino(int pid)
Determine if the PID is that of a neutrino.
Definition ParticleIdUtils.hh:217
bool isChargedLepton(int pid)
Determine if the PID is that of a charged lepton.
Definition ParticleIdUtils.hh:206
bool isLepton(int pid)
Definition ParticleIdUtils.hh:412
bool hasBottom(int pid)
Does this particle contain a bottom quark?
Definition ParticleIdUtils.hh:666
bool hasStrange(int pid)
Does this particle contain a strange quark?
Definition ParticleIdUtils.hh:658
bool hasCharm(int pid)
Does this particle contain a charm quark?
Definition ParticleIdUtils.hh:662
bool isBaryon(int pid)
Check to see if this is a valid baryon.
Definition ParticleIdUtils.hh:318
bool isMeson(int pid)
Check to see if this is a valid meson.
Definition ParticleIdUtils.hh:295
bool isHadron(int pid)
Definition ParticleIdUtils.hh:394
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
double p(const ParticleBase &p)
Unbound function access to p.
Definition ParticleBaseUtils.hh:653
Definition RivetHepMC.hh:52
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
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.