rivet is hosted by Hepforge, IPPP Durham

Project out jets found using the FastJet package jet algorithms. More...

#include <FastJets.hh>

List of all members.

Public Types

enum  JetAlgName {
  KT, CAM, SISCONE, ANTIKT,
  ATLASCONE, CMSCONE, CDFJETCLU, CDFMIDPOINT,
  D0ILCONE, JADE, DURHAM, TRACKJET
}
enum  MuonsStrategy { NO_MUONS, DECAY_MUONS, ALL_MUONS }
 Enum for the treatment of muons: whether to include all, some, or none in jet-finding. More...
enum  InvisiblesStrategy { NO_INVISIBLES, DECAY_INVISIBLES, ALL_INVISIBLES }
 Enum for the treatment of invisible particles: whether to include all, some, or none in jet-finding. More...
typedef Jet entity_type
typedef Jets collection_type

Public Member Functions

void reset ()
 Reset the projection. Jet def, etc. are unchanged.
void useJetArea (fastjet::AreaDefinition *adef)
 Use provided jet area definition.
void calc (const Particles &fsparticles, const Particles &tagparticles=Particles())
 Do the calculation locally (no caching).
virtual unique_ptr< Projectionclone () const =0
 Clone on the heap.
size_t numJets (const Cut &c=Cuts::open()) const
 Number of jets passing the provided Cut.
size_t size () const
 Number of jets (without cuts).
bool empty () const
 Whether the inclusive jet collection is empty.
collection_type entities () const
 Template-usable interface common to FinalState.
bool before (const Projection &p) const
virtual const std::set< PdgIdPairbeamPairs () const
virtual std::string name () const
 Get the name of the projection.
ProjectionaddPdgIdPair (PdgId beam1, PdgId beam2)
 Add a colliding beam pair.
LoggetLog () const
 Get a Log object based on the getName() property of the calling projection object.
void setName (const std::string &name)
 Used by derived classes to set their name.
void markAsOwned () const
Constructors etc.
 FastJets (const FinalState &fsp, const fastjet::JetDefinition &jdef, JetAlg::MuonsStrategy usemuons=JetAlg::ALL_MUONS, JetAlg::InvisiblesStrategy useinvis=JetAlg::NO_INVISIBLES, fastjet::AreaDefinition *adef=nullptr)
 FastJets (const FinalState &fsp, const fastjet::JetDefinition &jdef, fastjet::AreaDefinition *adef, JetAlg::MuonsStrategy usemuons=JetAlg::ALL_MUONS, JetAlg::InvisiblesStrategy useinvis=JetAlg::NO_INVISIBLES)
 FastJets (const FinalState &fsp, fastjet::JetAlgorithm type, fastjet::RecombinationScheme recom, double rparameter, JetAlg::MuonsStrategy usemuons=JetAlg::ALL_MUONS, JetAlg::InvisiblesStrategy useinvis=JetAlg::NO_INVISIBLES, fastjet::AreaDefinition *adef=nullptr)
 FastJets (const FinalState &fsp, fastjet::JetAlgorithm type, fastjet::RecombinationScheme recom, double rparameter, fastjet::AreaDefinition *adef, JetAlg::MuonsStrategy usemuons=JetAlg::ALL_MUONS, JetAlg::InvisiblesStrategy useinvis=JetAlg::NO_INVISIBLES)
 FastJets (const FinalState &fsp, fastjet::JetDefinition::Plugin *plugin, JetAlg::MuonsStrategy usemuons=JetAlg::ALL_MUONS, JetAlg::InvisiblesStrategy useinvis=JetAlg::NO_INVISIBLES, fastjet::AreaDefinition *adef=nullptr)
 Explicitly pass in an externally-constructed plugin.
 FastJets (const FinalState &fsp, fastjet::JetDefinition::Plugin *plugin, fastjet::AreaDefinition *adef, JetAlg::MuonsStrategy usemuons=JetAlg::ALL_MUONS, JetAlg::InvisiblesStrategy useinvis=JetAlg::NO_INVISIBLES)
 Explicitly pass in an externally-constructed plugin, with reordered args for easier specification of jet area definition.
 FastJets (const FinalState &fsp, JetAlgName alg, double rparameter, JetAlg::MuonsStrategy usemuons=JetAlg::ALL_MUONS, JetAlg::InvisiblesStrategy useinvis=JetAlg::NO_INVISIBLES, fastjet::AreaDefinition *adef=nullptr, double seed_threshold=1.0)
 Convenience constructor using Rivet enums for most common jet algs (including some plugins).
 DEFAULT_RIVET_PROJ_CLONE (FastJets)
 Clone on the heap.
Access to the jets
Jets _jets () const
 Get the jets (unordered) with pT > ptmin.
PseudoJets pseudoJets (double ptmin=0.0) const
 Get the pseudo jets (unordered).
PseudoJets pseudojets (double ptmin=0.0) const
 Alias.
PseudoJets pseudoJetsByPt (double ptmin=0.0) const
 Get the pseudo jets, ordered by $ p_T $.
PseudoJets pseudojetsByPt (double ptmin=0.0) const
 Alias.
PseudoJets pseudoJetsByE (double ptmin=0.0) const
 Get the pseudo jets, ordered by $ E $.
PseudoJets pseudojetsByE (double ptmin=0.0) const
 Alias.
PseudoJets pseudoJetsByRapidity (double ptmin=0.0) const
 Get the pseudo jets, ordered by rapidity.
PseudoJets pseudojetsByRapidity (double ptmin=0.0) const
 Alias.
Jet trimJet (const Jet &input, const fastjet::Filter &trimmer) const
 Trim (filter) a jet, keeping tag and constituent info in the resulting jet.
Access to the FastJet clustering objects such as jet def, area def, and cluster
const shared_ptr
< fastjet::ClusterSequence > 
clusterSeq () const
const shared_ptr
< fastjet::ClusterSequenceArea > 
clusterSeqArea () const
const fastjet::JetDefinition & jetDef () const
 Return the jet definition.
const shared_ptr
< fastjet::AreaDefinition > 
areaDef () const
 Return the area definition.
Control the treatment of muons and invisible particles

Since MC-based jet calibration (and/or particle flow) can add back in particles that weren't seen in calorimeters/trackers.

void useMuons (MuonsStrategy usemuons=ALL_MUONS)
 Include (some) muons in jet construction.
void useInvisibles (InvisiblesStrategy useinvis=DECAY_INVISIBLES)
 Include (some) invisible particles in jet construction.
void useInvisibles (bool useinvis)
 Include (some) invisible particles in jet construction.
Access to jet objects
virtual Jets jets (const Cut &c=Cuts::open()) const
template<typename F >
Jets jets (F sorter, const Cut &c=Cuts::open()) const
template<typename F >
Jets jets (const Cut &c, F sorter) const
Jets jetsByPt (const Cut &c=Cuts::open()) const
Old jet accessors
Deprecated:
Use the versions with Cut arguments
Jets jets (double ptmin, double ptmax=MAXDOUBLE, double rapmin=-MAXDOUBLE, double rapmax=MAXDOUBLE, RapScheme rapscheme=PSEUDORAPIDITY) const
Jets jetsByPt (double ptmin) const
Old sorted jet accessors
Deprecated:
Use the versions with sorter function arguments. These will be removed in Rivet v3
Jets jetsByP (const Cut &c=Cuts::open()) const
Jets jetsByE (const Cut &c=Cuts::open()) const
Jets jetsByEt (const Cut &c=Cuts::open()) const
Projection "getting" functions
std::set< ConstProjectionPtrgetProjections () const
 Get the contained projections, including recursion.
template<typename PROJ >
const PROJ & getProjection (const std::string &name) const
const ProjectiongetProjection (const std::string &name) const
template<typename PROJ >
const PROJ & get (const std::string &name) const
Projection applying functions
template<typename PROJ >
const PROJ & applyProjection (const Event &evt, const Projection &proj) const
 Apply the supplied projection on event evt.
template<typename PROJ >
const PROJ & applyProjection (const Event &evt, const PROJ &proj) const
 Apply the supplied projection on event evt.
template<typename PROJ >
const PROJ & applyProjection (const Event &evt, const std::string &name) const
template<typename PROJ >
const PROJ & apply (const Event &evt, const Projection &proj) const
template<typename PROJ >
const PROJ & apply (const Event &evt, const PROJ &proj) const
template<typename PROJ >
const PROJ & apply (const Event &evt, const std::string &name) const

Protected Member Functions

void project (const Event &e)
 Perform the projection on the Event.
int compare (const Projection &p) const
 Compare projections.
Cmp< ProjectionmkNamedPCmp (const Projection &otherparent, const std::string &pname) const
Cmp< ProjectionmkPCmp (const Projection &otherparent, const std::string &pname) const
ProjectionHandlergetProjHandler () const
 Get a reference to the ProjectionHandler for this thread.
const Projection_applyProjection (const Event &evt, const std::string &name) const
const Projection_applyProjection (const Event &evt, const Projection &proj) const
Projection registration functions
template<typename PROJ >
const PROJ & declareProjection (const PROJ &proj, const std::string &name)
 Register a contained projection.
template<typename PROJ >
const PROJ & declare (const PROJ &proj, const std::string &name)
 Register a contained projection (user-facing version)
template<typename PROJ >
const PROJ & addProjection (const PROJ &proj, const std::string &name)
 Register a contained projection (user-facing version)
const Projection_declareProjection (const Projection &proj, const std::string &name)
 Untemplated function to do the work...

Protected Attributes

MuonsStrategy _useMuons
 Flag to determine whether or not to exclude (some) muons from the would-be constituents.
InvisiblesStrategy _useInvisibles
 Flag to determine whether or not to exclude (some) invisible particles from the would-be constituents.
bool _allowProjReg
 Flag to forbid projection registration in analyses until the init phase.

Private Member Functions

void _initBase ()
void _initJdef (JetAlgName alg, double rparameter, double seed_threshold)
Jet _mkJet (const PseudoJet &pj) const
 Function to make Rivet::Jet from fastjet::PseudoJet, including constituent and tag info.

Private Attributes

fastjet::JetDefinition _jdef
 Jet definition.
std::shared_ptr
< fastjet::AreaDefinition > 
_adef
 Pointer to user-handled area definition.
std::shared_ptr
< fastjet::ClusterSequence > 
_cseq
 Cluster sequence.
std::shared_ptr
< fastjet::JetDefinition::Plugin > 
_plugin
 FastJet external plugin.
std::map< int, vector< double > > _yscales
 Map of vectors of y scales. This is mutable so we can use caching/lazy evaluation.
std::map< int, Particle_particles
 set of particles sorted by their PT2

Friends

class Event
 Event is a friend.
class Cmp< Projection >
 The Cmp specialization for Projection is a friend.

Detailed Description

Project out jets found using the FastJet package jet algorithms.

Definition at line 26 of file FastJets.hh.


Member Typedef Documentation

typedef Jets collection_type [inherited]

Definition at line 218 of file JetAlg.hh.

typedef Jet entity_type [inherited]

Definition at line 217 of file JetAlg.hh.


Member Enumeration Documentation

enum InvisiblesStrategy [inherited]

Enum for the treatment of invisible particles: whether to include all, some, or none in jet-finding.

Enumerator:
NO_INVISIBLES 
DECAY_INVISIBLES 
ALL_INVISIBLES 

Definition at line 22 of file JetAlg.hh.

enum JetAlgName

Wrapper enum for selected FastJet jet algorithms.

Todo:
Move to JetAlg and alias here?
Enumerator:
KT 
CAM 
SISCONE 
ANTIKT 
ATLASCONE 
CMSCONE 
CDFJETCLU 
CDFMIDPOINT 
D0ILCONE 
JADE 
DURHAM 
TRACKJET 

Definition at line 31 of file FastJets.hh.

enum MuonsStrategy [inherited]

Enum for the treatment of muons: whether to include all, some, or none in jet-finding.

Enumerator:
NO_MUONS 
DECAY_MUONS 
ALL_MUONS 

Definition at line 19 of file JetAlg.hh.


Constructor & Destructor Documentation

FastJets ( const FinalState fsp,
const fastjet::JetDefinition &  jdef,
JetAlg::MuonsStrategy  usemuons = JetAlg::ALL_MUONS,
JetAlg::InvisiblesStrategy  useinvis = JetAlg::NO_INVISIBLES,
fastjet::AreaDefinition *  adef = nullptr 
) [inline]

Constructor from a FastJet JetDefinition

Warning:
The AreaDefinition pointer must be heap-allocated: it will be stored/deleted via a shared_ptr.

Definition at line 44 of file FastJets.hh.

      : JetAlg(fsp, usemuons, useinvis), _jdef(jdef), _adef(adef)
    {
      _initBase();
    }
FastJets ( const FinalState fsp,
const fastjet::JetDefinition &  jdef,
fastjet::AreaDefinition *  adef,
JetAlg::MuonsStrategy  usemuons = JetAlg::ALL_MUONS,
JetAlg::InvisiblesStrategy  useinvis = JetAlg::NO_INVISIBLES 
) [inline]

JetDefinition-based constructor with reordered args for easier specification of jet area definition

Warning:
The AreaDefinition pointer must be heap-allocated: it will be stored/deleted via a shared_ptr.

Definition at line 57 of file FastJets.hh.

      : FastJets(fsp, jdef, usemuons, useinvis, adef)
    {    }
FastJets ( const FinalState fsp,
fastjet::JetAlgorithm  type,
fastjet::RecombinationScheme  recom,
double  rparameter,
JetAlg::MuonsStrategy  usemuons = JetAlg::ALL_MUONS,
JetAlg::InvisiblesStrategy  useinvis = JetAlg::NO_INVISIBLES,
fastjet::AreaDefinition *  adef = nullptr 
) [inline]

Native argument constructor, using FastJet alg/scheme enums.

Warning:
The AreaDefinition pointer must be heap-allocated: it will be stored/deleted via a shared_ptr.

Definition at line 68 of file FastJets.hh.

      : FastJets(fsp, fastjet::JetDefinition(type, rparameter, recom), usemuons, useinvis, adef)
    {    }
FastJets ( const FinalState fsp,
fastjet::JetAlgorithm  type,
fastjet::RecombinationScheme  recom,
double  rparameter,
fastjet::AreaDefinition *  adef,
JetAlg::MuonsStrategy  usemuons = JetAlg::ALL_MUONS,
JetAlg::InvisiblesStrategy  useinvis = JetAlg::NO_INVISIBLES 
) [inline]

Native argument constructor with reordered args for easier specification of jet area definition

Warning:
The AreaDefinition pointer must be heap-allocated: it will be stored/deleted via a shared_ptr.

Definition at line 80 of file FastJets.hh.

      : FastJets(fsp, type, recom, rparameter, usemuons, useinvis, adef)
    {    }
FastJets ( const FinalState fsp,
fastjet::JetDefinition::Plugin *  plugin,
JetAlg::MuonsStrategy  usemuons = JetAlg::ALL_MUONS,
JetAlg::InvisiblesStrategy  useinvis = JetAlg::NO_INVISIBLES,
fastjet::AreaDefinition *  adef = nullptr 
) [inline]

Explicitly pass in an externally-constructed plugin.

Warning:
Provided plugin and area definition pointers must be heap-allocated; Rivet will store/delete via a shared_ptr

Definition at line 92 of file FastJets.hh.

      : FastJets(fsp, fastjet::JetDefinition(plugin), usemuons, useinvis, adef)
    {
      _plugin.reset(plugin);
    }
FastJets ( const FinalState fsp,
fastjet::JetDefinition::Plugin *  plugin,
fastjet::AreaDefinition *  adef,
JetAlg::MuonsStrategy  usemuons = JetAlg::ALL_MUONS,
JetAlg::InvisiblesStrategy  useinvis = JetAlg::NO_INVISIBLES 
) [inline]

Explicitly pass in an externally-constructed plugin, with reordered args for easier specification of jet area definition.

Warning:
Provided plugin and area definition pointers must be heap-allocated; Rivet will store/delete via a shared_ptr

Definition at line 105 of file FastJets.hh.

      : FastJets(fsp, plugin, usemuons, useinvis, adef)
    {    }
FastJets ( const FinalState fsp,
JetAlgName  alg,
double  rparameter,
JetAlg::MuonsStrategy  usemuons = JetAlg::ALL_MUONS,
JetAlg::InvisiblesStrategy  useinvis = JetAlg::NO_INVISIBLES,
fastjet::AreaDefinition *  adef = nullptr,
double  seed_threshold = 1.0 
) [inline]

Convenience constructor using Rivet enums for most common jet algs (including some plugins).

For the built-in algs, E-scheme recombination is used. For full control of FastJet built-in jet algs, use the constructors from native-args or a plugin pointer.

Warning:
Provided area definition pointer must be heap-allocated; Rivet will store/delete via a shared_ptr

Definition at line 120 of file FastJets.hh.

      : JetAlg(fsp, usemuons, useinvis)
    {
      _initBase();
      _initJdef(alg, rparameter, seed_threshold);
    }

Member Function Documentation

const Projection & _applyProjection ( const Event evt,
const std::string &  name 
) const [protected, inherited]

Non-templated version of string-based applyProjection, to work around header dependency issue.

Definition at line 22 of file ProjectionApplier.cc.

                                                                                  {
    return evt.applyProjection(getProjection(name));
  }
const Projection & _applyProjection ( const Event evt,
const Projection proj 
) const [protected, inherited]

Non-templated version of proj-based applyProjection, to work around header dependency issue.

Definition at line 28 of file ProjectionApplier.cc.

                                                                                      {
    return evt.applyProjection(proj);
  }
const Projection & _declareProjection ( const Projection proj,
const std::string &  name 
) [protected, inherited]

Untemplated function to do the work...

Definition at line 34 of file ProjectionApplier.cc.

                                                                             {
    if (!_allowProjReg) {
      cerr << "Trying to register projection '"
           << proj.name() << "' before init phase in '" << this->name() << "'." << endl;
      exit(2);
    }
    const Projection& reg = getProjHandler().registerProjection(*this, proj, name);
    return reg;
  }
void _initBase ( ) [private]

Shared utility functions to implement constructor behaviour

Todo:
Replace with calls between constructors when C++11 available?

Definition at line 11 of file FastJets.cc.

                           {
    setName("FastJets");
    addProjection(HeavyHadrons(), "HFHadrons");
    addProjection(TauFinder(TauFinder::HADRONIC), "Taus");
  }
void _initJdef ( JetAlgName  alg,
double  rparameter,
double  seed_threshold 
) [private]

Definition at line 22 of file FastJets.cc.

                                                                                   {
    MSG_DEBUG("JetAlg = " << alg);
    MSG_DEBUG("R parameter = " << rparameter);
    MSG_DEBUG("Seed threshold = " << seed_threshold);
    if (alg == KT) {
      _jdef = fastjet::JetDefinition(fastjet::kt_algorithm, rparameter, fastjet::E_scheme);
    } else if (alg == CAM) {
      _jdef = fastjet::JetDefinition(fastjet::cambridge_algorithm, rparameter, fastjet::E_scheme);
    } else if (alg == ANTIKT) {
      _jdef = fastjet::JetDefinition(fastjet::antikt_algorithm, rparameter, fastjet::E_scheme);
    } else if (alg == DURHAM) {
      _jdef = fastjet::JetDefinition(fastjet::ee_kt_algorithm, fastjet::E_scheme);
    } else {
      // Plugins:
      if (alg == SISCONE) {
        const double OVERLAP_THRESHOLD = 0.75;
        _plugin.reset(new fastjet::SISConePlugin(rparameter, OVERLAP_THRESHOLD));
      // } else if (alg == PXCONE) {
      //   string msg = "PxCone currently not supported, since FastJet doesn't install it by default. ";
      //   msg += "Please notify the Rivet authors if this behaviour should be changed.";
      //   throw Error(msg);
      //  _plugin.reset(new fastjet::PxConePlugin(rparameter));
      } else if (alg == ATLASCONE) {
        const double OVERLAP_THRESHOLD = 0.5;
        _plugin.reset(new fastjet::ATLASConePlugin(rparameter, seed_threshold, OVERLAP_THRESHOLD));
      } else if (alg == CMSCONE) {
        _plugin.reset(new fastjet::CMSIterativeConePlugin(rparameter, seed_threshold));
      } else if (alg == CDFJETCLU) {
        const double OVERLAP_THRESHOLD = 0.75;
        _plugin.reset(new fastjet::CDFJetCluPlugin(rparameter, OVERLAP_THRESHOLD, seed_threshold));
      } else if (alg == CDFMIDPOINT) {
        const double OVERLAP_THRESHOLD = 0.5;
        _plugin.reset(new fastjet::CDFMidPointPlugin(rparameter, OVERLAP_THRESHOLD, seed_threshold));
      } else if (alg == D0ILCONE) {
        const double min_jet_Et = 6.0;
        _plugin.reset(new fastjet::D0RunIIConePlugin(rparameter, min_jet_Et));
      } else if (alg == JADE) {
        _plugin.reset(new fastjet::JadePlugin());
      } else if (alg == TRACKJET) {
        _plugin.reset(new fastjet::TrackJetPlugin(rparameter));
      }
      _jdef = fastjet::JetDefinition(_plugin.get());
    }
  }
Jets _jets ( ) const [virtual]

Get the jets (unordered) with pT > ptmin.

Todo:
Cache?

Implements JetAlg.

Definition at line 183 of file FastJets.cc.

                             {
    Jets rtn; rtn.reserve(pseudojets().size());
    foreach (const fastjet::PseudoJet& pj, pseudojets()) {
      rtn.push_back(_mkJet(pj));
    }
    /// @todo Cache?
    return rtn;
  }
Jet _mkJet ( const PseudoJet &  pj) const [private]

Function to make Rivet::Jet from fastjet::PseudoJet, including constituent and tag info.

Definition at line 205 of file FastJets.cc.

                                              {
    assert(clusterSeq());

    // Take the constituents from the cluster sequence, unless the jet was not
    // associated with the cluster sequence (could be the case for trimmed jets)
    const PseudoJets parts = (pj.associated_cluster_sequence() == clusterSeq().get())
      ? clusterSeq()->constituents(pj) : pj.constituents();

    vector<Particle> constituents, tags;
    constituents.reserve(parts.size());
    for (const fastjet::PseudoJet& p : parts) {
      map<int, Particle>::const_iterator found = _particles.find(p.user_index());
      // assert(found != _particles.end());
      if (found == _particles.end() && p.is_pure_ghost()) continue; //< Pure FJ ghosts are ok
      assert(found != _particles.end()); //< Anything else must be known
      assert(found->first != 0); //< All mapping IDs are pos-def (particles) or neg-def (tags)
      if (found->first > 0) constituents.push_back(found->second);
      else if (found->first < 0) tags.push_back(found->second);
    }

    return Jet(pj, constituents, tags);
  }
Projection& addPdgIdPair ( PdgId  beam1,
PdgId  beam2 
) [inline, inherited]

Add a colliding beam pair.

Definition at line 108 of file Projection.hh.

                                                       {
      _beamPairs.insert(PdgIdPair(beam1, beam2));
      return *this;
    }
const PROJ& addProjection ( const PROJ &  proj,
const std::string &  name 
) [inline, protected, inherited]

Register a contained projection (user-facing version)

Deprecated:
Use declareProjection() or declare()
Todo:
Add SFINAE to require that PROJ inherit from Projection

Definition at line 157 of file ProjectionApplier.hh.

{ return declareProjection(proj, name); }
const PROJ& apply ( const Event evt,
const Projection proj 
) const [inline, inherited]

Apply the supplied projection on event evt (user-facing alias).

Todo:
Add SFINAE to require that PROJ inherit from Projection

Definition at line 80 of file ProjectionApplier.hh.

{ return applyProjection<PROJ>(evt, proj); }
const PROJ& apply ( const Event evt,
const PROJ &  proj 
) const [inline, inherited]

Apply the supplied projection on event evt (user-facing alias).

Todo:
Add SFINAE to require that PROJ inherit from Projection

Definition at line 92 of file ProjectionApplier.hh.

{ return applyProjection<PROJ>(evt, proj); }
const PROJ& apply ( const Event evt,
const std::string &  name 
) const [inline, inherited]

Apply the supplied projection on event evt (user-facing alias).

Todo:
Add SFINAE to require that PROJ inherit from Projection

Definition at line 104 of file ProjectionApplier.hh.

{ return applyProjection<PROJ>(evt, name); }
const PROJ& applyProjection ( const Event evt,
const Projection proj 
) const [inline, inherited]

Apply the supplied projection on event evt.

Todo:
Add SFINAE to require that PROJ inherit from Projection

Definition at line 74 of file ProjectionApplier.hh.

                                                                                {
      return pcast<PROJ>(_applyProjection(evt, proj));
    }
const PROJ& applyProjection ( const Event evt,
const PROJ &  proj 
) const [inline, inherited]

Apply the supplied projection on event evt.

Todo:
Add SFINAE to require that PROJ inherit from Projection

Definition at line 86 of file ProjectionApplier.hh.

                                                                          {
      return pcast<PROJ>(_applyProjection(evt, proj));
    }
const PROJ& applyProjection ( const Event evt,
const std::string &  name 
) const [inline, inherited]

Apply the named projection on event evt.

Todo:
Add SFINAE to require that PROJ inherit from Projection

Definition at line 98 of file ProjectionApplier.hh.

                                                                               {
      return pcast<PROJ>(_applyProjection(evt, name));
    }
const shared_ptr<fastjet::AreaDefinition> areaDef ( ) const [inline]

Return the area definition.

Warning:
May be null!
Todo:
Care needed re. const shared_ptr<T> vs. shared_ptr<const T>

Definition at line 234 of file FastJets.hh.

                                                            {
      return _adef;
    }
const set< PdgIdPair > beamPairs ( ) const [virtual, inherited]

Return the allowed beam pairs on which this projection can operate, not including recursion. Derived classes should ensure that all contained projections are registered in the _projections set for the beam constraint chaining to work.

Todo:
Remove the beam constraints system from projections.

Definition at line 35 of file Projection.cc.

                                                   {
    set<PdgIdPair> ret = _beamPairs;
    set<ConstProjectionPtr> projs = getProjections();
    for (set<ConstProjectionPtr>::const_iterator ip = projs.begin(); ip != projs.end(); ++ip) {
      ConstProjectionPtr p = *ip;
      getLog() << Log::TRACE << "Proj addr = " << p << endl;
      if (p) ret = intersection(ret, p->beamPairs());
    }
    return ret;
  }
bool before ( const Projection p) const [inherited]

Determine whether this object should be ordered before the object p given as argument. If p is of a different class than this, the before() function of the corresponding type_info objects is used. Otherwise, if the objects are of the same class, the virtual compare(const Projection &) will be returned.

Definition at line 24 of file Projection.cc.

                                                   {
    const std::type_info& thisid = typeid(*this);
    const std::type_info& otherid = typeid(p);
    if (thisid == otherid) {
      return compare(p) < 0;
    } else {
      return thisid.before(otherid);
    }
  }
void calc ( const Particles fsparticles,
const Particles tagparticles = Particles() 
)

Do the calculation locally (no caching).

Todo:
Use FastJet3's UserInfo system

< Ghostify the momentum

Definition at line 136 of file FastJets.cc.

                                                                                 {
    _particles.clear();
    vector<fastjet::PseudoJet> pjs;

    MSG_DEBUG("Finding jets from " << fsparticles.size() << " input particles");

    /// @todo Use FastJet3's UserInfo system

    // Store 4 vector data about each particle into FastJet's PseudoJets
    int counter = 1;
    for (const Particle& p : fsparticles) {
      fastjet::PseudoJet pj = p;
      pj.set_user_index(counter);
      pjs.push_back(pj);
      _particles[counter] = p;
      counter += 1;
    }
    // And the same for ghost tagging particles (with negative user indices)
    counter = 1;
    for (const Particle& p : tagparticles) {
      fastjet::PseudoJet pj = p;
      pj *= 1e-20; ///< Ghostify the momentum
      pj.set_user_index(-counter);
      pjs.push_back(pj);
      _particles[-counter] = p;
      counter += 1;
    }

    // Choose cseq as basic or area-calculating
    if (_adef) {
      _cseq.reset(new fastjet::ClusterSequenceArea(pjs, _jdef, *_adef));
    } else {
      _cseq.reset(new fastjet::ClusterSequence(pjs, _jdef));
    }
    MSG_DEBUG("FastJet ClusterSequence constructed; Njets_tot = "
              << _cseq->inclusive_jets().size() << ", Njets_10 = "
              << _cseq->inclusive_jets(10*GeV).size()); //< only inefficient in debug mode
  }
virtual unique_ptr<Projection> clone ( ) const [pure virtual, inherited]

Clone on the heap.

Implements Projection.

const shared_ptr<fastjet::ClusterSequence> clusterSeq ( ) const [inline]

Return the cluster sequence.

Todo:
Care needed re. const shared_ptr<T> vs. shared_ptr<const T>

Definition at line 209 of file FastJets.hh.

                                                                {
      return _cseq;
    }
const shared_ptr<fastjet::ClusterSequenceArea> clusterSeqArea ( ) const [inline]

Return the area-enabled cluster sequence (if an area defn exists, otherwise returns a null ptr).

Todo:
Care needed re. const shared_ptr<T> vs. shared_ptr<const T>

Definition at line 218 of file FastJets.hh.

                                                                        {
      return areaDef() ? dynamic_pointer_cast<fastjet::ClusterSequenceArea>(_cseq) : nullptr;
    }
int compare ( const Projection p) const [protected, virtual]

Compare projections.

Implements JetAlg.

Definition at line 92 of file FastJets.cc.

                                                 {
    const FastJets& other = dynamic_cast<const FastJets&>(p);
    return \
      cmp(_useMuons, other._useMuons) ||
      cmp(_useInvisibles, other._useInvisibles) ||
      mkNamedPCmp(other, "FS") ||
      cmp(_jdef.jet_algorithm(), other._jdef.jet_algorithm()) ||
      cmp(_jdef.recombination_scheme(), other._jdef.recombination_scheme()) ||
      cmp(_jdef.plugin(), other._jdef.plugin()) ||
      cmp(_jdef.R(), other._jdef.R()) ||
      cmp(_adef, other._adef);
  }
const PROJ& declare ( const PROJ &  proj,
const std::string &  name 
) [inline, protected, inherited]

Register a contained projection (user-facing version)

Todo:
Add SFINAE to require that PROJ inherit from Projection

Definition at line 151 of file ProjectionApplier.hh.

{ return declareProjection(proj, name); }
const PROJ& declareProjection ( const PROJ &  proj,
const std::string &  name 
) [inline, protected, inherited]

Register a contained projection.

The type of the argument is used to instantiate a new projection internally: this new object is applied to events rather than the argument object. Hence you are advised to only use locally-scoped Projection objects in your Projection and Analysis constructors, and to avoid polymorphism (e.g. handling ConcreteProjection via a pointer or reference to type Projection) since this will screw up the internal type management.

Todo:
Add SFINAE to require that PROJ inherit from Projection

Definition at line 142 of file ProjectionApplier.hh.

                                                                           {
      const Projection& reg = _declareProjection(proj, name);
      const PROJ& rtn = dynamic_cast<const PROJ&>(reg);
      return rtn;
    }

Clone on the heap.

bool empty ( ) const [inline, inherited]

Whether the inclusive jet collection is empty.

Definition at line 212 of file JetAlg.hh.

{ return size() != 0; }
collection_type entities ( ) const [inline, inherited]

Template-usable interface common to FinalState.

Definition at line 221 of file JetAlg.hh.

{ return jets(); }
const PROJ& get ( const std::string &  name) const [inline, inherited]

Get the named projection, specifying return type via a template argument (user-facing alias).

Todo:
Add SFINAE to require that PROJ inherit from Projection

Definition at line 57 of file ProjectionApplier.hh.

{ return getProjection<PROJ>(name); }
Log& getLog ( ) const [inline, inherited]

Get a Log object based on the getName() property of the calling projection object.

Reimplemented from ProjectionApplier.

Definition at line 115 of file Projection.hh.

                        {
      string logname = "Rivet.Projection." + name();
      return Log::getLog(logname);
    }
const PROJ& getProjection ( const std::string &  name) const [inline, inherited]

Get the named projection, specifying return type via a template argument.

Todo:
Add SFINAE to require that PROJ inherit from Projection

Definition at line 50 of file ProjectionApplier.hh.

                                                           {
      const Projection& p = getProjHandler().getProjection(*this, name);
      return pcast<PROJ>(p);
    }
const Projection& getProjection ( const std::string &  name) const [inline, inherited]

Get the named projection (non-templated, so returns as a reference to a Projection base class).

Definition at line 61 of file ProjectionApplier.hh.

                                                                 {
      return getProjHandler().getProjection(*this, name);
    }
std::set<ConstProjectionPtr> getProjections ( ) const [inline, inherited]

Get the contained projections, including recursion.

Definition at line 43 of file ProjectionApplier.hh.

ProjectionHandler& getProjHandler ( ) const [inline, protected, inherited]

Get a reference to the ProjectionHandler for this thread.

Definition at line 122 of file ProjectionApplier.hh.

                                              {
      return _projhandler;
    }
const fastjet::JetDefinition& jetDef ( ) const [inline]

Return the jet definition.

Definition at line 226 of file FastJets.hh.

                                               {
      return _jdef;
    }
virtual Jets jets ( const Cut c = Cuts::open()) const [inline, virtual, inherited]

Get jets in no guaranteed order, with optional cuts on $ p_\perp $ and rapidity.

Note:
Returns a copy rather than a reference, due to cuts

Definition at line 85 of file JetAlg.hh.

                                                     {
      return filterBy(_jets(), c);
      // const Jets rawjets = _jets();
      // // Just return a copy of rawjets if the cut is open
      // if (c == Cuts::open()) return rawjets;
      // // If there is a non-trivial cut...
      // /// @todo Use an STL erase(remove_if) and lambda function for this
      // Jets rtn;
      // rtn.reserve(size());
      // foreach (const Jet& j, rawjets)
      //   if (c->accept(j)) rtn.push_back(j);
      // return rtn;
    }
Jets jets ( sorter,
const Cut c = Cuts::open() 
) const [inline, inherited]
Todo:
Want to add a general filtering function, but that clashes with the sorting functor... SFINAE?

Get the jets, ordered by supplied sorting function object, with optional cuts on $ p_\perp $ and rapidity.

Note:
Returns a copy rather than a reference, due to cuts and sorting
Todo:
Will the vector be efficiently std::move'd by value through this function chain?

Definition at line 106 of file JetAlg.hh.

                                                       {
      /// @todo Will the vector be efficiently std::move'd by value through this function chain?
      return sortBy(jets(c), sorter);
    }
Jets jets ( const Cut c,
sorter 
) const [inline, inherited]

Get the jets, ordered by supplied sorting function object, with optional cuts on $ p_\perp $ and rapidity.

Note:
Returns a copy rather than a reference, due to cuts and sorting
Todo:
Will the vector be efficiently std::move'd by value through this function chain?

Definition at line 114 of file JetAlg.hh.

                                            {
      /// @todo Will the vector be efficiently std::move'd by value through this function chain?
      return sortBy(jets(c), sorter);
    }
Jets jets ( double  ptmin,
double  ptmax = MAXDOUBLE,
double  rapmin = -MAXDOUBLE,
double  rapmax = MAXDOUBLE,
RapScheme  rapscheme = PSEUDORAPIDITY 
) const [inline, inherited]

Get jets in no guaranteed order, with optional cuts on $ p_\perp $ and rapidity.

Deprecated:
Use the version with a Cut argument
Note:
Returns a copy rather than a reference, due to cuts

Definition at line 173 of file JetAlg.hh.

                                                        {
      if (rapscheme == PSEUDORAPIDITY) {
        return jets((Cuts::pT >= ptmin) & (Cuts::pT < ptmax) & (Cuts::rapIn(rapmin, rapmax)));
      } else if (rapscheme == RAPIDITY) {
        return jets((Cuts::pT >= ptmin) & (Cuts::pT < ptmax) & (Cuts::etaIn(rapmin, rapmax)));
      }
      throw LogicError("Unknown rapidity scheme. This shouldn't be possible!");
    }
Jets jetsByE ( const Cut c = Cuts::open()) const [inline, inherited]

Get the jets, ordered by $ E $, with optional cuts on $ p_\perp $ and rapidity.

Note:
Returns a copy rather than a reference, due to cuts and sorting
Deprecated:
Use the version with a sorter function argument.

Definition at line 149 of file JetAlg.hh.

                                                {
      return jets(c, cmpMomByE);
    }
Jets jetsByEt ( const Cut c = Cuts::open()) const [inline, inherited]

Get the jets, ordered by $ E_T $, with optional cuts on $ p_\perp $ and rapidity.

Note:
Returns a copy rather than a reference, due to cuts and sorting
Deprecated:
Use the version with a sorter function argument.

Definition at line 157 of file JetAlg.hh.

                                                 {
      return jets(c, cmpMomByEt);
    }
Jets jetsByP ( const Cut c = Cuts::open()) const [inline, inherited]

Get the jets, ordered by $ |p| $, with optional cuts on $ p_\perp $ and rapidity.

Note:
Returns a copy rather than a reference, due to cuts and sorting
Deprecated:
Use the version with a sorter function argument.

Definition at line 141 of file JetAlg.hh.

                                                {
      return jets(c, cmpMomByP);
    }
Jets jetsByPt ( const Cut c = Cuts::open()) const [inline, inherited]

Get the jets, ordered by $ p_T $, with optional cuts.

Note:
Returns a copy rather than a reference, due to cuts and sorting

This is a very common use-case, so is available as syntatic sugar for jets(c, cmpMomByPt).

Todo:
The other sorted accessors should be removed in a cleanup.

Definition at line 126 of file JetAlg.hh.

                                                 {
      return jets(c, cmpMomByPt);
    }
Jets jetsByPt ( double  ptmin) const [inline, inherited]

Get the jets, ordered by $ p_T $, with a cut on $ p_\perp $.

Deprecated:
Use the version with a Cut argument
Note:
Returns a copy rather than a reference, due to cuts and sorting

This is a very common use-case, so is available as syntatic sugar for jets(Cuts::pT >= ptmin, cmpMomByPt).

Todo:
The other sorted accessors should be removed in a cleanup.

Definition at line 191 of file JetAlg.hh.

                                      {
      return jets(Cuts::pT >= ptmin, cmpMomByPt);
    }
void markAsOwned ( ) const [inline, inherited]

Mark object as owned by the _projhandler

Todo:
Huh? What's this for?

Definition at line 111 of file ProjectionApplier.hh.

{ _owned = true; }
Cmp< Projection > mkNamedPCmp ( const Projection otherparent,
const std::string &  pname 
) const [protected, inherited]

Shortcut to make a named Cmp<Projection> comparison with the *this object automatically passed as one of the parent projections.

Definition at line 47 of file Projection.cc.

                                                                                                  {
    return pcmp(*this, otherparent, pname);
  }
Cmp< Projection > mkPCmp ( const Projection otherparent,
const std::string &  pname 
) const [protected, inherited]

Shortcut to make a named Cmp<Projection> comparison with the *this object automatically passed as one of the parent projections.

Note:
Alias for mkNamedPCmp

Definition at line 51 of file Projection.cc.

                                                                                             {
    return pcmp(*this, otherparent, pname);
  }
virtual std::string name ( ) const [inline, virtual, inherited]

Get the name of the projection.

Implements ProjectionApplier.

Definition at line 102 of file Projection.hh.

                                   {
      return _name;
    }
size_t numJets ( const Cut c = Cuts::open()) const [inline, inherited]

Number of jets passing the provided Cut.

Definition at line 207 of file JetAlg.hh.

{ return jets(c).size(); }
void project ( const Event e) [protected, virtual]

Perform the projection on the Event.

Todo:
Allow the user to specify tag particle kinematic thresholds

Implements JetAlg.

Definition at line 114 of file FastJets.cc.

                                       {
    // Assemble final state particles
    const string fskey = (_useInvisibles == JetAlg::NO_INVISIBLES) ? "VFS" : "FS";
    Particles fsparticles = applyProjection<FinalState>(e, fskey).particles();
    // Remove prompt invisibles if needed (already done by VFS if using NO_INVISIBLES)
    if (_useInvisibles == JetAlg::DECAY_INVISIBLES)
      fsparticles.erase( std::remove_if(fsparticles.begin(), fsparticles.end(), isPromptInvisible), fsparticles.end() );
    // Remove prompt/all muons if needed
    if (_useMuons == JetAlg::DECAY_MUONS)
      fsparticles.erase( std::remove_if(fsparticles.begin(), fsparticles.end(), isPromptMuon), fsparticles.end() );
    else if (_useMuons == JetAlg::NO_MUONS)
      fsparticles.erase( std::remove_if(fsparticles.begin(), fsparticles.end(), isMuon), fsparticles.end() );

    // Tagging particles
    /// @todo Allow the user to specify tag particle kinematic thresholds
    const Particles chadrons = applyProjection<HeavyHadrons>(e, "HFHadrons").cHadrons();
    const Particles bhadrons = applyProjection<HeavyHadrons>(e, "HFHadrons").bHadrons();
    const Particles taus = applyProjection<FinalState>(e, "Taus").particles();
    calc(fsparticles, chadrons+bhadrons+taus);
  }
PseudoJets pseudoJets ( double  ptmin = 0.0) const

Get the pseudo jets (unordered).

Definition at line 200 of file FastJets.cc.

                                                    {
    return clusterSeq() ? clusterSeq()->inclusive_jets(ptmin) : PseudoJets();
  }
PseudoJets pseudojets ( double  ptmin = 0.0) const [inline]

Alias.

Definition at line 175 of file FastJets.hh.

{ return pseudoJets(ptmin); }
PseudoJets pseudoJetsByE ( double  ptmin = 0.0) const [inline]

Get the pseudo jets, ordered by $ E $.

Definition at line 185 of file FastJets.hh.

                                                     {
      return sorted_by_E(pseudoJets(ptmin));
    }
PseudoJets pseudojetsByE ( double  ptmin = 0.0) const [inline]

Alias.

Definition at line 189 of file FastJets.hh.

{ return pseudoJetsByE(ptmin); }
PseudoJets pseudoJetsByPt ( double  ptmin = 0.0) const [inline]

Get the pseudo jets, ordered by $ p_T $.

Definition at line 178 of file FastJets.hh.

                                                      {
      return sorted_by_pt(pseudoJets(ptmin));
    }
PseudoJets pseudojetsByPt ( double  ptmin = 0.0) const [inline]

Alias.

Definition at line 182 of file FastJets.hh.

{ return pseudoJetsByPt(ptmin); }
PseudoJets pseudoJetsByRapidity ( double  ptmin = 0.0) const [inline]

Get the pseudo jets, ordered by rapidity.

Definition at line 192 of file FastJets.hh.

                                                            {
      return sorted_by_rapidity(pseudoJets(ptmin));
    }
PseudoJets pseudojetsByRapidity ( double  ptmin = 0.0) const [inline]

Alias.

Definition at line 196 of file FastJets.hh.

{ return pseudoJetsByRapidity(ptmin); }
void reset ( ) [virtual]

Reset the projection. Jet def, etc. are unchanged.

Todo:
_cseq = fastjet::ClusterSequence();

Implements JetAlg.

Definition at line 176 of file FastJets.cc.

                       {
    _yscales.clear();
    _particles.clear();
    /// @todo _cseq = fastjet::ClusterSequence();
  }
void setName ( const std::string &  name) [inline, inherited]

Used by derived classes to set their name.

Definition at line 121 of file Projection.hh.

                                        {
      _name = name;
    }
size_t size ( ) const [inline, inherited]

Number of jets (without cuts).

Definition at line 210 of file JetAlg.hh.

{ return jets().size(); }
Jet trimJet ( const Jet input,
const fastjet::Filter &  trimmer 
) const

Trim (filter) a jet, keeping tag and constituent info in the resulting jet.

Definition at line 193 of file FastJets.cc.

                                                                          {
    assert(input.pseudojet().associated_cluster_sequence() == clusterSeq().get());
    PseudoJet pj = trimmer(input);
    return _mkJet(pj);
  }
void useInvisibles ( InvisiblesStrategy  useinvis = DECAY_INVISIBLES) [inline, inherited]

Include (some) invisible particles in jet construction.

The default behaviour is that jets are only constructed from visible particles. Some jet studies, including those from ATLAS, use a definition in which neutrinos from hadron decays are included via MC-based calibrations. Setting this flag to true avoids the automatic restriction to a VisibleFinalState.

Definition at line 61 of file JetAlg.hh.

                                                                     {
      _useInvisibles = useinvis;
    }
void useInvisibles ( bool  useinvis) [inline, inherited]

Include (some) invisible particles in jet construction.

The default behaviour is that jets are only constructed from visible particles. Some jet studies, including those from ATLAS, use a definition in which neutrinos from hadron decays are included via MC-based calibrations. Setting this flag to true avoids the automatic restriction to a VisibleFinalState.

Deprecated:
Use the enum-arg version instead. Will be removed in Rivet v3

Definition at line 73 of file JetAlg.hh.

void useJetArea ( fastjet::AreaDefinition *  adef) [inline]

Use provided jet area definition.

Warning:
The provided pointer must be heap-allocated: it will be stored/deleted via a shared_ptr.
Note:
Provide an adef null pointer to re-disable jet area calculation

Definition at line 161 of file FastJets.hh.

                                                 {
      _adef.reset(adef);
    }
void useMuons ( MuonsStrategy  usemuons = ALL_MUONS) [inline, inherited]

Include (some) muons in jet construction.

The default behaviour is that jets are only constructed from visible particles. Some jet studies, including those from ATLAS, use a definition in which neutrinos from hadron decays are included via MC-based calibrations. Setting this flag to true avoids the automatic restriction to a VisibleFinalState.

Definition at line 51 of file JetAlg.hh.

                                                    {
      _useMuons = usemuons;
    }

Friends And Related Function Documentation

friend class Cmp< Projection > [friend, inherited]

The Cmp specialization for Projection is a friend.

Definition at line 36 of file Projection.hh.

friend class Event [friend, inherited]

Event is a friend.

Definition at line 33 of file Projection.hh.


Member Data Documentation

std::shared_ptr<fastjet::AreaDefinition> _adef [private]

Pointer to user-handled area definition.

Definition at line 277 of file FastJets.hh.

bool _allowProjReg [protected, inherited]

Flag to forbid projection registration in analyses until the init phase.

Definition at line 176 of file ProjectionApplier.hh.

std::shared_ptr<fastjet::ClusterSequence> _cseq [private]

Cluster sequence.

Definition at line 280 of file FastJets.hh.

fastjet::JetDefinition _jdef [private]

Jet definition.

Definition at line 274 of file FastJets.hh.

std::map<int, Particle> _particles [private]

set of particles sorted by their PT2

Definition at line 290 of file FastJets.hh.

std::shared_ptr<fastjet::JetDefinition::Plugin> _plugin [private]

FastJet external plugin.

Definition at line 283 of file FastJets.hh.

InvisiblesStrategy _useInvisibles [protected, inherited]

Flag to determine whether or not to exclude (some) invisible particles from the would-be constituents.

Definition at line 242 of file JetAlg.hh.

MuonsStrategy _useMuons [protected, inherited]

Flag to determine whether or not to exclude (some) muons from the would-be constituents.

Definition at line 239 of file JetAlg.hh.

std::map<int, vector<double> > _yscales [mutable, private]

Map of vectors of y scales. This is mutable so we can use caching/lazy evaluation.

Definition at line 286 of file FastJets.hh.


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