rivet is hosted by Hepforge, IPPP Durham
Rivet  2.7.0
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/Math/LorentzTrans.hh"
11 // NOTE: Rivet/Tools/ParticleUtils.hh included at the end
12 #include "fastjet/PseudoJet.hh"
13 
14 namespace Rivet {
15 
16 
18  class Particle : public ParticleBase {
19  public:
20 
22 
23 
27  : ParticleBase(),
28  _original(nullptr), _id(PID::ANY), _isDirect{false,false}
29  { }
30 
32  Particle(PdgId pid, const FourMomentum& mom, const FourVector& pos=FourVector())
33  : ParticleBase(),
34  _original(nullptr), _id(pid),
35  _momentum(mom), _origin(pos),
36  _isDirect{false,false}
37  { }
38 
40  Particle(const GenParticle* gp)
41  : ParticleBase(),
42  _original(gp), _id(gp->pdg_id()),
43  _momentum(gp->momentum()),
44  _isDirect{false,false}
45  {
46  const GenVertex* vprod = gp->production_vertex();
47  if (vprod != nullptr) {
48  setOrigin(vprod->position().t(), vprod->position().x(), vprod->position().y(), vprod->position().z());
49  }
50  }
51 
53  Particle(const GenParticle& gp)
54  : Particle(&gp)
55  { }
56 
58 
59 
61 
62 
64  const FourMomentum& momentum() const {
65  return _momentum;
66  }
67 
70  _momentum = momentum;
71  return *this;
72  }
73 
75  Particle& setMomentum(double E, double px, double py, double pz) {
76  _momentum = FourMomentum(E, px, py, pz);
77  return *this;
78  }
79 
82 
83  //@
84 
85 
87 
88 
90  const FourVector& origin() const {
91  return _origin;
92  }
94  Particle& setOrigin(const FourVector& position) {
95  _origin = position;
96  return *this;
97  }
99  Particle& setOrigin(double t, double x, double y, double z) {
100  _origin = FourMomentum(t, x, y, z);
101  return *this;
102  }
103 
105 
106 
108 
109 
111  virtual fastjet::PseudoJet pseudojet() const {
112  return fastjet::PseudoJet(mom().px(), mom().py(), mom().pz(), mom().E());
113  }
114 
116  operator PseudoJet () const { return pseudojet(); }
117 
118 
120  const GenParticle* genParticle() const {
121  return _original;
122  }
123 
126  operator const GenParticle* () const { return genParticle(); }
127 
129 
130 
132 
133 
135  PdgId pid() const { return _id; }
137  PdgId abspid() const { return std::abs(_id); }
140  PdgId pdgId() const { return _id; }
141 
143 
144 
146 
147 
149  double charge() const { return PID::charge(pid()); }
150 
152  double abscharge() const { return PID::abscharge(pid()); }
153 
155  int charge3() const { return PID::charge3(pid()); }
156 
159  int threeCharge() const { return PID::threeCharge(pid()); }
160 
162  int abscharge3() const { return PID::abscharge3(pid()); }
163 
165  bool isCharged() const { return charge3() != 0; }
166 
168 
169 
171 
172 
174  bool isHadron() const { return PID::isHadron(pid()); }
175 
177  bool isMeson() const { return PID::isMeson(pid()); }
178 
180  bool isBaryon() const { return PID::isBaryon(pid()); }
181 
183  bool isLepton() const { return PID::isLepton(pid()); }
184 
186  bool isChargedLepton() const { return PID::isChargedLepton(pid()); }
187 
189  bool isNeutrino() const { return PID::isNeutrino(pid()); }
190 
192  bool hasBottom() const { return PID::hasBottom(pid()); }
193 
195  bool hasCharm() const { return PID::hasCharm(pid()); }
196 
197  // /// Does this (hadron) contain an s quark?
198  // bool hasStrange() const { return PID::hasStrange(pid()); }
199 
201  bool isVisible() const;
202 
204  bool isParton() const { return PID::isParton(pid()); }
205 
207 
208 
210 
211 
213  virtual void setConstituents(const Particles& cs, bool setmom=false);
214 
216  virtual void addConstituent(const Particle& c, bool addmom=false);
217 
219  virtual void addConstituents(const Particles& cs, bool addmom=false);
220 
221 
223  bool isComposite() const { return !constituents().empty(); }
224 
225 
230  const Particles& constituents() const { return _constituents; }
231 
234  const Particles constituents(const ParticleSorter& sorter) const {
235  return sortBy(constituents(), sorter);
236  }
237 
240  const Particles constituents(const Cut& c) const {
241  return filter_select(constituents(), c);
242  }
243 
246  const Particles constituents(const Cut& c, const ParticleSorter& sorter) const {
247  return sortBy(constituents(c), sorter);
248  }
249 
252  const Particles constituents(const ParticleSelector& selector) const {
253  return filter_select(constituents(), selector);
254  }
255 
258  const Particles constituents(const ParticleSelector& selector, const ParticleSorter& sorter) const {
259  return sortBy(constituents(selector), sorter);
260  }
261 
262 
265  Particles rawConstituents() const;
266 
269  const Particles rawConstituents(const ParticleSorter& sorter) const {
270  return sortBy(rawConstituents(), sorter);
271  }
272 
275  const Particles rawConstituents(const Cut& c) const {
276  return filter_select(rawConstituents(), c);
277  }
278 
281  const Particles rawConstituents(const Cut& c, const ParticleSorter& sorter) const {
282  return sortBy(rawConstituents(c), sorter);
283  }
284 
287  const Particles rawConstituents(const ParticleSelector& selector) const {
288  return filter_select(rawConstituents(), selector);
289  }
290 
293  const Particles rawConstituents(const ParticleSelector& selector, const ParticleSorter& sorter) const {
294  return sortBy(rawConstituents(selector), sorter);
295  }
296 
298 
299 
301 
302 
308  Particles parents(const Cut& c=Cuts::OPEN) const;
309 
315  Particles parents(const ParticleSelector& f) const {
316  return filter_select(parents(), f);
317  }
318 
324  bool hasParentWith(const ParticleSelector& f) const {
325  return !parents(f).empty();
326  }
332  bool hasParentWith(const Cut& c) const;
333 
339  bool hasParentWithout(const ParticleSelector& f) const {
340  return hasParentWith([&](const Particle& p){ return !f(p); });
341  }
347  bool hasParentWithout(const Cut& c) const;
348 
356  //DEPRECATED("Prefer e.g. hasParentWith(Cut::pid == 123)");
357  bool hasParent(PdgId pid) const;
358 
359 
360 
367  Particles ancestors(const Cut& c=Cuts::OPEN, bool only_physical=true) const;
368 
375  Particles ancestors(const ParticleSelector& f, bool only_physical=true) const {
376  return filter_select(ancestors(Cuts::OPEN, only_physical), f);
377  }
378 
384  bool hasAncestorWith(const ParticleSelector& f, bool only_physical=true) const {
385  return !ancestors(f, only_physical).empty();
386  }
392  bool hasAncestorWith(const Cut& c, bool only_physical=true) const;
393 
399  bool hasAncestorWithout(const ParticleSelector& f, bool only_physical=true) const {
400  return hasAncestorWith([&](const Particle& p){ return !f(p); }, only_physical);
401  }
407  bool hasAncestorWithout(const Cut& c, bool only_physical=true) const;
408 
416  //DEPRECATED("Prefer e.g. hasAncestorWith(Cut::pid == 123)");
417  bool hasAncestor(PdgId pid, bool only_physical=true) const;
418 
419 
425  bool fromBottom() const;
426 
435  bool fromCharm() const;
436 
437  // /// @brief Determine whether the particle is from a s-hadron decay
438  // ///
439  // /// @note If a hadron contains b or c quarks as well as strange it is
440  // /// considered a b or c hadron, but NOT a strange hadron.
441  // ///
442  // /// @note This question is valid in MC, but may not be perfectly answerable
443  // /// experimentally -- use this function with care when replicating
444  // /// experimental analyses!
445  // bool fromStrange() const;
446 
452  bool fromHadron() const;
453 
459  bool fromTau(bool prompt_taus_only=false) const;
460 
466  bool fromPromptTau() const { return fromTau(true); }
467 
473  bool fromHadronicTau(bool prompt_taus_only=false) const;
474 
485  bool fromDecay() const { return fromHadron() || fromPromptTau(); }
486 
497  bool isDirect(bool allow_from_direct_tau=false, bool allow_from_direct_mu=false) const;
498 
500  bool isPrompt(bool allow_from_prompt_tau=false, bool allow_from_prompt_mu=false) const {
501  return isDirect(allow_from_prompt_tau, allow_from_prompt_mu);
502  }
503 
505 
506 
508 
509 
511  bool isStable() const;
512 
514 
515 
517  Particles children(const Cut& c=Cuts::OPEN) const;
518 
520  Particles children(const ParticleSelector& f) const {
521  return filter_select(children(), f);
522  }
523 
529  bool hasChildWith(const ParticleSelector& f) const {
530  return !children(f).empty();
531  }
537  bool hasChildWith(const Cut& c) const;
538 
544  bool hasChildWithout(const ParticleSelector& f) const {
545  return hasChildWith([&](const Particle& p){ return !f(p); });
546  }
552  bool hasChildWithout(const Cut& c) const;
553 
554 
556  Particles allDescendants(const Cut& c=Cuts::OPEN, bool remove_duplicates=true) const;
557 
559  Particles allDescendants(const ParticleSelector& f, bool remove_duplicates=true) const {
560  return filter_select(allDescendants(Cuts::OPEN, remove_duplicates), f);
561  }
562 
568  bool hasDescendantWith(const ParticleSelector& f, bool remove_duplicates=true) const {
569  return !allDescendants(f, remove_duplicates).empty();
570  }
576  bool hasDescendantWith(const Cut& c, bool remove_duplicates=true) const;
577 
583  bool hasDescendantWithout(const ParticleSelector& f, bool remove_duplicates=true) const {
584  return hasDescendantWith([&](const Particle& p){ return !f(p); }, remove_duplicates);
585  }
591  bool hasDescendantWithout(const Cut& c, bool remove_duplicates=true) const;
592 
593 
598  Particles stableDescendants(const Cut& c=Cuts::OPEN) const;
599 
601  Particles stableDescendants(const ParticleSelector& f) const {
602  return filter_select(stableDescendants(), f);
603  }
604 
610  bool hasStableDescendantWith(const ParticleSelector& f) const {
611  return !stableDescendants(f).empty();
612  }
618  bool hasStableDescendantWith(const Cut& c) const;
619 
625  bool hasStableDescendantWithout(const ParticleSelector& f) const {
626  return hasStableDescendantWith([&](const Particle& p){ return !f(p); });
627  }
633  bool hasStableDescendantWithout(const Cut& c) const;
634 
635 
637  double flightLength() const;
638 
640 
641 
643 
644 
646  inline bool isFirstWith(const ParticleSelector& f) const {
647  if (!f(*this)) return false; //< This doesn't even meet f, let alone being the last to do so
648  if (any(parents(), f)) return false; //< If a direct parent has this property, this isn't the first
649  return true;
650  }
651 
653  inline bool isFirstWithout(const ParticleSelector& f) const {
654  return isFirstWith([&](const Particle& p){ return !f(p); });
655  }
656 
658  inline bool isLastWith(const ParticleSelector& f) const {
659  if (!f(*this)) return false; //< This doesn't even meet f, let alone being the last to do so
660  if (any(children(), f)) return false; //< If a child has this property, this isn't the last
661  return true;
662  }
663 
665  inline bool isLastWithout(const ParticleSelector& f) const {
666  return isLastWith([&](const Particle& p){ return !f(p); });
667  }
668 
670 
671 
672  protected:
673 
675  const GenParticle* _original;
676 
678  Particles _constituents;
679 
681  PdgId _id;
682 
684  FourMomentum _momentum;
685 
687  FourVector _origin;
688 
691  mutable std::pair<bool,bool> _isDirect;
692 
693  };
694 
695 
697 
698 
700  std::ostream& operator << (std::ostream& os, const Particle& p);
701 
703  std::ostream& operator << (std::ostream& os, const ParticlePair& pp);
704 
706 
707 
708 }
709 
710 
711 #include "Rivet/Tools/ParticleUtils.hh"
712 
713 #endif
Definition: ALICE_2010_I880049.cc:13
Particle(PdgId pid, const FourMomentum &mom, const FourVector &pos=FourVector())
Constructor without GenParticle.
Definition: Particle.hh:32
const Particles constituents(const ParticleSelector &selector) const
Direct constituents of this particle, filtered by a selection functor.
Definition: Particle.hh:252
virtual void addConstituents(const Particles &cs, bool addmom=false)
Add direct constituents to this particle.
Definition: Particle.cc:21
PdgId pdgId() const
Definition: Particle.hh:140
double flightLength() const
Flight length (divide by mm or cm to get the appropriate units)
Definition: Particle.cc:174
double pz() const
z component of momentum.
Definition: ParticleBase.hh:124
Particles children(const ParticleSelector &f) const
Get a list of the direct descendants from the current particle (with selector function) ...
Definition: Particle.hh:520
bool fromHadron() const
Determine whether the particle is from a hadron decay.
Definition: Particle.cc:229
Particles rawConstituents() const
Fundamental constituents of this particle.
Definition: Particle.cc:29
bool hasParent(PdgId pid) const
Definition: Particle.cc:184
Particles ancestors(const Cut &c=Cuts::OPEN, bool only_physical=true) const
Definition: Particle.cc:64
bool fromCharm() const
Determine whether the particle is from a c-hadron decay.
Definition: Particle.cc:223
const FourMomentum & mom() const
Get equivalent single momentum four-vector (const) (alias).
Definition: ParticleBase.hh:39
virtual void addConstituent(const Particle &c, bool addmom=false)
Add a single direct constituent to this particle.
Definition: Particle.cc:15
const Particles rawConstituents(const ParticleSorter &sorter) const
Fundamental constituents of this particle, sorted by a functor.
Definition: Particle.hh:269
Particle(const GenParticle &gp)
Constructor from a HepMC GenParticle reference.
Definition: Particle.hh:53
bool any(const CONTAINER &c)
Return true if x is true for any x in container c, otherwise false.
Definition: Utils.hh:287
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:653
const Particles constituents(const Cut &c, const ParticleSorter &sorter) const
Direct constituents of this particle, sorted by a functor.
Definition: Particle.hh:246
bool isParton() const
Is this a parton? (Hopefully not very often... fiducial FTW)
Definition: Particle.hh:204
bool hasDescendantWith(const ParticleSelector &f, bool remove_duplicates=true) const
Definition: Particle.hh:568
bool isStable() const
Whether this particle is stable according to the generator.
Definition: Particle.cc:57
int threeCharge() const
Definition: Particle.hh:159
PdgId pid() const
This Particle&#39;s PDG ID code.
Definition: Particle.hh:135
bool hasAncestorWith(const ParticleSelector &f, bool only_physical=true) const
Definition: Particle.hh:384
Particles children(const Cut &c=Cuts::OPEN) const
Get a list of the direct descendants from the current particle (with optional selection Cut) ...
Definition: Particle.cc:104
Particle representation, either from a HepMC::GenEvent or reconstructed.
Definition: Particle.hh:18
const Particles constituents(const Cut &c) const
Direct constituents of this particle, filtered by a Cut.
Definition: Particle.hh:240
Particle & transformBy(const LorentzTransform &lt)
Apply an active Lorentz transform to this particle.
Definition: Particle.cc:37
bool isHadron() const
Is this a hadron?
Definition: Particle.hh:174
Particle & setMomentum(double E, double px, double py, double pz)
Set the momentum via components.
Definition: Particle.hh:75
Particles parents(const ParticleSelector &f) const
Definition: Particle.hh:315
bool isCharged() const
Is this Particle charged?
Definition: Particle.hh:165
bool hasStableDescendantWith(const ParticleSelector &f) const
Definition: Particle.hh:610
Particle & setOrigin(double t, double x, double y, double z)
Set the origin position via components.
Definition: Particle.hh:99
bool hasChildWith(const ParticleSelector &f) const
Definition: Particle.hh:529
int abscharge3() const
Three times the absolute charge of this Particle (i.e. integer multiple of smallest quark charge)...
Definition: Particle.hh:162
double E() const
Get the energy directly.
Definition: ParticleBase.hh:51
const Particles & constituents() const
Direct constituents of this particle, returned by reference.
Definition: Particle.hh:230
Particles ancestors(const ParticleSelector &f, bool only_physical=true) const
Definition: Particle.hh:375
Particles stableDescendants(const ParticleSelector &f) const
Get a list of all the stable descendants from the current particle (with selector function) ...
Definition: Particle.hh:601
bool hasStableDescendantWithout(const ParticleSelector &f) const
Definition: Particle.hh:625
const GenParticle * genParticle() const
Get a const pointer to the original GenParticle.
Definition: Particle.hh:120
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:143
bool isChargedLepton() const
Is this a charged lepton?
Definition: Particle.hh:186
double charge() const
The charge of this Particle.
Definition: Particle.hh:149
bool fromPromptTau() const
Determine whether the particle is from a prompt tau decay.
Definition: Particle.hh:466
const Particles constituents(const ParticleSorter &sorter) const
Direct constituents of this particle, sorted by a functor.
Definition: Particle.hh:234
const Particles rawConstituents(const Cut &c) const
Fundamental constituents of this particle, filtered by a Cut.
Definition: Particle.hh:275
const FourVector & origin() const
The origin position.
Definition: Particle.hh:90
Object implementing Lorentz transform calculations and boosts.
Definition: LorentzTrans.hh:21
bool isBaryon() const
Is this a baryon?
Definition: Particle.hh:180
bool isVisible() const
Is this particle potentially visible in a detector?
Definition: Particle.cc:43
bool fromHadronicTau(bool prompt_taus_only=false) const
Determine whether the particle is from a tau which decayed hadronically.
Definition: Particle.cc:242
Specialisation of VectorN to a general (non-momentum) Lorentz 4-vector.
Definition: Vector4.hh:22
bool isNeutrino() const
Is this a neutrino?
Definition: Particle.hh:189
const Particles rawConstituents(const ParticleSelector &selector, const ParticleSorter &sorter) const
Fundamental constituents of this particle, filtered and sorted by functors.
Definition: Particle.hh:293
const FourMomentum & momentum() const
The momentum.
Definition: Particle.hh:64
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:559
Base class for particle-like things like Particle and Jet.
Definition: ParticleBase.hh:13
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:658
bool hasDescendantWithout(const ParticleSelector &f, bool remove_duplicates=true) const
Definition: Particle.hh:583
bool isMeson() const
Is this a meson?
Definition: Particle.hh:177
virtual void setConstituents(const Particles &cs, bool setmom=false)
Set direct constituents of this particle.
Definition: Particle.cc:9
bool isComposite() const
Determine if this Particle is a composite of other Rivet Particles.
Definition: Particle.hh:223
int charge3() const
Three times the charge of this Particle (i.e. integer multiple of smallest quark charge).
Definition: Particle.hh:155
Particle()
Definition: Particle.hh:26
Particle(const GenParticle *gp)
Constructor from a HepMC GenParticle pointer.
Definition: Particle.hh:40
MOMS & sortBy(MOMS &pbs, const CMP &cmp)
Sort a container of momenta by cmp and return by reference for non-const inputs.
Definition: Vector4.hh:1431
bool hasAncestorWithout(const ParticleSelector &f, bool only_physical=true) const
Definition: Particle.hh:399
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:665
double px() const
x component of momentum.
Definition: ParticleBase.hh:120
bool isLepton() const
Is this a lepton?
Definition: Particle.hh:183
const Particles rawConstituents(const ParticleSelector &selector) const
Fundamental constituents of this particle, filtered by a selection functor.
Definition: Particle.hh:287
const Particles rawConstituents(const Cut &c, const ParticleSorter &sorter) const
Fundamental constituents of this particle, sorted by a functor.
Definition: Particle.hh:281
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:646
Particle & setMomentum(const FourMomentum &momentum)
Set the momentum.
Definition: Particle.hh:69
bool hasChildWithout(const ParticleSelector &f) const
Definition: Particle.hh:544
double abscharge() const
The absolute charge of this Particle.
Definition: Particle.hh:152
PdgId abspid() const
Absolute value of the PDG ID code.
Definition: Particle.hh:137
bool fromTau(bool prompt_taus_only=false) const
Determine whether the particle is from a tau decay.
Definition: Particle.cc:235
bool hasBottom() const
Does this (hadron) contain a b quark?
Definition: Particle.hh:192
bool isDirect(bool allow_from_direct_tau=false, bool allow_from_direct_mu=false) const
Shorthand definition of &#39;promptness&#39; based on set definition flags.
Definition: Particle.cc:249
bool hasParentWithout(const ParticleSelector &f) const
Definition: Particle.hh:339
Particles stableDescendants(const Cut &c=Cuts::OPEN) const
Definition: Particle.cc:155
const Particles constituents(const ParticleSelector &selector, const ParticleSorter &sorter) const
Direct constituents of this particle, filtered and sorted by functors.
Definition: Particle.hh:258
Specialized version of the FourVector with momentum/energy functionality.
Definition: Vector4.hh:301
bool hasCharm() const
Does this (hadron) contain a c quark?
Definition: Particle.hh:195
double p() const
Get the 3-momentum magnitude directly.
Definition: ParticleBase.hh:110
bool isPrompt(bool allow_from_prompt_tau=false, bool allow_from_prompt_mu=false) const
Alias for isDirect.
Definition: Particle.hh:500
Particles parents(const Cut &c=Cuts::OPEN) const
Definition: Particle.cc:87
Particle & setOrigin(const FourVector &position)
Set the origin position.
Definition: Particle.hh:94
double py() const
y component of momentum.
Definition: ParticleBase.hh:122
bool hasAncestor(PdgId pid, bool only_physical=true) const
Definition: Particle.cc:193
virtual fastjet::PseudoJet pseudojet() const
Converter to FastJet3 PseudoJet.
Definition: Particle.hh:111
bool fromDecay() const
Determine whether the particle is from a hadron or tau decay.
Definition: Particle.hh:485
bool hasParentWith(const ParticleSelector &f) const
Definition: Particle.hh:324
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) ...
Definition: Particle.cc:127
bool fromBottom() const
Determine whether the particle is from a b-hadron decay.
Definition: Particle.cc:217