rivet is hosted by Hepforge, IPPP Durham
Rivet 4.0.0

Find jets using jet algorithms via the FastJet package. More...

#include <FastJets.hh>

Inheritance diagram for Rivet::FastJets:
Rivet::JetFinder Rivet::Projection Rivet::ProjectionApplier

Public Types

typedef Jet entity_type
 
typedef Jets collection_type
 

Public Member Functions

void reset ()
 Reset the projection. Jet def, etc. are unchanged.
 
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 size () const
 Count the jets.
 
size_t size (const Cut &c) const
 Count the jets after a Cut is applied.
 
size_t size (const JetSelector &s) const
 Count the jets after a selection functor is applied.
 
bool empty () const
 Is this jet finder empty?
 
bool empty (const Cut &c) const
 Is this jet finder empty after a Cut is applied?
 
bool empty (const JetSelector &s) const
 Is this jet finder empty after a selection functor is applied?
 
collection_type entities () const
 Template-usable interface common to FinalState.
 
virtual std::string name () const
 Get the name of the projection.
 
bool valid () const
 Get the state of the projetion.
 
bool failed () const
 Get the state of the projetion.
 
void markAsOwned () const
 Mark this object as owned by a proj-handler.
 
Constructors etc.
 FastJets (const FinalState &fsp, const fastjet::JetDefinition &jdef, JetMuons usemuons=JetMuons::ALL, JetInvisibles useinvis=JetInvisibles::NONE, fastjet::AreaDefinition *adef=nullptr)
 
 FastJets (const FinalState &fsp, const fastjet::JetDefinition &jdef, fastjet::AreaDefinition *adef, JetMuons usemuons=JetMuons::ALL, JetInvisibles useinvis=JetInvisibles::NONE)
 
 FastJets (const FinalState &fsp, fastjet::JetAlgorithm type, fastjet::RecombinationScheme recom, double rparameter, JetMuons usemuons=JetMuons::ALL, JetInvisibles useinvis=JetInvisibles::NONE, fastjet::AreaDefinition *adef=nullptr)
 
 FastJets (const FinalState &fsp, fastjet::JetAlgorithm type, fastjet::RecombinationScheme recom, double rparameter, fastjet::AreaDefinition *adef, JetMuons usemuons=JetMuons::ALL, JetInvisibles useinvis=JetInvisibles::NONE)
 
 FastJets (const FinalState &fsp, fastjet::JetDefinition::Plugin *plugin, JetMuons usemuons=JetMuons::ALL, JetInvisibles useinvis=JetInvisibles::NONE, fastjet::AreaDefinition *adef=nullptr)
 Explicitly pass in an externally-constructed plugin.
 
 FastJets (const FinalState &fsp, fastjet::JetDefinition::Plugin *plugin, fastjet::AreaDefinition *adef, JetMuons usemuons=JetMuons::ALL, JetInvisibles useinvis=JetInvisibles::NONE)
 Explicitly pass in an externally-constructed plugin, with reordered args for easier specification of jet area definition.
 
 FastJets (const FinalState &fsp, JetAlg alg, double rparameter, JetMuons usemuons=JetMuons::ALL, JetInvisibles useinvis=JetInvisibles::NONE, fastjet::AreaDefinition *adef=nullptr, double seed_threshold=1.0)
 Convenience constructor using Rivet enums for most common jet algs (including some plugins).
 
 RIVET_DEFAULT_PROJ_CLONE (FastJets)
 Clone on the heap.
 
Jet-area calculations
void useJetArea (fastjet::AreaDefinition *adef)
 Use provided jet area definition.
 
void clearJetArea ()
 Don't calculate a jet area.
 
Jet grooming
void addTrf (fastjet::Transformer *trf)
 Add a grooming transformer (base class of fastjet::Filter, etc.)
 
template<typename TRFS , typename TRF = typename TRFS::value_type>
std::enable_if< Derefable< TRF >::value, void >::type addTrfs (const TRFS &trfs)
 Add a list of grooming transformers.
 
void clearTrfs ()
 Don't apply any jet transformers.
 
Access to the jets
PseudoJets pseudojets (double ptmin=0.0) const
 
PseudoJets pseudojetsByPt (double ptmin=0.0) const
 Get the pseudo jets, ordered by \( p_T \).
 
PseudoJets pseudojetsByE (double ptmin=0.0) const
 Get the pseudo jets, ordered by \( E \).
 
PseudoJets pseudojetsByRapidity (double ptmin=0.0) const
 Get the pseudo jets, ordered by rapidity.
 
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 (JetMuons usemuons=JetMuons::ALL)
 Include (some) muons in jet construction.
 
void useInvisibles (JetInvisibles useinvis=JetInvisibles::DECAY)
 Include (some) invisible particles in jet construction.
 
Access to jet objects
virtual Jets jets (const Cut &c=Cuts::open()) const
 
virtual Jets jets (const JetSelector &selector) const
 
Jets jets (const Cut &c, const JetSorter &sorter) const
 
Jets jets (const JetSorter &sorter, const Cut &c=Cuts::open()) const
 
Jets jets (const JetSelector &selector, const JetSorter &sorter) const
 
Jets jets (const JetSorter &sorter, const JetSelector selector) const
 
Jets jetsByPt (const Cut &c=Cuts::open()) const
 
Jets jetsByPt (const JetSelector &selector) const
 
Projection operation and comparison
bool before (const Projection &p) const
 
Projection "getting" functions
std::set< ConstProjectionPtr > getProjections () const
 Get the contained projections, including recursion.
 
std::set< ConstProjectionPtr > getImmediateChildProjections () const
 Get the contained projections, excluding recursion.
 
bool hasProjection (const std::string &name) const
 Does this applier have a projection registered under the name name?
 
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
 
template<typename PROJ >
const PROJ & getProjectionFromDeclQueue (const std::string name) const
 
Projection applying functions
template<typename PROJ = Projection>
std::enable_if_t< std::is_base_of< Projection, PROJ >::value, const PROJ & > apply (const Event &evt, const Projection &proj) const
 Apply the supplied projection on event evt.
 
template<typename PROJ = Projection>
std::enable_if_t< std::is_base_of< Projection, PROJ >::value, const PROJ & > apply (const Event &evt, const PROJ &proj) const
 Apply the supplied projection on event evt (user-facing alias).
 
template<typename PROJ = Projection>
std::enable_if_t< std::is_base_of< Projection, PROJ >::value, const PROJ & > apply (const Event &evt, const std::string &name) const
 Apply the supplied projection on event evt (user-facing alias).
 
template<typename PROJ = Projection>
std::enable_if_t< std::is_base_of< Projection, PROJ >::value, const PROJ & > apply (const std::string &name, const Event &evt) const
 Apply the supplied projection on event evt (convenience arg-reordering alias).
 

Static Public Member Functions

Static helper functions for FastJet interaction, with tagging
static PseudoJets mkClusterInputs (const Particles &fsparticles, const Particles &tagparticles=Particles())
 Make PseudoJets for input to a ClusterSequence, with user_index codes for constituent- and tag-particle linking.
 
static Jet mkJet (const PseudoJet &pj, const Particles &fsparticles, const Particles &tagparticles=Particles())
 Make a Rivet Jet from a PseudoJet holding a user_index code for lookup of Rivet fsparticle or tagparticle links.
 
static Jets mkJets (const PseudoJets &pjs, const Particles &fsparticles, const Particles &tagparticles=Particles())
 Convert a whole list of PseudoJets to a list of Jets, with mkJet-style unpacking.
 

Protected Member Functions

void project (const Event &e)
 Perform the projection on the Event.
 
CmpState compare (const Projection &p) const
 Compare projections.
 
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 fail ()
 Set the projection in an unvalid state.
 
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.
 
void setProjectionHandler (ProjectionHandler &projectionHandler) const
 
Projection registration functions
template<typename PROJ >
const PROJ & declare (const PROJ &proj, const std::string &name) const
 Register a contained projection (user-facing version)
 
template<typename PROJ >
const PROJ & declare (const std::string &name, const PROJ &proj) const
 Register a contained projection (user-facing, arg-reordered version)
 

Detailed Description

Find jets using jet algorithms via the FastJet package.

Constructor & Destructor Documentation

◆ FastJets() [1/7]

Rivet::FastJets::FastJets ( const FinalState fsp,
const fastjet::JetDefinition &  jdef,
JetMuons  usemuons = JetMuons::ALL,
JetInvisibles  useinvis = JetInvisibles::NONE,
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.

◆ FastJets() [2/7]

Rivet::FastJets::FastJets ( const FinalState fsp,
const fastjet::JetDefinition &  jdef,
fastjet::AreaDefinition *  adef,
JetMuons  usemuons = JetMuons::ALL,
JetInvisibles  useinvis = JetInvisibles::NONE 
)
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.

◆ FastJets() [3/7]

Rivet::FastJets::FastJets ( const FinalState fsp,
fastjet::JetAlgorithm  type,
fastjet::RecombinationScheme  recom,
double  rparameter,
JetMuons  usemuons = JetMuons::ALL,
JetInvisibles  useinvis = JetInvisibles::NONE,
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.

◆ FastJets() [4/7]

Rivet::FastJets::FastJets ( const FinalState fsp,
fastjet::JetAlgorithm  type,
fastjet::RecombinationScheme  recom,
double  rparameter,
fastjet::AreaDefinition *  adef,
JetMuons  usemuons = JetMuons::ALL,
JetInvisibles  useinvis = JetInvisibles::NONE 
)
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.

◆ FastJets() [5/7]

Rivet::FastJets::FastJets ( const FinalState fsp,
fastjet::JetDefinition::Plugin *  plugin,
JetMuons  usemuons = JetMuons::ALL,
JetInvisibles  useinvis = JetInvisibles::NONE,
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

References reset().

◆ FastJets() [6/7]

Rivet::FastJets::FastJets ( const FinalState fsp,
fastjet::JetDefinition::Plugin *  plugin,
fastjet::AreaDefinition *  adef,
JetMuons  usemuons = JetMuons::ALL,
JetInvisibles  useinvis = JetInvisibles::NONE 
)
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

◆ FastJets() [7/7]

Rivet::FastJets::FastJets ( const FinalState fsp,
JetAlg  alg,
double  rparameter,
JetMuons  usemuons = JetMuons::ALL,
JetInvisibles  useinvis = JetInvisibles::NONE,
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

Member Function Documentation

◆ addTrf()

void Rivet::FastJets::addTrf ( fastjet::Transformer *  trf)
inline

Add a grooming transformer (base class of fastjet::Filter, etc.)

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

Referenced by addTrfs().

◆ addTrfs()

template<typename TRFS , typename TRF = typename TRFS::value_type>
std::enable_if< Derefable< TRF >::value, void >::type Rivet::FastJets::addTrfs ( const TRFS &  trfs)
inline

Add a list of grooming transformers.

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

References addTrf().

◆ apply()

template<typename PROJ = Projection>
std::enable_if_t< std::is_base_of< Projection, PROJ >::value, const PROJ & > Rivet::ProjectionApplier::apply ( const Event evt,
const Projection proj 
) const
inlineinherited

Apply the supplied projection on event evt.

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

Referenced by Rivet::ALICE::V0Trigger< MODE >::project().

◆ areaDef()

const shared_ptr< fastjet::AreaDefinition > Rivet::FastJets::areaDef ( ) const
inline

Return the area definition.

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

Referenced by clusterSeqArea().

◆ before()

bool Rivet::Projection::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.

◆ clone()

virtual unique_ptr< Projection > Rivet::JetFinder::clone ( ) const
pure virtualinherited

Clone on the heap.

Implements Rivet::Projection.

◆ clusterSeq()

const shared_ptr< fastjet::ClusterSequence > Rivet::FastJets::clusterSeq ( ) const
inline

Return the cluster sequence.

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

◆ clusterSeqArea()

const shared_ptr< fastjet::ClusterSequenceArea > Rivet::FastJets::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>

References areaDef().

◆ compare()

CmpState Rivet::FastJets::compare ( const Projection p) const
protectedvirtual

Compare projections.

Implements Rivet::JetFinder.

◆ declare() [1/2]

template<typename PROJ >
const PROJ & Rivet::ProjectionApplier::declare ( const PROJ &  proj,
const std::string &  name 
) const
inlineprotectedinherited

Register a contained projection (user-facing version)

Todo:
Add SFINAE to require that PROJ inherit from Projection

Referenced by Rivet::CentralEtHCM::CentralEtHCM(), Rivet::CentralityEstimator::CentralityEstimator(), Rivet::ChargedLeptons::ChargedLeptons(), Rivet::ALICE::CLMultiplicity< INNER >::CLMultiplicity(), Rivet::DISDiffHadron::DISDiffHadron(), Rivet::DISFinalState::DISFinalState(), Rivet::DISKinematics::DISKinematics(), Rivet::DISLepton::DISLepton(), Rivet::EventMixingBase::EventMixingBase(), Rivet::GammaGammaKinematics::GammaGammaKinematics(), Rivet::GammaGammaLeptons::GammaGammaLeptons(), Rivet::GammaGammaLeptons::GammaGammaLeptons(), Rivet::GeneratedCentrality::GeneratedCentrality(), Rivet::HadronicFinalState::HadronicFinalState(), Rivet::HeavyHadrons::HeavyHadrons(), Rivet::Hemispheres::Hemispheres(), Rivet::InvisibleFinalState::InvisibleFinalState(), Rivet::LeadingParticlesFinalState::LeadingParticlesFinalState(), Rivet::LossyFinalState< FILTER >::LossyFinalState(), Rivet::LossyFinalState< FILTER >::LossyFinalState(), Rivet::MC_pPbMinBiasTrigger::MC_pPbMinBiasTrigger(), Rivet::MC_SumETFwdPbCentrality::MC_SumETFwdPbCentrality(), Rivet::ATLAS::MinBiasTrigger::MinBiasTrigger(), Rivet::MissingMomentum::MissingMomentum(), Rivet::NeutralFinalState::NeutralFinalState(), Rivet::NeutralFinalState::NeutralFinalState(), Rivet::NonHadronicFinalState::NonHadronicFinalState(), Rivet::ParisiTensor::ParisiTensor(), Rivet::PercentileProjection::PercentileProjection(), Rivet::PrimaryHadrons::PrimaryHadrons(), Rivet::PrimaryHadrons::PrimaryHadrons(), Rivet::SmearedMET::SmearedMET(), Rivet::SmearedMET::SmearedMET(), Rivet::Spherocity::Spherocity(), Rivet::ATLAS::SumET_PB_Centrality::SumET_PB_Centrality(), Rivet::ATLAS::SumET_PBPB_Centrality::SumET_PBPB_Centrality(), Rivet::TauFinder::TauFinder(), Rivet::TriggerCDFRun0Run1::TriggerCDFRun0Run1(), Rivet::TriggerCDFRun2::TriggerCDFRun2(), Rivet::UndressBeamLeptons::UndressBeamLeptons(), Rivet::ALICE::V0AndTrigger::V0AndTrigger(), Rivet::ALICE::V0Trigger< MODE >::V0Trigger(), Rivet::VetoedFinalState::VetoedFinalState(), Rivet::VisibleFinalState::VisibleFinalState(), Rivet::VisibleFinalState::VisibleFinalState(), Rivet::CentralityProjection::add(), Rivet::CentralityBinner< T, MDist >::setProjection(), and Rivet::VetoedFinalState::vetoFinalState().

◆ declare() [2/2]

template<typename PROJ >
const PROJ & Rivet::ProjectionApplier::declare ( const std::string &  name,
const PROJ &  proj 
) const
inlineprotectedinherited

Register a contained projection (user-facing, arg-reordered version)

Todo:
Add SFINAE to require that PROJ inherit from Projection

◆ get()

template<typename PROJ >
const PROJ & Rivet::ProjectionApplier::get ( const std::string &  name) const
inlineinherited

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

Todo:
Add SFINAE to require that PROJ inherit from Projection

◆ getProjection() [1/2]

template<typename PROJ >
const PROJ & Rivet::ProjectionApplier::getProjection ( const std::string &  name) const
inlineinherited

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

Todo:
Add SFINAE to require that PROJ inherit from Projection

References Rivet::ProjectionHandler::getProjection(), and Rivet::ProjectionApplier::getProjHandler().

Referenced by Rivet::CentralityProjection::compare(), Rivet::pcmp(), Rivet::pcmp(), Rivet::pcmp(), and Rivet::pcmp().

◆ getProjection() [2/2]

const Projection & Rivet::ProjectionApplier::getProjection ( const std::string &  name) const
inlineinherited

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

References Rivet::ProjectionHandler::getProjection(), and Rivet::ProjectionApplier::getProjHandler().

◆ getProjectionFromDeclQueue()

template<typename PROJ >
const PROJ & Rivet::ProjectionApplier::getProjectionFromDeclQueue ( const std::string  name) const
inlineinherited

Get a named projection from this projection appliers declqueue TODO for TP: Recursion?

References MSG_ERROR.

◆ jets() [1/6]

Jets Rivet::JetFinder::jets ( const Cut &  c,
const JetSorter sorter 
) const
inlineinherited

Get the jets with a Cut applied, and ordered by supplied sorting functor

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?

References Rivet::JetFinder::jets(), and Rivet::sortBy().

◆ jets() [2/6]

virtual Jets Rivet::JetFinder::jets ( const Cut &  c = Cuts::open()) const
inlinevirtualinherited

◆ jets() [3/6]

virtual Jets Rivet::JetFinder::jets ( const JetSelector selector) const
inlinevirtualinherited

Get jets in no guaranteed order, with a selection functor

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

References Rivet::select().

◆ jets() [4/6]

Jets Rivet::JetFinder::jets ( const JetSelector selector,
const JetSorter sorter 
) const
inlineinherited

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?

References Rivet::JetFinder::jets(), and Rivet::sortBy().

◆ jets() [5/6]

Jets Rivet::JetFinder::jets ( const JetSorter sorter,
const Cut &  c = Cuts::open() 
) const
inlineinherited

Get the jets, ordered by supplied sorting functor, with an optional Cut

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?

References Rivet::JetFinder::jets().

◆ jets() [6/6]

Jets Rivet::JetFinder::jets ( const JetSorter sorter,
const JetSelector  selector 
) const
inlineinherited

Get the jets, ordered by supplied sorting functor and with a selection functor applied

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

References Rivet::JetFinder::jets().

◆ jetsByPt() [1/2]

Jets Rivet::JetFinder::jetsByPt ( const Cut &  c = Cuts::open()) const
inlineinherited

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).

References Rivet::cmpMomByPt(), and Rivet::JetFinder::jets().

◆ jetsByPt() [2/2]

Jets Rivet::JetFinder::jetsByPt ( const JetSelector selector) const
inlineinherited

Get the jets, ordered by \( p_T \), with cuts via a selection functor.

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).

References Rivet::cmpMomByPt(), and Rivet::JetFinder::jets().

◆ mkNamedPCmp()

◆ mkPCmp()

Cmp< Projection > Rivet::Projection::mkPCmp ( const Projection otherparent,
const std::string &  pname 
) const
protectedinherited

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

Note
Alias for mkNamedPCmp

Referenced by Rivet::ALICE::PrimaryParticles::compare(), Rivet::SmearedJets::compare(), Rivet::SmearedMET::compare(), Rivet::SmearedParticles::compare(), and Rivet::Correlators::compare().

◆ name()

virtual std::string Rivet::Projection::name ( ) const
inlinevirtualinherited

◆ project()

void Rivet::FastJets::project ( const Event e)
protectedvirtual

Perform the projection on the Event.

Implements Rivet::JetFinder.

◆ pseudojets()

PseudoJets Rivet::FastJets::pseudojets ( double  ptmin = 0.0) const

Get the pseudo jets (unordered).

Deprecated:
Use pseudojets

Referenced by pseudojetsByE(), pseudojetsByPt(), and pseudojetsByRapidity().

◆ reset()

void Rivet::FastJets::reset ( )
virtual

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

Implements Rivet::JetFinder.

Referenced by FastJets().

◆ setProjectionHandler()

void Rivet::ProjectionApplier::setProjectionHandler ( ProjectionHandler projectionHandler) const
protectedinherited
Todo:
AB: Add Doxygen comment, follow surrounding coding style

◆ useInvisibles()

void Rivet::JetFinder::useInvisibles ( JetInvisibles  useinvis = JetInvisibles::DECAY)
inlineinherited

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.

◆ useJetArea()

void Rivet::FastJets::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

◆ useMuons()

void Rivet::JetFinder::useMuons ( JetMuons  usemuons = JetMuons::ALL)
inlineinherited

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.


The documentation for this class was generated from the following file:
  • /Users/chrisg/software/rivet/include/Rivet/Projections/FastJets.hh