rivet is hosted by Hepforge, IPPP Durham
Rivet 3.1.6
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
24 // typedef std::vector<Particle> Particles;
25 class Particles : public std::vector<Particle> {
26 public:
27 using base = std::vector<Particle>; //< using-declarations don't like template syntax
28 using base::base; //< import base-class constructors
29 Particles();
30 Particles(const std::vector<Particle>& vps);
31 FourMomenta moms() const;
32 PseudoJets pseudojets() const;
33 operator FourMomenta () const { return moms(); }
34 operator PseudoJets () const { return pseudojets(); }
35 Particles& operator += (const Particle& p);
36 Particles& operator += (const Particles& ps);
37 };
38
39 Particles operator + (const Particles& a, const Particles& b);
40
42 typedef std::pair<Particle, Particle> ParticlePair;
43
45
46
47
49
50
51
53 class Particle : public ParticleBase {
54 public:
55
58
62 : ParticleBase(),
63 _original(nullptr), _id(PID::ANY), _isDirect(4, std::make_pair(false,false))
64 { }
65
67 Particle(PdgId pid, const FourMomentum& mom, const FourVector& pos=FourVector(), ConstGenParticlePtr gp=nullptr)
68 : ParticleBase(),
69 _original(gp), _id(pid),
70 _momentum(mom), _origin(pos),
71 _isDirect(4, std::make_pair(false,false))
72 { }
73
75 Particle(PdgId pid, const FourMomentum& mom, ConstGenParticlePtr gp, const FourVector& pos=FourVector())
76 : Particle(pid, mom, pos, gp)
77 { }
78
80 Particle(ConstGenParticlePtr gp)
81 : ParticleBase(),
82 _original(gp), _id(gp->pdg_id()),
83 _momentum(gp->momentum()),
84 _isDirect(4, std::make_pair(false,false))
85 {
86 ConstGenVertexPtr vprod = gp->production_vertex();
87 if (vprod != nullptr) {
88 setOrigin(vprod->position().t(), vprod->position().x(), vprod->position().y(), vprod->position().z());
89 }
90 }
91
93 Particle(const RivetHepMC::GenParticle& gp)
94 : Particle(HepMCUtils::getParticlePtr(gp))
95 { }
96
98
101
103 const FourMomentum& momentum() const {
104 return _momentum;
105 }
106
109 _momentum = momentum;
110 return *this;
111 }
112
114 Particle& setMomentum(double E, double px, double py, double pz) {
115 _momentum = FourMomentum(E, px, py, pz);
116 return *this;
117 }
118
121
122 //@
123
124
127
129 const FourVector& origin() const {
130 return _origin;
131 }
133 Particle& setOrigin(const FourVector& position) {
134 _origin = position;
135 return *this;
136 }
138 Particle& setOrigin(double t, double x, double y, double z) {
139 _origin = FourMomentum(t, x, y, z);
140 return *this;
141 }
142
144
147
150 const FourVector& v0 = origin();
152 const double rho0 = origin().perp() / sin(this->phi() - origin().phi());
153 const double phi0 = M_PI/2 - this->phi();
154 const double x0 = rho0 * cos(phi0);
155 const double y0 = rho0 * sin(phi0);
156 const double z0 = origin().z() - v0.perp()/tan(this->theta());
157 return Vector3(x0, y0, z0);
158 }
159
161
162
165
167 virtual fastjet::PseudoJet pseudojet() const {
168 return fastjet::PseudoJet(mom().px(), mom().py(), mom().pz(), mom().E());
169 }
170
172 operator PseudoJet () const { return pseudojet(); }
173
174
176 Particle& setGenParticle(ConstGenParticlePtr gp) {
177 _original = gp;
178 return *this;
179 }
180
182 ConstGenParticlePtr genParticle() const {
183 return _original;
184 }
185
188 explicit operator ConstGenParticlePtr () const { return genParticle(); }
189
191
192
195
197 PdgId pid() const { return _id; }
199 PdgId abspid() const { return std::abs(_id); }
200
202
203
206
208 double charge() const { return PID::charge(pid()); }
209
211 double abscharge() const { return PID::abscharge(pid()); }
212
214 int charge3() const { return PID::charge3(pid()); }
215
217 int abscharge3() const { return PID::abscharge3(pid()); }
218
220 bool isCharged() const { return charge3() != 0; }
221
223
224
227
229 bool isHadron() const { return PID::isHadron(pid()); }
230
232 bool isMeson() const { return PID::isMeson(pid()); }
233
235 bool isBaryon() const { return PID::isBaryon(pid()); }
236
238 bool isLepton() const { return PID::isLepton(pid()); }
239
241 bool isChargedLepton() const { return PID::isChargedLepton(pid()); }
242
244 bool isNeutrino() const { return PID::isNeutrino(pid()); }
245
247 bool hasBottom() const { return PID::hasBottom(pid()); }
248
250 bool hasCharm() const { return PID::hasCharm(pid()); }
251
252 // /// Does this (hadron) contain an s quark?
253 // bool hasStrange() const { return PID::hasStrange(pid()); }
254
256 bool isVisible() const;
257
259 bool isParton() const { return PID::isParton(pid()); }
260
262
263
266
268 virtual void setConstituents(const Particles& cs, bool setmom=false);
269
271 virtual void addConstituent(const Particle& c, bool addmom=false);
272
274 virtual void addConstituents(const Particles& cs, bool addmom=false);
275
276
278 bool isComposite() const { return !constituents().empty(); }
279
280
285 const Particles& constituents() const { return _constituents; }
286
289 const Particles constituents(const ParticleSorter& sorter) const {
290 return sortBy(constituents(), sorter);
291 }
292
295 const Particles constituents(const Cut& c) const {
296 return filter_select(constituents(), c);
297 }
298
301 const Particles constituents(const Cut& c, const ParticleSorter& sorter) const {
302 return sortBy(constituents(c), sorter);
303 }
304
307 const Particles constituents(const ParticleSelector& selector) const {
308 return filter_select(constituents(), selector);
309 }
310
313 const Particles constituents(const ParticleSelector& selector, const ParticleSorter& sorter) const {
314 return sortBy(constituents(selector), sorter);
315 }
316
317
321
324 const Particles rawConstituents(const ParticleSorter& sorter) const {
325 return sortBy(rawConstituents(), sorter);
326 }
327
330 const Particles rawConstituents(const Cut& c) const {
331 return filter_select(rawConstituents(), c);
332 }
333
336 const Particles rawConstituents(const Cut& c, const ParticleSorter& sorter) const {
337 return sortBy(rawConstituents(c), sorter);
338 }
339
342 const Particles rawConstituents(const ParticleSelector& selector) const {
343 return filter_select(rawConstituents(), selector);
344 }
345
348 const Particles rawConstituents(const ParticleSelector& selector, const ParticleSorter& sorter) const {
349 return sortBy(rawConstituents(selector), sorter);
350 }
351
353
354
357
363 Particles parents(const Cut& c=Cuts::OPEN) const;
364
370 Particles parents(const ParticleSelector& f) const {
371 return filter_select(parents(), f);
372 }
373
379 bool hasParentWith(const ParticleSelector& f) const {
380 return !parents(f).empty();
381 }
387 bool hasParentWith(const Cut& c) const;
388
394 bool hasParentWithout(const ParticleSelector& f) const {
395 return hasParentWith([&](const Particle& p){ return !f(p); });
396 }
402 bool hasParentWithout(const Cut& c) const;
403
411 bool hasParent(PdgId pid) const;
412
413
414
421 Particles ancestors(const Cut& c=Cuts::OPEN, bool only_physical=true) const;
422
429 Particles ancestors(const ParticleSelector& f, bool only_physical=true) const {
430 return filter_select(ancestors(Cuts::OPEN, only_physical), f);
431 }
432
438 bool hasAncestorWith(const ParticleSelector& f, bool only_physical=true) const {
439 return !ancestors(f, only_physical).empty();
440 }
446 bool hasAncestorWith(const Cut& c, bool only_physical=true) const;
447
453 bool hasAncestorWithout(const ParticleSelector& f, bool only_physical=true) const {
454 return hasAncestorWith([&](const Particle& p){ return !f(p); }, only_physical);
455 }
461 bool hasAncestorWithout(const Cut& c, bool only_physical=true) const;
462
470 bool hasAncestor(PdgId pid, bool only_physical=true) const;
471
472
478 bool fromBottom() const;
479
488 bool fromCharm() const;
489
490 // /// @brief Determine whether the particle is from a s-hadron decay
491 // ///
492 // /// @note If a hadron contains b or c quarks as well as strange it is
493 // /// considered a b or c hadron, but NOT a strange hadron.
494 // ///
495 // /// @note This question is valid in MC, but may not be perfectly answerable
496 // /// experimentally -- use this function with care when replicating
497 // /// experimental analyses!
498 // bool fromStrange() const;
499
505 bool fromHadron() const;
506
512 bool fromTau(bool prompt_taus_only=false) const;
513
519 bool fromPromptTau() const { return fromTau(true); }
520
526 bool fromHadronicTau(bool prompt_taus_only=false) const;
527
537 DEPRECATED("Too vague: use fromHadron() || fromPromptTau(), or isDirect()")
538 bool fromDecay() const { return fromHadron() || fromPromptTau(); }
539
550 bool isDirect(bool allow_from_direct_tau=false, bool allow_from_direct_mu=false) const;
551
553 bool isPrompt(bool allow_from_prompt_tau=false, bool allow_from_prompt_mu=false) const {
554 return isDirect(allow_from_prompt_tau, allow_from_prompt_mu);
555 }
556
558
559
562
564 bool isStable() const;
565
567
568
570 Particles children(const Cut& c=Cuts::OPEN) const;
571
573 Particles children(const ParticleSelector& f) const {
574 return filter_select(children(), f);
575 }
576
582 bool hasChildWith(const ParticleSelector& f) const {
583 return !children(f).empty();
584 }
590 bool hasChildWith(const Cut& c) const;
591
597 bool hasChildWithout(const ParticleSelector& f) const {
598 return hasChildWith([&](const Particle& p){ return !f(p); });
599 }
605 bool hasChildWithout(const Cut& c) const;
606
607
609 Particles allDescendants(const Cut& c=Cuts::OPEN, bool remove_duplicates=true) const;
610
612 Particles allDescendants(const ParticleSelector& f, bool remove_duplicates=true) const {
613 return filter_select(allDescendants(Cuts::OPEN, remove_duplicates), f);
614 }
615
621 bool hasDescendantWith(const ParticleSelector& f, bool remove_duplicates=true) const {
622 return !allDescendants(f, remove_duplicates).empty();
623 }
629 bool hasDescendantWith(const Cut& c, bool remove_duplicates=true) const;
630
636 bool hasDescendantWithout(const ParticleSelector& f, bool remove_duplicates=true) const {
637 return hasDescendantWith([&](const Particle& p){ return !f(p); }, remove_duplicates);
638 }
644 bool hasDescendantWithout(const Cut& c, bool remove_duplicates=true) const;
645
646
651 Particles stableDescendants(const Cut& c=Cuts::OPEN) const;
652
654 Particles stableDescendants(const ParticleSelector& f) const {
655 return filter_select(stableDescendants(), f);
656 }
657
663 bool hasStableDescendantWith(const ParticleSelector& f) const {
664 return !stableDescendants(f).empty();
665 }
671 bool hasStableDescendantWith(const Cut& c) const;
672
678 bool hasStableDescendantWithout(const ParticleSelector& f) const {
679 return hasStableDescendantWith([&](const Particle& p){ return !f(p); });
680 }
686 bool hasStableDescendantWithout(const Cut& c) const;
687
688
692 double flightLength() const;
693
695
696
699
701 inline bool isFirstWith(const ParticleSelector& f) const {
702 if (!f(*this)) return false; //< This doesn't even meet f, let alone being the last to do so
703 if (any(parents(), f)) return false; //< If a direct parent has this property, this isn't the first
704 return true;
705 }
706
708 inline bool isFirstWithout(const ParticleSelector& f) const {
709 return isFirstWith([&](const Particle& p){ return !f(p); });
710 }
711
713 inline bool isLastWith(const ParticleSelector& f) const {
714 if (!f(*this)) return false; //< This doesn't even meet f, let alone being the last to do so
715 if (any(children(), f)) return false; //< If a child has this property, this isn't the last
716 return true;
717 }
718
720 inline bool isLastWithout(const ParticleSelector& f) const {
721 return isLastWith([&](const Particle& p){ return !f(p); });
722 }
723
725
726
729
733 bool isSame(const Particle& other) const {
734 if (pid() != other.pid()) return false;
735 if (!isZero((mom() - other.mom()).mod())) return false;
736 if (!isZero((origin() - other.origin()).mod())) return false;
737 return true;
738 }
739
741
742
743 protected:
744
746 ConstGenParticlePtr _original;
747
749 Particles _constituents;
750
752 PdgId _id;
753
755 FourMomentum _momentum;
756
758 FourVector _origin;
759
762 mutable std::vector<std::pair<bool,bool> > _isDirect;
763
764 };
765
766
769
771 std::ostream& operator << (std::ostream& os, const Particle& p);
772
774 std::ostream& operator << (std::ostream& os, const ParticlePair& pp);
775
777
778
779}
780
781
782#include "Rivet/Tools/ParticleUtils.hh"
783
784#endif
Specialized version of the FourVector with momentum/energy functionality.
Definition: Vector4.hh:306
Specialisation of VectorN to a general (non-momentum) Lorentz 4-vector.
Definition: Vector4.hh:27
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 E() const
Get the energy directly.
Definition: ParticleBase.hh:51
Particle representation, either from a HepMC::GenEvent or reconstructed.
Definition: Particle.hh:53
Particle & setMomentum(double E, double px, double py, double pz)
Set the momentum via components.
Definition: Particle.hh:114
bool isHadron() const
Is this a hadron?
Definition: Particle.hh:229
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:612
const Particles constituents(const ParticleSelector &selector, const ParticleSorter &sorter) const
Direct constituents of this particle, filtered and sorted by functors.
Definition: Particle.hh:313
const Particles rawConstituents(const ParticleSorter &sorter) const
Fundamental constituents of this particle, sorted by a functor.
Definition: Particle.hh:324
bool hasStableDescendantWithout(const ParticleSelector &f) const
Definition: Particle.hh:678
virtual fastjet::PseudoJet pseudojet() const
Converter to FastJet3 PseudoJet.
Definition: Particle.hh:167
bool hasChildWith(const ParticleSelector &f) const
Definition: Particle.hh:582
Vector3 closestApproach() const
Find the point of closest approach to the primary vertex.
Definition: Particle.hh:149
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:553
bool hasAncestorWith(const ParticleSelector &f, bool only_physical=true) const
Definition: Particle.hh:438
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:301
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:244
bool isSame(const Particle &other) const
Definition: Particle.hh:733
const Particles rawConstituents(const ParticleSelector &selector, const ParticleSorter &sorter) const
Fundamental constituents of this particle, filtered and sorted by functors.
Definition: Particle.hh:348
Particle & setOrigin(const FourVector &position)
Set the origin position.
Definition: Particle.hh:133
Particle()
Definition: Particle.hh:61
const Particles & constituents() const
Direct constituents of this particle, returned by reference.
Definition: Particle.hh:285
bool hasChildWithout(const Cut &c) const
bool hasParentWith(const ParticleSelector &f) const
Definition: Particle.hh:379
bool isLepton() const
Is this a lepton?
Definition: Particle.hh:238
bool isChargedLepton() const
Is this a charged lepton?
Definition: Particle.hh:241
bool hasParentWith(const Cut &c) const
bool hasBottom() const
Does this (hadron) contain a b quark?
Definition: Particle.hh:247
DEPRECATED("Too vague: use fromHadron() || fromPromptTau(), or isDirect()") bool fromDecay() const
Determine whether the particle is from a hadron or tau decay.
Definition: Particle.hh:537
bool hasChildWithout(const ParticleSelector &f) const
Definition: Particle.hh:597
bool isCharged() const
Is this Particle charged?
Definition: Particle.hh:220
bool isMeson() const
Is this a meson?
Definition: Particle.hh:232
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:713
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:307
double abscharge() const
The absolute charge of this Particle.
Definition: Particle.hh:211
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:701
bool hasDescendantWith(const ParticleSelector &f, bool remove_duplicates=true) const
Definition: Particle.hh:621
bool hasStableDescendantWith(const Cut &c) const
bool hasParent(PdgId pid) const
PdgId pid() const
This Particle's PDG ID code.
Definition: Particle.hh:197
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:720
Particles ancestors(const ParticleSelector &f, bool only_physical=true) const
Definition: Particle.hh:429
Particle(const RivetHepMC::GenParticle &gp)
Constructor from a HepMC GenParticle reference.
Definition: Particle.hh:93
const Particles constituents(const Cut &c) const
Direct constituents of this particle, filtered by a Cut.
Definition: Particle.hh:295
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:342
double charge() const
The charge of this Particle.
Definition: Particle.hh:208
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:176
bool hasParentWithout(const ParticleSelector &f) const
Definition: Particle.hh:394
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:217
Particles stableDescendants(const ParticleSelector &f) const
Get a list of all the stable descendants from the current particle (with selector function)
Definition: Particle.hh:654
bool isComposite() const
Determine if this Particle is a composite of other Rivet Particles.
Definition: Particle.hh:278
Particle & setOrigin(double t, double x, double y, double z)
Set the origin position via components.
Definition: Particle.hh:138
bool fromPromptTau() const
Determine whether the particle is from a prompt tau decay.
Definition: Particle.hh:519
bool hasCharm() const
Does this (hadron) contain a c quark?
Definition: Particle.hh:250
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:708
bool isStable() const
Whether this particle is stable according to the generator.
const FourMomentum & momentum() const
The momentum.
Definition: Particle.hh:103
Particle(ConstGenParticlePtr gp)
Constructor from a HepMC GenParticle pointer.
Definition: Particle.hh:80
const Particles rawConstituents(const Cut &c, const ParticleSorter &sorter) const
Fundamental constituents of this particle, sorted by a functor.
Definition: Particle.hh:336
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:573
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:75
PdgId abspid() const
Absolute value of the PDG ID code.
Definition: Particle.hh:199
const Particles rawConstituents(const Cut &c) const
Fundamental constituents of this particle, filtered by a Cut.
Definition: Particle.hh:330
Particles parents(const ParticleSelector &f) const
Definition: Particle.hh:370
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:636
const Particles constituents(const ParticleSorter &sorter) const
Direct constituents of this particle, sorted by a functor.
Definition: Particle.hh:289
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:235
Particles stableDescendants(const Cut &c=Cuts::OPEN) const
ConstGenParticlePtr genParticle() const
Get a const pointer to the original GenParticle.
Definition: Particle.hh:182
bool hasStableDescendantWith(const ParticleSelector &f) const
Definition: Particle.hh:663
bool isParton() const
Is this a parton? (Hopefully not very often... fiducial FTW)
Definition: Particle.hh:259
bool hasAncestor(PdgId pid, bool only_physical=true) const
int charge3() const
Three times the charge of this Particle (i.e. integer multiple of smallest quark charge).
Definition: Particle.hh:214
bool hasAncestorWithout(const ParticleSelector &f, bool only_physical=true) const
Definition: Particle.hh:453
Particle(PdgId pid, const FourMomentum &mom, const FourVector &pos=FourVector(), ConstGenParticlePtr gp=nullptr)
Constructor from PID and momentum.
Definition: Particle.hh:67
Particle & setMomentum(const FourMomentum &momentum)
Set the momentum.
Definition: Particle.hh:108
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:129
Specialised vector of Particle objects.
Definition: Particle.hh:25
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:334
Jets filter_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:865
double charge(int pid)
Return the EM charge (as floating point)
Definition: ParticleIdUtils.hh:868
double abscharge(int pid)
Return the EM charge (as floating point)
Definition: ParticleIdUtils.hh:871
int charge3(int pid)
Three times the EM charge (as integer)
Definition: ParticleIdUtils.hh:800
bool isParton(int pid)
Determine if the PID is that of a parton (quark or gluon)
Definition: ParticleIdUtils.hh:132
bool isNeutrino(int pid)
Determine if the PID is that of a neutrino.
Definition: ParticleIdUtils.hh:167
bool isChargedLepton(int pid)
Determine if the PID is that of a charged lepton.
Definition: ParticleIdUtils.hh:158
bool isLepton(int pid)
Definition: ParticleIdUtils.hh:340
bool hasBottom(int pid)
Does this particle contain a bottom quark?
Definition: ParticleIdUtils.hh:581
bool hasCharm(int pid)
Does this particle contain a charm quark?
Definition: ParticleIdUtils.hh:579
bool isBaryon(int pid)
Check to see if this is a valid baryon.
Definition: ParticleIdUtils.hh:260
bool isMeson(int pid)
Check to see if this is a valid meson.
Definition: ParticleIdUtils.hh:237
bool isHadron(int pid)
Definition: ParticleIdUtils.hh:322
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:1441
double perp() const
Synonym for polarRadius.
Definition: Vector4.hh:113
double p(const ParticleBase &p)
Unbound function access to p.
Definition: ParticleBaseUtils.hh:653
Definition: MC_Cent_pPb.hh:10
std::ostream & operator<<(std::ostream &os, const AnalysisInfo &ai)
Stream an AnalysisInfo as a text description.
Definition: AnalysisInfo.hh:362
std::pair< Particle, Particle > ParticlePair
Typedef for a pair of Particle objects.
Definition: Particle.hh:42
std::vector< PseudoJet > PseudoJets
Definition: RivetFastJet.hh:30
std::enable_if< std::is_floating_point< NUM >::value, bool >::type isZero(NUM val, double tolerance=1e-8)
Compare a number to zero.
Definition: MathUtils.hh:24
STL namespace.