FastJet 3.0.2
|
Implementation of the spherical version of the SISCone algorithm (plugin for fastjet v2.1 upwards) More...
#include <fastjet/SISConeSphericalPlugin.hh>
Inherits fastjet::SISConeBasePlugin.
Public Types | |
enum | SplitMergeScale { SM_E, SM_Etilde } |
enum for the different split-merge scale choices; Note that order _must_ be the same as in siscone More... | |
Public Member Functions | |
SISConeSphericalPlugin (double cone_radius_in, double overlap_threshold_in, int n_pass_max_in=0, double protojet_Emin_in=0.0, bool caching_in=false, SplitMergeScale split_merge_scale_in=SM_Etilde, double split_merge_stopping_scale_in=0.0) | |
Main constructor for the SISConeSpherical Plugin class. | |
double | protojet_Emin () const |
minimum energy for a protojet to be considered in the split-merge step of the algorithm | |
double | protojet_or_ghost_Emin () const |
return the scale to be passed to SISCone as the protojet_Emin | |
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_use_E_weighted_splitting () const |
indicate if the splittings are done using the anti-kt distance | |
void | set_split_merge_use_E_weighted_splitting (bool val) |
virtual bool | supports_ghosted_passive_areas () const |
overload the default as we don't provide support for passive areas. | |
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: | |
Protected Member Functions | |
virtual void | reset_stored_plugin () const |
Implementation of the spherical version of the SISCone algorithm (plugin for fastjet v2.1 upwards)
SISConeSphericalPlugin 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.
This is the version of SISCone using spherical coordinates. Compared to the original cylindrical version:
For moderate values of R the problem should not be too severe (or may even be absent for some values of the overlap parameter), however the user should be aware of the issue.
The default split-merge scale may change at a later date to resolve this issue.
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 94 of file SISConeSphericalPlugin.hh.
enum for the different split-merge scale choices; Note that order _must_ be the same as in siscone
SM_E |
Energy (IR unsafe with momentum conservation) |
SM_Etilde |
sum_{i jet} E_i [1+sin^2(theta_iJ)] |
Definition at line 99 of file SISConeSphericalPlugin.hh.
double fastjet::SISConeSphericalPlugin::protojet_or_ghost_Emin | ( | ) | const [inline] |
return the scale to be passed to SISCone as the protojet_Emin
-- 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 133 of file SISConeSphericalPlugin.hh.
virtual bool fastjet::SISConeSphericalPlugin::supports_ghosted_passive_areas | ( | ) | const [inline, virtual] |
overload the default as we don't provide support for passive areas.
Reimplemented from fastjet::JetDefinition::Plugin.
Definition at line 148 of file SISConeSphericalPlugin.hh.
void fastjet::SISConeSphericalPlugin::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:
Implements fastjet::JetDefinition::Plugin.
Definition at line 75 of file SISConeSphericalPlugin.cc.