Rivet API documentation
![]() |
Rivet 4.1.3
|
FastJets.hh
260 static PseudoJets mkClusterInputs(const Particles& fsparticles, const Particles& tagparticles=Particles());
262 static Jet mkJet(const PseudoJet& pj, const Particles& fsparticles, const Particles& tagparticles=Particles());
264 static Jets mkJets(const PseudoJets& pjs, const Particles& fsparticles=Particles(), const Particles& tagparticles=Particles());
311 static CONTAINER reclusterJets(const CONTAINER &jetsIn, const fastjet::JetDefinition &jDef, const fastjet::Filter &filter){
345 static CONTAINER reclusterJets(const CONTAINER &jetsIn, const FJPluginPtr &jetAlg, const fastjet::Filter &filter){
384 static CONTAINER reclusterJets(const CONTAINER &jetsIn, const fastjet::Filter &filter, Args&&... args){
401 for ( auto const &[key, jetsIn] : jetsMap ) rtn[key] = reclusterJets(jetsIn, std::forward<Args>(args)...);
411 for ( auto const &[key, jetsIn] : jetsMap ) rtn[key] = reclusterJets<JETALG>(jetsIn, std::forward<Args>(args)...);
603 static const std::map<JetAlg, std::pair<fastjet::JetAlgorithm, fastjet::RecombinationScheme>> jetAlgMap;
610 template<typename X> struct mapJetAlg2Plugin<JetAlg::SISCONE, X>{ using type = fastjet::SISConePlugin; };
611 template<typename X> struct mapJetAlg2Plugin<JetAlg::PXCONE, X>{ using type = Rivet::PxConePlugin; };
612 template<typename X> struct mapJetAlg2Plugin<JetAlg::CDFJETCLU, X>{ using type = fastjet::CDFJetCluPlugin; };
613 template<typename X> struct mapJetAlg2Plugin<JetAlg::CDFMIDPOINT, X>{ using type = fastjet::CDFMidPointPlugin; };
614 template<typename X> struct mapJetAlg2Plugin<JetAlg::D0ILCONE, X>{ using type = fastjet::D0RunIIConePlugin; };
615 template<typename X> struct mapJetAlg2Plugin<JetAlg::JADE, X>{ using type = fastjet::JadePlugin; };
616 template<typename X> struct mapJetAlg2Plugin<JetAlg::TRACKJET, X>{ using type = fastjet::TrackJetPlugin; };
617 template<typename X> struct mapJetAlg2Plugin<JetAlg::VARIABLER, X>{ using type = fastjet::contrib::VariableRPlugin; };
Representation of a HepMC event, and enabler of Projection caching.
Definition Event.hh:22
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 tagpart...
FastJets(const FinalState &fsp, FJPluginPtr 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 ...
Definition FastJets.hh:173
static PJetsParts reclusterJetsParts(const Jets &jetsIn, const fastjet::JetDefinition &jDef)
Aux function for reclustering.
PseudoJets pseudojetsByE(double ptmin=0.0) const
Get the pseudo jets, ordered by .
Definition FastJets.hh:502
FastJets(const FinalState &fsp, JetAlg alg, double rparameter=-1, JetMuons usemuons=JetMuons::ALL, JetInvisibles useinvis=JetInvisibles::NONE, fastjet::AreaDefinition *adef=nullptr)
Convenience constructor using Rivet enums for most common jet algs (including some plugins).
Definition FastJets.hh:200
void clearTrfs()
Don't apply any jet transformers.
Definition FastJets.hh:479
FastJets(const FinalState &fsp, const fastjet::JetDefinition &jdef, fastjet::AreaDefinition *adef, JetMuons usemuons=JetMuons::ALL, JetInvisibles useinvis=JetInvisibles::NONE)
Definition FastJets.hh:71
std::enable_if< Derefable< FILTER >::value, void >::type addFilters(const FILTERS &filters)
Add a list of grooming filters.
Definition FastJets.hh:474
PseudoJets pseudojetsByRapidity(double ptmin=0.0) const
Get the pseudo jets, ordered by rapidity.
Definition FastJets.hh:507
FastJets(const FinalState &fsp, const fastjet::JetDefinition &jdef, JetMuons usemuons=JetMuons::ALL, JetInvisibles useinvis=JetInvisibles::NONE, fastjet::AreaDefinition *adef=nullptr)
Definition FastJets.hh:42
FastJets(const FinalState &fsp, const fastjet::JetDefinition &jdef, const Cut &c, JetMuons usemuons=JetMuons::ALL, JetInvisibles useinvis=JetInvisibles::NONE, fastjet::AreaDefinition *adef=nullptr)
Definition FastJets.hh:56
FastJets(const FinalState &fsp, fastjet::JetAlgorithm type, fastjet::RecombinationScheme recom, double rparameter, JetMuons usemuons=JetMuons::ALL, JetInvisibles useinvis=JetInvisibles::NONE, fastjet::AreaDefinition *adef=nullptr)
Definition FastJets.hh:96
FastJets(const FinalState &fsp, fastjet::JetAlgorithm type, fastjet::RecombinationScheme recom, double rparameter, const Cut &c, JetMuons usemuons=JetMuons::ALL, JetInvisibles useinvis=JetInvisibles::NONE, fastjet::AreaDefinition *adef=nullptr)
Definition FastJets.hh:108
std::enable_if< Derefable< TRF >::value, void >::type addTrfs(const TRFS &trfs)
Add a list of grooming transformers.
Definition FastJets.hh:458
FastJets(const FinalState &fsp, const fastjet::JetDefinition &jdef, fastjet::AreaDefinition *adef, const Cut &c, JetMuons usemuons=JetMuons::ALL, JetInvisibles useinvis=JetInvisibles::NONE)
Definition FastJets.hh:83
void clearJetArea()
Don't calculate a jet area.
Definition FastJets.hh:434
const shared_ptr< fastjet::ClusterSequence > clusterSeq() const
Definition FastJets.hh:519
void addFilter(fastjet::Filter *filter)
Add a grooming filter.
Definition FastJets.hh:465
static FJPluginPtr mkPlugin(JetAlg alg, double rparameter=-1)
Non-templated version only allowing to set rparameter.
FastJets(const FinalState &fsp, FJPluginPtr plugin, const Cut &c, JetMuons usemuons=JetMuons::ALL, JetInvisibles useinvis=JetInvisibles::NONE, fastjet::AreaDefinition *adef=nullptr)
Explicitly pass in an externally-constructed plugin.
Definition FastJets.hh:159
FastJets(const FinalState &fsp, FJPluginPtr plugin, fastjet::AreaDefinition *adef, const Cut &c, JetMuons usemuons=JetMuons::ALL, JetInvisibles useinvis=JetInvisibles::NONE)
Explicitly pass in an externally-constructed plugin, with reordered args for easier specification of ...
Definition FastJets.hh:184
PseudoJets pseudojets(double ptmin=0.0) const
FastJets(const FinalState &fsp, JetAlg alg, double rparameter, const Cut &c, JetMuons usemuons=JetMuons::ALL, JetInvisibles useinvis=JetInvisibles::NONE, fastjet::AreaDefinition *adef=nullptr)
Convenience constructor using Cut argument and Rivet enums for most common jet algs (including some p...
Definition FastJets.hh:219
FastJets(const FinalState &fsp, FJPluginPtr plugin, JetMuons usemuons=JetMuons::ALL, JetInvisibles useinvis=JetInvisibles::NONE, fastjet::AreaDefinition *adef=nullptr)
Explicitly pass in an externally-constructed plugin.
Definition FastJets.hh:146
static fastjet::JetDefinition mkJetDef(JetAlg alg, double rparameter)
Make a FastJet JetDefinition according to enum JetAlg.
static Jets mkTaggedJets(const Jets &jetsIn, const PJetsParts &pJetsParts)
Aux function for reclustering.
void calc(const Particles &fsparticles, const Particles &tagparticles=Particles())
Do the calculation locally (no caching).
typename mapJetAlg2Plugin< JETALG, void >::type mapJetAlg2Plugin_t
Make usage of "map" a bit more convenient.
Definition FastJets.hh:621
static const std::map< JetAlg, std::pair< fastjet::JetAlgorithm, fastjet::RecombinationScheme > > jetAlgMap
Map of Rivet::JetAlg to targeted fastjet::JetAlgorithm and fastjet::RecombinationScheme.
Definition FastJets.hh:603
const shared_ptr< fastjet::ClusterSequenceArea > clusterSeqArea() const
Definition FastJets.hh:525
void useJetArea(fastjet::AreaDefinition *adef)
Use provided jet area definition.
Definition FastJets.hh:429
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-partic...
const fastjet::JetDefinition & jetDef() const
Return the jet definition.
Definition FastJets.hh:530
void addTrf(fastjet::Transformer *trf)
Add a grooming transformer (base class of fastjet::Filter, etc.).
Definition FastJets.hh:448
static Jets mkJets(const PseudoJets &pjs, const Particles &fsparticles=Particles(), const Particles &tagparticles=Particles())
Convert a whole list of PseudoJets to a list of Jets, with mkJet-style unpacking.
PseudoJets pseudojetsByPt(double ptmin=0.0) const
Get the pseudo jets, ordered by .
Definition FastJets.hh:497
FastJets(const FinalState &fsp, fastjet::JetAlgorithm type, fastjet::RecombinationScheme recom, double rparameter, fastjet::AreaDefinition *adef, JetMuons usemuons=JetMuons::ALL, JetInvisibles useinvis=JetInvisibles::NONE)
Definition FastJets.hh:121
const shared_ptr< fastjet::AreaDefinition > areaDef() const
Return the area definition.
Definition FastJets.hh:538
FastJets(const FinalState &fsp, fastjet::JetAlgorithm type, fastjet::RecombinationScheme recom, double rparameter, fastjet::AreaDefinition *adef, const Cut &c, JetMuons usemuons=JetMuons::ALL, JetInvisibles useinvis=JetInvisibles::NONE)
Definition FastJets.hh:133
Project out all final-state particles in an event. Probably the most important projection in Rivet!
Definition FinalState.hh:12
JetFinder(const FinalState &fs, JetMuons usemuons=JetMuons::ALL, JetInvisibles useinvis=JetInvisibles::NONE)
Constructor.
const PROJ & get(const std::string &name) const
Definition ProjectionApplier.hh:83
PseudoJet & ifilterPseudoJets(PseudoJet &pj, const fastjet::Filter &filter)
Apply given FastJet::Filter.
static CONTAINER reclusterJets(const CONTAINER &jetsIn, const fastjet::JetDefinition &jDef, const fastjet::Filter &filter)
Recluster Rivet::Jets.
Definition FastJets.hh:311
static CONTAINER reclusterJets(const CONTAINER &jetsIn, Args &&... args)
Apply the.
Definition FastJets.hh:360
static CONTAINER reclusterJets(const CONTAINER &jetsIn, const fastjet::Filter &filter, Args &&... args)
Apply the.
Definition FastJets.hh:384
static CONTAINER reclusterJets(const CONTAINER &jetsIn, const fastjet::JetDefinition &jDef)
Recluster Rivet::Jets.
Definition FastJets.hh:291
static CONTAINER reclusterJets(const CONTAINER &jetsIn, const FJPluginPtr &jetAlg)
Apply the.
Definition FastJets.hh:327
static CONTAINER reclusterJets(const CONTAINER &jetsIn, const FJPluginPtr &jetAlg, const fastjet::Filter &filter)
Apply the.
Definition FastJets.hh:345
static std::map< T, U > reclusterJets(const std::map< T, U > &jetsMap, Args &&... args)
Apply the.
Definition FastJets.hh:409
static std::map< T, U > reclusterJets(const std::map< T, U > &jetsMap, Args &&... args)
Apply the.
Definition FastJets.hh:399
#define MSG_DEBUG(x)
Debug messaging, not enabled by default, using MSG_LVL.
Definition Logging.hh:182
double p(const ParticleBase &p)
Unbound function access to p.
Definition ParticleBaseUtils.hh:653
Definition MC_CENT_PPB_Projections.hh:10
JetInvisibles
Enum for the treatment of invisible particles: whether to include all, some, or none in jet-finding.
Definition JetFinder.hh:18
JetMuons
Enum for the treatment of muons: whether to include all, some, or none in jet-finding.
Definition JetFinder.hh:15
Definition FastJets.hh:607
Generated on for Rivet by
