FastJet 3.4.0
Loading...
Searching...
No Matches
fastjet::SISConePlugin Class Reference

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

#include <fastjet/SISConePlugin.hh>

Inheritance diagram for fastjet::SISConePlugin:
Collaboration diagram for fastjet::SISConePlugin:

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.
 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).
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.
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.
void set_split_merge_use_pt_weighted_splitting (bool val)
virtual std::string description () const
 plugin description
virtual void run_clustering (ClusterSequence &) const
 really do the clustering work
Public Member Functions inherited from fastjet::SISConeBasePlugin
 SISConeBasePlugin ()
 default ctor
 SISConeBasePlugin (const SISConeBasePlugin &plugin)
 copy constructor
void set_progressive_removal (bool progressive_removal_in=true)
 set whether to use SISCone with progressive removal instead of the default split_merge step.
bool progressive_removal () const
 returns true if progressive_removal is enabled
double cone_radius () const
 the cone radius
double overlap_threshold () const
 Fraction of overlap energy in a jet above which jets are merged and below which jets are split.
int n_pass_max () const
 the maximum number of passes of stable-cone searching (<=0 is same as infinity).
void set_split_merge_stopping_scale (double scale)
 set the "split_merge_stopping_scale": if the scale variable for all protojets is below this, then stop the split-merge procedure and keep only those jets found so far.
double split_merge_stopping_scale ()
 return the value of the split_merge_stopping_scale (see set_split_merge_stopping_scale(...) for description)
void set_use_jet_def_recombiner (bool choice)
 allow the user to decide if one uses the jet_def's own recombination scheme
bool use_jet_def_recombiner () const
 indicate if the jet_def's recombination scheme is being used
bool caching () const
 indicates whether caching is turned on or not.
virtual double R () const
 the plugin mechanism's standard way of accessing the jet radius
virtual bool supports_ghosted_passive_areas () const
 return true since there is specific support for the measurement of passive areas, in the sense that areas determined from all particles below the ghost separation scale will be a passive area.
virtual void set_ghost_separation_scale (double scale) const
 set the ghost separation scale for passive area determinations just in the next run (strictly speaking that makes the routine a non const, so related internal info must be stored as a mutable)
virtual double ghost_separation_scale () const
void set_user_scale (const UserScaleBase *user_scale_in)
 set a user-defined scale for stable-cone ordering in progressive removal
const UserScaleBaseuser_scale () const
 returns the user-defined scale in use (0 if none)
Public Member Functions inherited from fastjet::JetDefinition::Plugin
virtual bool exclusive_sequence_meaningful () const
 if this returns false then a warning will be given whenever the user requests "exclusive" jets from the cluster sequence
virtual bool is_spherical () const
 returns true if the plugin implements an algorithm intended for use on a spherical geometry (e.g.
virtual ~Plugin ()
 a destructor to be replaced if necessary in derived classes...

Protected Member Functions

virtual void reset_stored_plugin () const
 call the re-clustering itself

Additional Inherited Members

Protected Attributes inherited from fastjet::SISConeBasePlugin
double _cone_radius
double _overlap_threshold
int _n_pass_max
bool _caching
double _split_merge_stopping_scale
bool _use_jet_def_recombiner
bool _progressive_removal
double _ghost_sep_scale
const UserScaleBase_user_scale

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:

  • the cone_radius (this should be self-explanatory!)
  • the overlap_threshold is the parameter which dictates how much two jets must overlap (pt_overlap/min(pt1,pt2)) if they are to be merged
  • Not all particles are in stable cones in the first round of searching for stable cones; one can therefore optionally have the the jet finder carry out additional passes of searching for stable cones among particles that were in no stable cone in previous passes — the maximum number of passes carried out is n_pass_max. If this is zero then additional passes are carried out until no new stable cones are found.
  • Protojet ptmin: protojets that are below this ptmin (default = 0) are discarded before each iteration of the split-merge loop.

One parameter governs some internal algorithmic shortcuts:

  • if "caching" is turned on then the last event clustered by siscone is stored – if the current event is identical and the cone_radius and n_pass_mass are identical, then the only part of the clustering that needs to be rerun is the split-merge part, leading to significant speed gains; there is a small (O(N) storage and speed) penalty for caching, so it should be kept off (default) if only a single overlap_threshold is used.

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

◆ SplitMergeScale

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 = \sum_{i in jet} |p_{ti}|, should be IR safe in all cases

Definition at line 77 of file SISConePlugin.hh.

Constructor & Destructor Documentation

◆ SISConePlugin() [1/3]

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.

◆ SISConePlugin() [2/3]

fastjet::SISConePlugin::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 )
inline

Backwards compatible constructor for the SISCone Plugin class.

Definition at line 117 of file SISConePlugin.hh.

◆ SISConePlugin() [3/3]

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

◆ protojet_ptmin()

double fastjet::SISConePlugin::protojet_ptmin ( ) const
inline

minimum pt for a protojet to be considered in the split-merge step of the algorithm

Definition at line 152 of file SISConePlugin.hh.

◆ protojet_or_ghost_ptmin()

double fastjet::SISConePlugin::protojet_or_ghost_ptmin ( ) const
inline

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

Definition at line 158 of file SISConePlugin.hh.

◆ split_merge_scale()

SplitMergeScale fastjet::SISConePlugin::split_merge_scale ( ) const
inline

indicates scale used in split-merge

Definition at line 162 of file SISConePlugin.hh.

◆ set_split_merge_scale()

void fastjet::SISConePlugin::set_split_merge_scale ( SplitMergeScale sms)
inline

sets scale used in split-merge

Definition at line 164 of file SISConePlugin.hh.

◆ split_merge_on_transverse_mass()

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.

◆ set_split_merge_on_transverse_mass()

void fastjet::SISConePlugin::set_split_merge_on_transverse_mass ( bool val)
inline

Definition at line 169 of file SISConePlugin.hh.

◆ split_merge_use_pt_weighted_splitting()

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.

◆ set_split_merge_use_pt_weighted_splitting()

void fastjet::SISConePlugin::set_split_merge_use_pt_weighted_splitting ( bool val)
inline

Definition at line 175 of file SISConePlugin.hh.

◆ description()

string fastjet::SISConePlugin::description ( ) const
virtual

plugin description

Implements fastjet::SISConeBasePlugin.

Definition at line 88 of file SISConePlugin.cc.

◆ run_clustering()

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

really do the clustering work

Implements fastjet::SISConeBasePlugin.

Definition at line 139 of file SISConePlugin.cc.

◆ reset_stored_plugin()

void fastjet::SISConePlugin::reset_stored_plugin ( ) const
protectedvirtual

call the re-clustering itself

Implements fastjet::SISConeBasePlugin.

Definition at line 297 of file SISConePlugin.cc.


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