rivet is hosted by Hepforge, IPPP Durham
Jet Class Reference

Representation of a clustered jet of particles. More...

#include <Jet.hh>

List of all members.

Public Member Functions

Constructors
 Jet (const fastjet::PseudoJet &pj, const Particles &particles=Particles(), const Particles &tags=Particles())
 Constructor from a FastJet PseudoJet, with optional full particle constituents information.
 Jet (const FourMomentum &pjet, const Particles &particles=Particles(), const Particles &tags=Particles())
 Set the jet data, with optional full particle information.
 Jet (const Particles &particles, const FourMomentum &pjet)
 Jet ()
 Default constructor -- only for STL storability.
Access jet constituents
size_t size () const
 Number of particles in this jet.
Particlesparticles ()
 Get the particles in this jet.
const Particlesparticles () const
 Get the particles in this jet (const version)
const Particles particles (const Cut &c) const
 Get the particles in this jet which pass a cut (const)
Particlesconstituents ()
 Get the particles in this jet (FastJet-like alias)
const Particlesconstituents () const
 Get the particles in this jet (FastJet-like alias, const version)
const Particles constituents (const Cut &c) const
 Get the particles in this jet which pass a cut (FastJet-like alias, const)
bool containsParticle (const Particle &particle) const
 Check whether this jet contains a particular particle.
bool containsPID (const Particle &particle) const
 Nicer alias for containsParticleId.
bool containsParticleId (PdgId pid) const
 Check whether this jet contains a certain particle type.
bool containsPID (PdgId pid) const
 Nicer alias for containsParticleId.
bool containsParticleId (const vector< PdgId > &pids) const
 Check whether this jet contains at least one of certain particle types.
bool containsPID (const vector< PdgId > &pids) const
 Nicer alias for containsParticleId.
Tagging
Note:
General sources of tag particles are planned. The default jet finding adds b-hadron, c-hadron, and tau tags by ghost association.
Particlestags ()
 Particles which have been tag-matched to this jet.
const Particlestags () const
 Particles which have been tag-matched to this jet (const version)
template<typename FN >
Particles tags (const FN &f) const
 Particles which have been tag-matched to this jet _and_ pass a selector function.
Particles tags (const Cut &c) const
 Particles which have been tag-matched to this jet _and_ pass a Cut.
Particles bTags (const Cut &c=Cuts::open()) const
 b particles which have been tag-matched to this jet (and pass an optional Cut)
template<typename FN >
Particles bTags (const FN &f) const
 b particles which have been tag-matched to this jet _and_ pass a selector function
bool bTagged (const Cut &c=Cuts::open()) const
 Does this jet have at least one b-tag (that passes an optional Cut)?
template<typename FN >
Particles bTagged (const FN &f) const
 Does this jet have at least one b-tag (that passes the supplied selector function)?
Particles cTags (const Cut &c=Cuts::open()) const
 c (and not b) particles which have been tag-matched to this jet (and pass an optional Cut)
template<typename FN >
Particles cTags (const FN &f) const
 c (and not b) particles which have been tag-matched to this jet and pass a selector function
bool cTagged (const Cut &c=Cuts::open()) const
 Does this jet have at least one c-tag (that passes an optional Cut)?
template<typename FN >
bool cTagged (const FN &f) const
 Does this jet have at least one c-tag (that passes the supplied selector function)?
Particles tauTags (const Cut &c=Cuts::open()) const
 Tau particles which have been tag-matched to this jet (and pass an optional Cut)
template<typename FN >
Particles tauTags (const FN &f) const
 Tau particles which have been tag-matched to this jet and pass a selector function.
bool tauTagged (const Cut &c=Cuts::open()) const
 Does this jet have at least one tau-tag (that passes an optional Cut)?
template<typename FN >
bool tauTagged (const FN &f) const
 Does this jet have at least one tau-tag (that passes the supplied selector function)?
bool containsBottom (bool include_decay_products=true) const
 Check whether this jet contains a bottom-flavoured hadron.
bool containsCharm (bool include_decay_products=true) const
 Check whether this jet contains a charm-flavoured hadron.
Effective jet 4-vector properties
const FourMomentummomentum () const
 Get equivalent single momentum four-vector.
JettransformBy (const LorentzTransform &lt)
double totalEnergy () const
 Get the total energy of this jet.
double neutralEnergy () const
 Get the energy carried in this jet by neutral particles.
double hadronicEnergy () const
 Get the energy carried in this jet by hadrons.
Interaction with FastJet
const fastjet::PseudoJet & pseudojet () const
 Access the internal FastJet3 PseudoJet (as a const reference)
 operator const fastjet::PseudoJet & () const
 Cast operator to FastJet3 PseudoJet (as a const reference)
Set the jet constituents and properties
JetsetState (const fastjet::PseudoJet &pj, const Particles &particles=Particles(), const Particles &tags=Particles())
 Set the jet data from a FastJet PseudoJet, with optional particle constituents and tags lists.
JetsetState (const FourMomentum &mom, const Particles &particles, const Particles &tags=Particles())
 Set all the jet data, with optional full particle constituent and tag information.
JetsetState (const Particles &particles, const FourMomentum &mom)
JetsetParticles (const Particles &particles)
 Set the particles collection with full particle information.
JetsetConstituents (const Particles &particles)
Jetclear ()
 Reset this jet as empty.
Effective momentum accessors
const FourMomentummom () const
 Get equivalent single momentum four-vector (const) (alias).
 operator const FourMomentum & () const
 Cast operator for conversion to FourMomentum.
Convenience access to the effective 4-vector properties
double E () const
 Get the energy directly.
double energy () const
 Get the energy directly (alias).
double E2 () const
 Get the energy-squared.
double energy2 () const
 Get the energy-squared (alias).
double pt () const
 Get the $ p_T $ directly.
double pT () const
 Get the $ p_T $ directly (alias).
double perp () const
 Get the $ p_T $ directly (alias).
double pt2 () const
 Get the $ p_T^2 $ directly.
double pT2 () const
 Get the $ p_T^2 $ directly (alias).
double perp2 () const
 Get the $ p_T^2 $ directly (alias).
double Et () const
 Get the $ E_T $ directly.
double Et2 () const
 Get the $ E_T^2 $ directly.
double mass () const
 Get the mass directly.
double mass2 () const
 Get the mass**2 directly.
double pseudorapidity () const
 Get the $ \eta $ directly.
double eta () const
 Get the $ \eta $ directly (alias).
double abspseudorapidity () const
 Get the $ |\eta| $ directly.
double abseta () const
 Get the $ |\eta| $ directly (alias).
double rapidity () const
 Get the $ y $ directly.
double rap () const
 Get the $ y $ directly (alias).
double absrapidity () const
 Get the $ |y| $ directly.
double absrap () const
 Get the $ |y| $ directly (alias).
double azimuthalAngle (const PhiMapping mapping=ZERO_2PI) const
 Azimuthal angle $ \phi $.
double phi (const PhiMapping mapping=ZERO_2PI) const
 Get the $ \phi $ directly.
Vector3 p3 () const
 Get the 3-momentum directly.
double p () const
 Get the 3-momentum magnitude directly.
double p2 () const
 Get the 3-momentum magnitude-squared directly.
double px () const
 x component of momentum.
double py () const
 y component of momentum.
double pz () const
 z component of momentum.
double px2 () const
 x component of momentum, squared.
double py2 () const
 y component of momentum, squared.
double pz2 () const
 z component of momentum, squared.
double polarAngle () const
 Angle subtended by the 3-vector and the z-axis.
double theta () const
 Synonym for polarAngle.
double angle (const ParticleBase &v) const
 Angle between this vector and another.
double angle (const FourVector &v) const
 Angle between this vector and another.
double angle (const Vector3 &v3) const
 Angle between this vector and another (3-vector)

Private Attributes

fastjet::PseudoJet _pseudojet
 FJ3 PseudoJet member to unify PseudoJet and Jet.
Particles _particles
Particles _tags
 Particles used to tag this jet (can be anything, but c and b hadrons are the most common)
FourMomentum _momentum
 Effective jet 4-vector (just for caching)

Detailed Description

Representation of a clustered jet of particles.

Definition at line 18 of file Jet.hh.


Constructor & Destructor Documentation

Jet ( const fastjet::PseudoJet &  pj,
const Particles particles = Particles(),
const Particles tags = Particles() 
) [inline]

Constructor from a FastJet PseudoJet, with optional full particle constituents information.

Definition at line 25 of file Jet.hh.

                                                                                                               {
      setState(pj, particles, tags);
    }
Jet ( const FourMomentum pjet,
const Particles particles = Particles(),
const Particles tags = Particles() 
) [inline]

Set the jet data, with optional full particle information.

Definition at line 30 of file Jet.hh.

                                                                                                             {
      setState(pjet, particles, tags);
    }
Jet ( const Particles particles,
const FourMomentum pjet 
) [inline]

Set all the jet data, with full particle information.

Deprecated:
Prefer the form where the 4-vec comes first and the particles list is optional.

Definition at line 37 of file Jet.hh.

                                                              {
      setState(pjet, particles);
    }
Jet ( ) [inline]

Default constructor -- only for STL storability.

Definition at line 42 of file Jet.hh.

{ clear(); }

Member Function Documentation

double abseta ( ) const [inline, inherited]

Get the $ |\eta| $ directly (alias).

Definition at line 81 of file ParticleBase.hh.

{ return momentum().abseta(); }
double abspseudorapidity ( ) const [inline, inherited]

Get the $ |\eta| $ directly.

Definition at line 79 of file ParticleBase.hh.

{ return momentum().abspseudorapidity(); }
double absrap ( ) const [inline, inherited]

Get the $ |y| $ directly (alias).

Definition at line 90 of file ParticleBase.hh.

{ return momentum().absrap(); }
double absrapidity ( ) const [inline, inherited]

Get the $ |y| $ directly.

Definition at line 88 of file ParticleBase.hh.

{ return momentum().absrapidity(); }
double angle ( const ParticleBase v) const [inline, inherited]

Angle between this vector and another.

Definition at line 125 of file ParticleBase.hh.

{ return momentum().angle(v.momentum()); }
double angle ( const FourVector v) const [inline, inherited]

Angle between this vector and another.

Definition at line 127 of file ParticleBase.hh.

{ return momentum().angle(v); }
double angle ( const Vector3 v3) const [inline, inherited]

Angle between this vector and another (3-vector)

Definition at line 129 of file ParticleBase.hh.

{ return momentum().angle(v3); }
double azimuthalAngle ( const PhiMapping  mapping = ZERO_2PI) const [inline, inherited]

Azimuthal angle $ \phi $.

Definition at line 93 of file ParticleBase.hh.

{ return momentum().azimuthalAngle(mapping); }
bool bTagged ( const Cut c = Cuts::open()) const [inline]

Does this jet have at least one b-tag (that passes an optional Cut)?

Definition at line 115 of file Jet.hh.

{ return !bTags(c).empty(); }
Particles bTagged ( const FN &  f) const [inline]

Does this jet have at least one b-tag (that passes the supplied selector function)?

Definition at line 118 of file Jet.hh.

{ return !bTags(f).empty(); }
Particles bTags ( const Cut c = Cuts::open()) const

b particles which have been tag-matched to this jet (and pass an optional Cut)

The default jet finding adds b-hadron tags by ghost association.

Definition at line 166 of file Jet.cc.

                                         {
    Particles rtn;
    for (const Particle& tp : tags()) {
      if (hasBottom(tp) && c->accept(tp)) rtn.push_back(tp);
    }
    return rtn;
  }
Particles bTags ( const FN &  f) const [inline]

b particles which have been tag-matched to this jet _and_ pass a selector function

Definition at line 112 of file Jet.hh.

{ return filter_select(bTags(), f); }
Jet & clear ( )

Reset this jet as empty.

Definition at line 10 of file Jet.cc.

                  {
    _momentum = FourMomentum();
    _pseudojet.reset(0,0,0,0);
    _particles.clear();
    return *this;
  }
Particles& constituents ( ) [inline]

Get the particles in this jet (FastJet-like alias)

Definition at line 61 of file Jet.hh.

{ return particles(); }
const Particles& constituents ( ) const [inline]

Get the particles in this jet (FastJet-like alias, const version)

Definition at line 63 of file Jet.hh.

{ return particles(); }
const Particles constituents ( const Cut c) const [inline]

Get the particles in this jet which pass a cut (FastJet-like alias, const)

Definition at line 65 of file Jet.hh.

{ return particles(c); }
bool containsBottom ( bool  include_decay_products = true) const

Check whether this jet contains a bottom-flavoured hadron.

Deprecated:
The bTags() or bTagged() function is probably what you want for tagging. This one ignores the tags() list and draws conclusions based directly on the jet constituents; the other gives a much better match to typical experimental methods.
Note:
The decision is made by first trying to find a bottom-flavoured particle in the particles list. Most likely this will fail unless bottom hadrons are set stable. If include_decay_products is true (the default), a fallback is attempted, using the post-hadronization ancestor history of all constituents.

Definition at line 143 of file Jet.cc.

                                                            {
    foreach (const Particle& p, particles()) {
      const PdgId pid = p.pid();
      if (abs(pid) == PID::BQUARK) return true;
      if (PID::isHadron(pid) && PID::hasBottom(pid)) return true;
      if (include_decay_products) {
        const HepMC::GenVertex* gv = p.genParticle()->production_vertex();
        if (gv) {
          foreach (const GenParticle* pi, Rivet::particles(gv, HepMC::ancestors)) {
            const PdgId pid2 = pi->pdg_id();
            if (PID::isHadron(pid2) && PID::hasBottom(pid2)) return true;
          }
        }
      }
    }
    return false;
  }
bool containsCharm ( bool  include_decay_products = true) const

Check whether this jet contains a charm-flavoured hadron.

Deprecated:
The cTags() or cTagged() function is probably what you want for tagging. This one ignores the tags() list and draws conclusions based directly on the jet constituents; the other gives a much better match to typical experimental methods.
Note:
The decision is made by first trying to find a charm-flavoured particle in the particles list. Most likely this will fail unless charmed hadrons are set stable. If include_decay_products is true (the default), a fallback is attempted, using the post-hadronization ancestor history of all constituents.

Definition at line 124 of file Jet.cc.

                                                           {
    foreach (const Particle& p, particles()) {
      const PdgId pid = p.pid();
      if (abs(pid) == PID::CQUARK) return true;
      if (PID::isHadron(pid) && PID::hasCharm(pid)) return true;
      if (include_decay_products) {
        const HepMC::GenVertex* gv = p.genParticle()->production_vertex();
        if (gv) {
          foreach (const GenParticle* pi, Rivet::particles(gv, HepMC::ancestors)) {
            const PdgId pid2 = pi->pdg_id();
            if (PID::isHadron(pid2) && PID::hasCharm(pid2)) return true;
          }
        }
      }
    }
    return false;
  }
bool containsParticle ( const Particle particle) const

Check whether this jet contains a particular particle.

Definition at line 61 of file Jet.cc.

                                                           {
    const int barcode = particle.genParticle()->barcode();
    foreach (const Particle& p, particles()) {
      if (p.genParticle()->barcode() == barcode) return true;
    }
    return false;
  }
bool containsParticleId ( PdgId  pid) const

Check whether this jet contains a certain particle type.

Definition at line 70 of file Jet.cc.

                                              {
    foreach (const Particle& p, particles()) {
      if (p.pid() == pid) return true;
    }
    return false;
  }
bool containsParticleId ( const vector< PdgId > &  pids) const

Check whether this jet contains at least one of certain particle types.

Definition at line 78 of file Jet.cc.

                                                              {
    foreach (const Particle& p, particles()) {
      foreach (PdgId pid, pids) {
        if (p.pid() == pid) return true;
      }
    }
    return false;
  }
bool containsPID ( const Particle particle) const [inline]

Nicer alias for containsParticleId.

Definition at line 70 of file Jet.hh.

{ return containsParticle(particle); }
bool containsPID ( PdgId  pid) const [inline]

Nicer alias for containsParticleId.

Definition at line 75 of file Jet.hh.

{ return containsParticleId(pid); }
bool containsPID ( const vector< PdgId > &  pids) const [inline]

Nicer alias for containsParticleId.

Definition at line 80 of file Jet.hh.

{ return containsParticleId(pids); }
bool cTagged ( const Cut c = Cuts::open()) const [inline]

Does this jet have at least one c-tag (that passes an optional Cut)?

Definition at line 130 of file Jet.hh.

{ return !cTags(c).empty(); }
bool cTagged ( const FN &  f) const [inline]

Does this jet have at least one c-tag (that passes the supplied selector function)?

Definition at line 133 of file Jet.hh.

{ return !cTags(f).empty(); }
Particles cTags ( const Cut c = Cuts::open()) const

c (and not b) particles which have been tag-matched to this jet (and pass an optional Cut)

The default jet finding adds c-hadron tags by ghost association.

Todo:
Is making b and c tags exclusive the right thing to do?

Definition at line 174 of file Jet.cc.

                                         {
    Particles rtn;
    for (const Particle& tp : tags()) {
      /// @todo Is making b and c tags exclusive the right thing to do?
      if (hasCharm(tp) && !hasBottom(tp) && c->accept(tp)) rtn.push_back(tp);
    }
    return rtn;
  }
Particles cTags ( const FN &  f) const [inline]

c (and not b) particles which have been tag-matched to this jet and pass a selector function

Definition at line 127 of file Jet.hh.

{ return filter_select(cTags(), f); }
double E ( ) const [inline, inherited]

Get the energy directly.

Definition at line 41 of file ParticleBase.hh.

{ return momentum().E(); }
double E2 ( ) const [inline, inherited]

Get the energy-squared.

Definition at line 46 of file ParticleBase.hh.

{ return momentum().E2(); }
double energy ( ) const [inline, inherited]

Get the energy directly (alias).

Definition at line 43 of file ParticleBase.hh.

{ return momentum().E(); }
double energy2 ( ) const [inline, inherited]

Get the energy-squared (alias).

Definition at line 48 of file ParticleBase.hh.

{ return momentum().E2(); }
double Et ( ) const [inline, inherited]

Get the $ E_T $ directly.

Definition at line 65 of file ParticleBase.hh.

{ return momentum().Et(); }
double Et2 ( ) const [inline, inherited]

Get the $ E_T^2 $ directly.

Definition at line 67 of file ParticleBase.hh.

{ return momentum().Et2(); }
double eta ( ) const [inline, inherited]

Get the $ \eta $ directly (alias).

Definition at line 77 of file ParticleBase.hh.

{ return momentum().eta(); }
double hadronicEnergy ( ) const

Get the energy carried in this jet by hadrons.

Definition at line 112 of file Jet.cc.

                                   {
    double e_hadr = 0.0;
    foreach (const Particle& p, particles()) {
      const PdgId pid = p.pid();
      if (PID::isHadron(pid)) {
        e_hadr += p.E();
      }
    }
    return e_hadr;
  }
double mass ( ) const [inline, inherited]

Get the mass directly.

Definition at line 70 of file ParticleBase.hh.

{ return momentum().mass(); }
double mass2 ( ) const [inline, inherited]

Get the mass**2 directly.

Definition at line 72 of file ParticleBase.hh.

{ return momentum().mass2(); }
const FourMomentum& mom ( ) const [inline, inherited]

Get equivalent single momentum four-vector (const) (alias).

Definition at line 29 of file ParticleBase.hh.

{ return momentum(); };
const FourMomentum& momentum ( ) const [inline, virtual]

Get equivalent single momentum four-vector.

Implements ParticleBase.

Definition at line 188 of file Jet.hh.

{ return _momentum; }
double neutralEnergy ( ) const

Get the energy carried in this jet by neutral particles.

Definition at line 100 of file Jet.cc.

                                  {
    double e_neutral = 0.0;
    foreach (const Particle& p, particles()) {
      const PdgId pid = p.pid();
      if (PID::threeCharge(pid) == 0) {
        e_neutral += p.E();
      }
    }
    return e_neutral;
  }
operator const fastjet::PseudoJet & ( ) const [inline]

Cast operator to FastJet3 PseudoJet (as a const reference)

Definition at line 214 of file Jet.hh.

{ return pseudojet(); }
operator const FourMomentum & ( ) const [inline, inherited]

Cast operator for conversion to FourMomentum.

Definition at line 32 of file ParticleBase.hh.

{ return momentum(); }
double p ( ) const [inline, inherited]

Get the 3-momentum magnitude directly.

Definition at line 101 of file ParticleBase.hh.

{ return momentum().p(); }
double p2 ( ) const [inline, inherited]

Get the 3-momentum magnitude-squared directly.

Definition at line 103 of file ParticleBase.hh.

{ return momentum().p2(); }
Vector3 p3 ( ) const [inline, inherited]

Get the 3-momentum directly.

Definition at line 99 of file ParticleBase.hh.

{ return momentum().vector3(); }
Particles& particles ( ) [inline]

Get the particles in this jet.

Definition at line 54 of file Jet.hh.

{ return _particles; }
const Particles& particles ( ) const [inline]

Get the particles in this jet (const version)

Definition at line 56 of file Jet.hh.

{ return _particles; }
const Particles particles ( const Cut c) const [inline]

Get the particles in this jet which pass a cut (const)

Definition at line 58 of file Jet.hh.

{ return filterBy(_particles, c); }
double perp ( ) const [inline, inherited]

Get the $ p_T $ directly (alias).

Definition at line 55 of file ParticleBase.hh.

{ return pt(); }
double perp2 ( ) const [inline, inherited]

Get the $ p_T^2 $ directly (alias).

Definition at line 62 of file ParticleBase.hh.

{ return pt2(); }
double phi ( const PhiMapping  mapping = ZERO_2PI) const [inline, inherited]

Get the $ \phi $ directly.

Definition at line 95 of file ParticleBase.hh.

{ return momentum().phi(mapping); }
double polarAngle ( ) const [inline, inherited]

Angle subtended by the 3-vector and the z-axis.

Definition at line 120 of file ParticleBase.hh.

{ return momentum().polarAngle(); }
const fastjet::PseudoJet& pseudojet ( ) const [inline]

Access the internal FastJet3 PseudoJet (as a const reference)

Definition at line 211 of file Jet.hh.

{ return _pseudojet; }
double pseudorapidity ( ) const [inline, inherited]

Get the $ \eta $ directly.

Definition at line 75 of file ParticleBase.hh.

{ return momentum().eta(); }
double pt ( ) const [inline, inherited]

Get the $ p_T $ directly.

Definition at line 51 of file ParticleBase.hh.

{ return momentum().pt(); }
double pT ( ) const [inline, inherited]

Get the $ p_T $ directly (alias).

Definition at line 53 of file ParticleBase.hh.

{ return pt(); }
double pt2 ( ) const [inline, inherited]

Get the $ p_T^2 $ directly.

Definition at line 58 of file ParticleBase.hh.

{ return momentum().pt2(); }
double pT2 ( ) const [inline, inherited]

Get the $ p_T^2 $ directly (alias).

Definition at line 60 of file ParticleBase.hh.

{ return pt2(); }
double px ( ) const [inline, inherited]

x component of momentum.

Definition at line 106 of file ParticleBase.hh.

{ return momentum().x(); }
double px2 ( ) const [inline, inherited]

x component of momentum, squared.

Definition at line 113 of file ParticleBase.hh.

{ return momentum().x2(); }
double py ( ) const [inline, inherited]

y component of momentum.

Definition at line 108 of file ParticleBase.hh.

{ return momentum().y(); }
double py2 ( ) const [inline, inherited]

y component of momentum, squared.

Definition at line 115 of file ParticleBase.hh.

{ return momentum().y2(); }
double pz ( ) const [inline, inherited]

z component of momentum.

Definition at line 110 of file ParticleBase.hh.

{ return momentum().z(); }
double pz2 ( ) const [inline, inherited]

z component of momentum, squared.

Definition at line 117 of file ParticleBase.hh.

{ return momentum().z2(); }
double rap ( ) const [inline, inherited]

Get the $ y $ directly (alias).

Definition at line 86 of file ParticleBase.hh.

{ return momentum().rapidity(); }
double rapidity ( ) const [inline, inherited]

Get the $ y $ directly.

Definition at line 84 of file ParticleBase.hh.

{ return momentum().rapidity(); }
Jet& setConstituents ( const Particles particles) [inline]

Definition at line 239 of file Jet.hh.

{ return setParticles(particles); }
Jet & setParticles ( const Particles particles)

Set the particles collection with full particle information.

If set, this overrides particle info extracted from the PseudoJet

Definition at line 55 of file Jet.cc.

                                                   {
    _particles = particles;
    return *this;
  }
Jet & setState ( const fastjet::PseudoJet &  pj,
const Particles particles = Particles(),
const Particles tags = Particles() 
)

Set the jet data from a FastJet PseudoJet, with optional particle constituents and tags lists.

Note:
The particles() list will be extracted from PseudoJet constituents by default, making use of an attached user info if one is found.

Definition at line 28 of file Jet.cc.

                                                                                                  {
    clear();
    _pseudojet = pj;
    _momentum = FourMomentum(pj.e(), pj.px(), pj.py(), pj.pz());
    _particles = particles;
    _tags = tags;
    // if (_particles.empty()) {
    //   foreach (const fastjet::PseudoJet pjc, _pseudojet.constituents()) {
    //     // If there is no attached user info, we can't create a meaningful particle, so skip
    //     if (!pjc.has_user_info<RivetFJInfo>()) continue;
    //     const RivetFJInfo& fjinfo = pjc.user_info<RivetFJInfo>();
    //     // Don't add ghosts to the particles list
    //     if (fjinfo.isGhost) continue;
    //     // Otherwise construct a Particle from the PseudoJet, preferably from an associated GenParticle
    //     ?if (fjinfo.genParticle != NULL) {
    //       _particles.push_back(Particle(fjinfo.genParticle));
    //     } else {
    //       if (fjinfo.pid == 0) continue; // skip if there is a null PID entry in the FJ info
    //       const FourMomentum pjcmom(pjc.e(), pjc.px(), pjc.py(), pjc.pz());
    //       _particles.push_back(Particle(fjinfo.pid, pjcmom));
    //     }
    //   }
    // }
    return *this;
  }
Jet & setState ( const FourMomentum mom,
const Particles particles,
const Particles tags = Particles() 
)

Set all the jet data, with optional full particle constituent and tag information.

Definition at line 18 of file Jet.cc.

                                                                                               {
    clear();
    _momentum = mom;
    _pseudojet = fastjet::PseudoJet(mom.px(), mom.py(), mom.pz(), mom.E());
    _particles = particles;
    _tags = tags;
    return *this;
  }
Jet& setState ( const Particles particles,
const FourMomentum mom 
) [inline]
Deprecated:
Prefer the 4-mom first-arg versions. Remove in Rivet v3

Definition at line 233 of file Jet.hh.

{ return setState(mom, particles); }
size_t size ( ) const [inline]

Number of particles in this jet.

Definition at line 51 of file Jet.hh.

{ return _particles.size(); }
Particles& tags ( ) [inline]

Particles which have been tag-matched to this jet.

Definition at line 92 of file Jet.hh.

{ return _tags; }
const Particles& tags ( ) const [inline]

Particles which have been tag-matched to this jet (const version)

Definition at line 94 of file Jet.hh.

{ return _tags; }
Particles tags ( const FN &  f) const [inline]

Particles which have been tag-matched to this jet _and_ pass a selector function.

Note:
Note the less efficient return by value, due to the filtering.

Definition at line 99 of file Jet.hh.

{ return filter_select(tags(), f); }
Particles tags ( const Cut c) const

Particles which have been tag-matched to this jet _and_ pass a Cut.

Note:
Note the less efficient return by value, due to the cut-pass filtering.

Definition at line 162 of file Jet.cc.

                                        {
    return filter_select(tags(), c);
  }
bool tauTagged ( const Cut c = Cuts::open()) const [inline]

Does this jet have at least one tau-tag (that passes an optional Cut)?

Definition at line 145 of file Jet.hh.

{ return !tauTags(c).empty(); }
bool tauTagged ( const FN &  f) const [inline]

Does this jet have at least one tau-tag (that passes the supplied selector function)?

Definition at line 148 of file Jet.hh.

{ return !tauTags(f).empty(); }
Particles tauTags ( const Cut c = Cuts::open()) const

Tau particles which have been tag-matched to this jet (and pass an optional Cut)

The default jet finding adds tau tags by ghost association.

Definition at line 183 of file Jet.cc.

                                           {
    Particles rtn;
    for (const Particle& tp : tags()) {
      if (isTau(tp) && c->accept(tp)) rtn.push_back(tp);
    }
    return rtn;
  }
Particles tauTags ( const FN &  f) const [inline]

Tau particles which have been tag-matched to this jet and pass a selector function.

Definition at line 142 of file Jet.hh.

{ return filter_select(tauTags(), f); }
double theta ( ) const [inline, inherited]

Synonym for polarAngle.

Definition at line 122 of file ParticleBase.hh.

{ return momentum().theta(); }
double totalEnergy ( ) const [inline]

Get the total energy of this jet.

Definition at line 196 of file Jet.hh.

{ return momentum().E(); }
Jet & transformBy ( const LorentzTransform lt)

Apply an active Lorentz transform to this jet

Note:
The Rivet jet momentum, constituent particles, and tag particles will be modified.
Warning:
The FastJet cluster sequence and pseudojets will not be modified: don't use them after transformation!
Todo:
Jet::containsMatch(Matcher m) { ... if m(pid) return true; ... }

Definition at line 91 of file Jet.cc.

                                                  {
    _momentum = lt.transform(_momentum);
    for (Particle& p : _particles) p.transformBy(lt);
    for (Particle& t : _tags) t.transformBy(lt);
    _pseudojet.reset(_momentum.px(), _momentum.py(), _momentum.pz(), _momentum.E()); //< lose ClusterSeq etc.
    return *this;
  }

Member Data Documentation

FourMomentum _momentum [mutable, private]

Effective jet 4-vector (just for caching)

Definition at line 260 of file Jet.hh.

Particles _particles [private]

Full constituent particle information. (Filled from PseudoJet if possible.)

Todo:
Make these mutable or similar? Add a flag to force a cache rebuild?

Definition at line 254 of file Jet.hh.

fastjet::PseudoJet _pseudojet [private]

FJ3 PseudoJet member to unify PseudoJet and Jet.

Definition at line 250 of file Jet.hh.

Particles _tags [private]

Particles used to tag this jet (can be anything, but c and b hadrons are the most common)

Definition at line 257 of file Jet.hh.


The documentation for this class was generated from the following files: