1#ifndef __SISCONEBASEPLUGIN_HH__
2#define __SISCONEBASEPLUGIN_HH__
4#include "fastjet/JetDefinition.hh"
5#include "fastjet/ClusterSequence.hh"
13FASTJET_BEGIN_NAMESPACE
44 _use_jet_def_recombiner =
false;
45 set_progressive_removal(
false);
68 _progressive_removal = progressive_removal_in;
92 _split_merge_stopping_scale = scale;}
108 virtual double R()
const {
return cone_radius();}
122 _ghost_sep_scale = scale;
125 virtual double ghost_separation_scale()
const {
126 return _ghost_sep_scale;
199 double _cone_radius, _overlap_threshold;
202 double _split_merge_stopping_scale;
203 bool _use_jet_def_recombiner;
204 bool _progressive_removal;
206 mutable double _ghost_sep_scale;
239 const std::vector<PseudoJet> &
stable_cones()
const {
return _protocones;}
242 const std::vector<PseudoJet> &
protocones()
const {
return _protocones;}
251 std::ostringstream ostr;
252 ostr <<
"This SISCone clustering found " << protocones().size()
253 <<
" stable protocones";
263 std::vector<PseudoJet> _protocones;
264 std::vector<int> _pass;
265 double _most_ambiguous_split;
270inline SISConeBaseExtras::~SISConeBaseExtras(){}
293 return "PseudoJet wrapping a siscone jet from a stable cone";
304 std::vector<PseudoJet> constits;
305 constits.reserve(size());
306 for (
unsigned int i=0; i<size();i++)
307 constits.push_back(constituent(i));
315 virtual unsigned int size()
const = 0;
323 return _cs.jets()[constituent_index(i)];
354 virtual unsigned int size()
const{
361 return _jet.contents[i];
base class providing interface for a generic function of a PseudoJet
a class that allows a user to introduce their own "plugin" jet finder
Contains any information related to the clustering that should be directly accessible to PseudoJet.
Class to contain pseudojets, including minimal information of use to jet-clustering routines.
int cluster_hist_index() const
return the cluster_hist_index, intended to be used by clustering routines.
template class derived from UserScaleBase::StryctureType that works for both SISCone jet classes impl...
virtual double ordering_var2() const
returns the sm_var2 (signed ordering variable squared) for this stable cone
virtual ~UserScaleBaseStructureType()
empty virtual dtor
virtual int constituent_index(unsigned int i) const
returns the index (in the original particle list) of the ith constituent
const Tjet & _jet
a reference to the internal jet in SISCone
virtual unsigned int size() const
returns the number of constituents
the structure that allows to store the information contained into a siscone::Cjet (built internally i...
virtual bool has_constituents() const
this structure has constituents
virtual unsigned int size() const =0
returns the number of constituents
virtual std::vector< PseudoJet > constituents(const PseudoJet &) const
retrieve the constituents
virtual std::string description() const
the textual descripotion
virtual int constituent_index(unsigned int i) const =0
returns the index (in the original particle list) of the ith constituent
StructureType(const ClusterSequence &cs)
base ctor (constructed from a ClusterSequence tin order to have access to the initial particles
virtual double ordering_var2() const =0
returns the sm_var2 (signed ordering variable squared) for this stable cone
const ClusterSequence & _cs
a reference to the CS (for access to the particles)
const PseudoJet & constituent(unsigned int i) const
returns the ith constituent (as a PseusoJet)
virtual ~StructureType()
empty virtual dtor
base class for user-defined ordering of stable cones (used for prorgessive removal)
virtual double result(const PseudoJet &jet) const =0
returns the scale associated with a given jet
virtual ~UserScaleBase()
empty virtual dtor
double split_merge_stopping_scale()
return the value of the split_merge_stopping_scale (see set_split_merge_stopping_scale(....
virtual void set_ghost_separation_scale(double scale) const
set the ghost separation scale for passive area determinations just in the next run (strictly speakin...
void set_user_scale(const UserScaleBase *user_scale_in)
set a user-defined scale for stable-cone ordering in progressive removal
bool use_jet_def_recombiner() const
indicate if the jet_def's recombination scheme is being used
void set_use_jet_def_recombiner(bool choice)
allow the user to decide if one uses the jet_def's own recombination scheme
double cone_radius() const
the cone radius
bool caching() const
indicates whether caching is turned on or not.
SISConeBasePlugin()
default ctor
void set_progressive_removal(bool progressive_removal_in=true)
set whether to use SISCone with progressive removal instead of the default split_merge step.
const UserScaleBase * user_scale() const
returns the user-defined scale in use (0 if none)
SISConeBasePlugin(const SISConeBasePlugin &plugin)
copy constructor
virtual void run_clustering(ClusterSequence &) const =0
really do the clustering work
virtual bool supports_ghosted_passive_areas() const
return true since there is specific support for the measurement of passive areas, in the sense that a...
virtual double R() const
the plugin mechanism's standard way of accessing the jet 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).
virtual std::string description() const =0
plugin description
void set_split_merge_stopping_scale(double scale)
set the "split_merge_stopping_scale": if the scale variable for all protojets is below this,...
bool progressive_removal() const
returns true if progressive_removal is enabled
virtual void reset_stored_plugin() const =0
call the re-clustering itself