FastJet  3.1.3
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Modules Pages
Public Types | Public Member Functions | Protected Member Functions | List of all members
fastjet::SISConePlugin Class Reference

Implementation of the SISCone algorithm (plugin for fastjet v2.1 upwards) More...

#include <fastjet/SISConePlugin.hh>

Inherits fastjet::SISConeBasePlugin.

Public Types

enum  SplitMergeScale { SM_pt, SM_Et, SM_mt, SM_pttilde }
 enum for the different split-merge scale choices; Note that order must be the same as in siscone More...
 

Public Member Functions

 SISConePlugin (double cone_radius_in, double overlap_threshold_in, int n_pass_max_in=0, double protojet_ptmin_in=0.0, bool caching_in=false, SplitMergeScale split_merge_scale_in=SM_pttilde, double split_merge_stopping_scale_in=0.0)
 Main constructor for the SISCone Plugin class. More...
 
 SISConePlugin (double cone_radius_in, double overlap_threshold_in, int n_pass_max_in, double protojet_ptmin_in, bool caching_in, bool split_merge_on_transverse_mass_in)
 Backwards compatible constructor for the SISCone Plugin class.
 
 SISConePlugin (double cone_radius_in, double overlap_threshold_in, int n_pass_max_in, bool caching_in)
 backwards compatible constructor for the SISCone Plugin class (avoid using this in future). More...
 
double protojet_ptmin () const
 minimum pt for a protojet to be considered in the split-merge step of the algorithm
 
double protojet_or_ghost_ptmin () const
 return the scale to be passed to SISCone as the protojet_ptmin – if we have a ghost separation scale that is above the protojet_ptmin, then the ghost_separation_scale becomes the relevant one to use here
 
SplitMergeScale split_merge_scale () const
 indicates scale used in split-merge
 
void set_split_merge_scale (SplitMergeScale sms)
 sets scale used in split-merge
 
bool split_merge_on_transverse_mass () const
 indicates whether the split-merge orders on transverse mass or not. More...
 
void set_split_merge_on_transverse_mass (bool val)
 
bool split_merge_use_pt_weighted_splitting () const
 indicates whether the split-merge orders on transverse mass or not. More...
 
void set_split_merge_use_pt_weighted_splitting (bool val)
 
virtual std::string description () const
 return a textual description of the jet-definition implemented in this plugin
 
virtual void run_clustering (ClusterSequence &) const
 given a ClusterSequence that has been filled up with initial particles, the following function should fill up the rest of the ClusterSequence, using the following member functions of ClusterSequence: More...
 

Protected Member Functions

virtual void reset_stored_plugin () const
 

Detailed Description

Implementation of the SISCone algorithm (plugin for fastjet v2.1 upwards)

SISConePlugin is a plugin for fastjet (v2.1 upwards) that provides an interface to the seedless infrared safe cone jet finder by Gregory Soyez and Gavin Salam.

SISCone uses geometrical techniques to exhaustively consider all possible distinct cones. It then finds out which ones are stable and sends the result to the Tevatron Run-II type split-merge procedure for overlapping cones.

Four parameters govern the "physics" of the algorithm:

One parameter governs some internal algorithmic shortcuts:

The final jets can be accessed by requestion the inclusive_jets(...) from the ClusterSequence object. Note that these PseudoJets have their user_index() set to the index of the pass in which they were found (first pass = 0). NB: This does not currently work for jets that consist of a single particle.

For further information on the details of the algorithm see the SISCone paper, arXiv:0704.0292 [JHEP 0705:086,2007].

For documentation about the implementation, see the siscone/doc/html/index.html file.

Definition at line 72 of file SISConePlugin.hh.

Member Enumeration Documentation

enum for the different split-merge scale choices; Note that order must be the same as in siscone

Enumerator
SM_pt 

transverse momentum (E-scheme), IR unsafe

SM_Et 

transverse energy (E-scheme), not long.

boost invariant original run-II choice [may not be implemented]

SM_mt 

transverse mass (E-scheme), IR safe except in decays of two identical narrow heavy particles

SM_pttilde 

pt-scheme pt = {i in jet} |p_{ti}|, should be IR safe in all cases

Definition at line 77 of file SISConePlugin.hh.

Constructor & Destructor Documentation

fastjet::SISConePlugin::SISConePlugin ( double  cone_radius_in,
double  overlap_threshold_in,
int  n_pass_max_in = 0,
double  protojet_ptmin_in = 0.0,
bool  caching_in = false,
SplitMergeScale  split_merge_scale_in = SM_pttilde,
double  split_merge_stopping_scale_in = 0.0 
)
inline

Main constructor for the SISCone Plugin class.

Note: wrt version prior to 2.4 this constructor differs in that a the default value has been removed for overlap_threshold. The former has been removed because the old default of 0.5 was found to be unsuitable in high-noise environments; so the user should now explicitly think about the value for this – we recommend 0.75.

Definition at line 97 of file SISConePlugin.hh.

fastjet::SISConePlugin::SISConePlugin ( double  cone_radius_in,
double  overlap_threshold_in,
int  n_pass_max_in,
bool  caching_in 
)
inline

backwards compatible constructor for the SISCone Plugin class (avoid using this in future).

Definition at line 135 of file SISConePlugin.hh.

Member Function Documentation

bool fastjet::SISConePlugin::split_merge_on_transverse_mass ( ) const
inline

indicates whether the split-merge orders on transverse mass or not.

retained for backwards compatibility with 2.1.0b3

Definition at line 168 of file SISConePlugin.hh.

bool fastjet::SISConePlugin::split_merge_use_pt_weighted_splitting ( ) const
inline

indicates whether the split-merge orders on transverse mass or not.

retained for backwards compatibility with 2.1.0b3

Definition at line 174 of file SISConePlugin.hh.

void fastjet::SISConePlugin::run_clustering ( ClusterSequence ) const
virtual

given a ClusterSequence that has been filled up with initial particles, the following function should fill up the rest of the ClusterSequence, using the following member functions of ClusterSequence:

  • plugin_do_ij_recombination(...)
  • plugin_do_iB_recombination(...)

Implements fastjet::JetDefinition::Plugin.

Definition at line 139 of file SISConePlugin.cc.


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