FastJet 3.4.1
SISConeBasePlugin.hh
1#ifndef __SISCONEBASEPLUGIN_HH__
2#define __SISCONEBASEPLUGIN_HH__
3
4#include "fastjet/JetDefinition.hh"
5#include "fastjet/ClusterSequence.hh"
6#include <vector>
7#include <memory>
8#include <cmath>
9
10#include <sstream>
11
12// questionable whether this should be in fastjet namespace or not...
13FASTJET_BEGIN_NAMESPACE // defined in fastjet/internal/base.hh
14
15//----------------------------------------------------------------------
16//
17/// \if internal_doc
18/// @ingroup internal
19/// \class SISConeBasePlugin
20/// Implementation of the SISCone algorithm, base class (plugin for fastjet v2.1 upwards)
21///
22/// SISConeBasePlugin is a plugin for fastjet (v2.1 upwards) that
23/// provides a base interface to SISCone-type cone jet finder by
24/// Gregory Soyez and Gavin Salam.
25///
26/// This is a purely virtual class that needs to be overloaded
27/// for the specific implementations of SISCone (i.e. regular or
28/// spherical as of July 16th 2008).
29///
30/// any derived plugin MUST overload the following methods:
31/// description()
32/// run_siscone_clustering()
33/// reset_stored_plugin()
34///
35/// For further details, see the derived plugins or
36/// http://projects.hepforge.com/siscone
37///
38/// \endif
39//
41public:
42 /// default ctor
44 _use_jet_def_recombiner = false;
45 set_progressive_removal(false);
46 }
47
48 /// copy constructor
50 *this = plugin;
51 }
52
53 /// set whether to use SISCone with progressive removal instead of
54 /// the default split_merge step.
55 ///
56 /// If progressive removal is enabled, the following SISCone
57 /// variables are not used:
58 ///
59 /// - overlap_threshold
60 /// - caching
61 /// - split_merge_stopping_scale
62 ///
63 /// The split_merge_scale choice is reinterpreted as the ordering
64 /// variable for progressive removal. It is also possible for the
65 /// user to supply his/her own function for the scale that orders
66 /// progressive removal, with set_user_scale(...)
67 void set_progressive_removal(bool progressive_removal_in=true){
68 _progressive_removal = progressive_removal_in;
69 }
70
71 /// returns true if progressive_removal is enabled
72 bool progressive_removal() const{ return _progressive_removal;}
73
74 /// the cone radius
75 double cone_radius () const {return _cone_radius ;}
76
77 /// Fraction of overlap energy in a jet above which jets are merged
78 /// and below which jets are split.
79 double overlap_threshold () const {return _overlap_threshold ;}
80
81 /// the maximum number of passes of stable-cone searching (<=0 is same
82 /// as infinity).
83 int n_pass_max () const {return _n_pass_max ;}
84
85 /// set the "split_merge_stopping_scale": if the scale variable for
86 /// all protojets is below this, then stop the split-merge procedure
87 /// and keep only those jets found so far. This is useful in
88 /// determination of areas of hard jets because it can be used to
89 /// avoid running the split-merging on the pure ghost-part of the
90 /// event.
91 void set_split_merge_stopping_scale(double scale) {
92 _split_merge_stopping_scale = scale;}
93
94 /// return the value of the split_merge_stopping_scale (see
95 /// set_split_merge_stopping_scale(...) for description)
96 double split_merge_stopping_scale() {return _split_merge_stopping_scale;}
97
98 /// allow the user to decide if one uses the jet_def's own recombination scheme
99 void set_use_jet_def_recombiner(bool choice) {_use_jet_def_recombiner = choice;}
100
101 /// indicate if the jet_def's recombination scheme is being used
102 bool use_jet_def_recombiner() const {return _use_jet_def_recombiner;}
103
104 /// indicates whether caching is turned on or not.
105 bool caching() const {return _caching ;}
106
107 /// the plugin mechanism's standard way of accessing the jet radius
108 virtual double R() const {return cone_radius();}
109
110 /// return true since there is specific support for the measurement
111 /// of passive areas, in the sense that areas determined from all
112 /// particles below the ghost separation scale will be a passive
113 /// area.
114 virtual bool supports_ghosted_passive_areas() const {
115 return true;
116 }
117
118 /// set the ghost separation scale for passive area determinations
119 /// _just_ in the next run (strictly speaking that makes the routine
120 /// a non const, so related internal info must be stored as a mutable)
121 virtual void set_ghost_separation_scale(double scale) const {
122 _ghost_sep_scale = scale;
123 }
124
125 virtual double ghost_separation_scale() const {
126 return _ghost_sep_scale;
127 }
128
129 // user-defined scale for progressive removal
130 //------------------------------------------------------------
131
132 /// \class UserScaleBase
133 /// base class for user-defined ordering of stable cones (used for
134 /// prorgessive removal)
135 ///
136 /// derived classes have to implement the () operator that returns
137 /// the scale associated with a given jet.
138 ///
139 /// It is also highly recommended to implement the is_larger()
140 /// method whenever possible, in order to avoid rounding issues
141 /// known to lead to possible infrared unsafeties.
142 ///
143 /// The jets that are passed to this class will carry the structure
144 /// of type SISConePlugin::StructureType which allows to retreive
145 /// easily the following information:
146 ///
147 /// vector<PseudoJet> constituents = jet.constituents();
148 /// unsigned int n_constituents = jet.structure_of<SISConePlugin::UserScaleBase>().size();
149 /// int index = jet.structure_of<SISConePlugin::UserScaleBase>().constituent_index(index i);
150 /// const PseudoJet & p = jet.structure_of<SISConePlugin::UserScaleBase>().constituent(index i);
151 /// double scalar_pt = jet.structure_of<SISConePlugin::UserScaleBase>().pt_tilde();
152 ///
153 /// see SISConePlugin::StructureType below for further details
154 class UserScaleBase : public FunctionOfPseudoJet<double>{
155 public:
156 /// empty virtual dtor
157 virtual ~UserScaleBase(){}
158
159 /// returns the scale associated with a given jet
160 ///
161 /// "progressive removal" iteratively removes the stable cone with
162 /// the largest scale
163 virtual double result(const PseudoJet & jet) const = 0;
164
165 /// returns true when the scale associated with jet a is larger than
166 /// the scale associated with jet b
167 ///
168 /// By default this does a simple direct comparison but it can be
169 /// overloaded for higher precision [recommended if possible]
170 virtual bool is_larger(const PseudoJet & a, const PseudoJet & b) const;
171
172 class StructureType; // defined below
173 };
174
175 // template class derived from UserScaleBase::StryctureType that
176 // works for both SISCone jet classes
177 // implemented below
178 template<class Tjet>
180
181 /// set a user-defined scale for stable-cone ordering in
182 /// progressive removal
183 void set_user_scale(const UserScaleBase *user_scale_in){ _user_scale = user_scale_in;}
184
185 /// returns the user-defined scale in use (0 if none)
186 const UserScaleBase * user_scale() const{ return _user_scale;}
187
188
189 // the things that one MUST overload required by base class
190 //---------------------------------------------------------
191
192 /// plugin description
193 virtual std::string description () const =0;
194
195 /// really do the clustering work
196 virtual void run_clustering(ClusterSequence &) const = 0;
197
198protected:
199 double _cone_radius, _overlap_threshold;
200 int _n_pass_max;
201 bool _caching;//, _split_merge_on_transverse_mass;
202 double _split_merge_stopping_scale;
203 bool _use_jet_def_recombiner;
204 bool _progressive_removal;
205
206 mutable double _ghost_sep_scale;
207
208 // the part that HAS to be overloaded
209 /// call the re-clustering itself
210 virtual void reset_stored_plugin() const =0;
211
212 const UserScaleBase * _user_scale;
213
214};
215
216
217//======================================================================
218/// @ingroup extra_info
219/// \class SISConeBaseExtras
220/// Class that provides extra information about a SISCone clustering
221///
222/// This is only the base class that the "regular" and "spherical"
223/// implementations of SISCone will have to overload. The only thing
224/// that needs to be done for the derived classes is to define
225/// '_jet_def_plugin', implement
226/// jet_def_plugin();
227/// and add the corresponding plugin class as a friend
229public:
230
231 /// constructor
232 // it just initialises the pass information
233 SISConeBaseExtras(int nparticles) : _pass(nparticles*2,-1) {}
234
235 /// purely virtual destructor
236 inline virtual ~SISConeBaseExtras()=0;
237
238 /// returns a reference to the vector of stable cones (aka protocones)
239 const std::vector<PseudoJet> & stable_cones() const {return _protocones;}
240
241 /// an old name for getting the vector of stable cones (aka protocones)
242 const std::vector<PseudoJet> & protocones() const {return _protocones;}
243
244 /// return the # of the pass at which a given jet was found; will
245 /// return -1 if the pass is invalid
246 int pass(const PseudoJet & jet) const {return _pass[jet.cluster_hist_index()];}
247
248 /// return a brief summary of the contents of the extras object
249 /// (specifically, the number of protocones.
250 std::string description() const{
251 std::ostringstream ostr;
252 ostr << "This SISCone clustering found " << protocones().size()
253 << " stable protocones";
254 return ostr.str();
255 };
256
257 /// return the smallest difference in squared distance encountered
258 /// during splitting between a particle and two overlapping
259 /// protojets.
260 inline double most_ambiguous_split() const {return _most_ambiguous_split;}
261
262protected:
263 std::vector<PseudoJet> _protocones;
264 std::vector<int> _pass;
265 double _most_ambiguous_split;
266 const SISConeBasePlugin * _jet_def_plugin;
267};
268
269/// give the destructor its required implementation
270inline SISConeBaseExtras::~SISConeBaseExtras(){}
271
272//----------------------------------------------------------------------
273// implementation of the structure type associated with the UserScaleBase class
274
275/// \class SISConeBasePlugin::UserScaleBase::StructureType
276/// the structure that allows to store the information contained
277/// into a siscone::Cjet (built internally in SISCone from a stable
278/// cone) into a PseudoJet
280public:
281 /// base ctor (constructed from a ClusterSequence tin order to have
282 /// access to the initial particles
284 : _cs(cs){}
285
286 /// empty virtual dtor
287 virtual ~StructureType(){}
288
289 //--------------------------------------------------
290 // members inherited from the base class
291 /// the textual descripotion
292 virtual std::string description() const{
293 return "PseudoJet wrapping a siscone jet from a stable cone";
294 }
295
296 /// this structure has constituents
297 virtual bool has_constituents() const {return true;}
298
299 /// retrieve the constituents
300 ///
301 /// if you simply need to iterate over the constituents, it will be
302 /// faster to access them via constituent(i)
303 virtual std::vector<PseudoJet> constituents(const PseudoJet & /*reference*/) const{
304 std::vector<PseudoJet> constits;
305 constits.reserve(size());
306 for (unsigned int i=0; i<size();i++)
307 constits.push_back(constituent(i));
308 return constits;
309 }
310
311 //--------------------------------------------------
312 // additional information relevant for this structure
313
314 /// returns the number of constituents
315 virtual unsigned int size() const = 0;
316
317 /// returns the index (in the original particle list) of the ith
318 /// constituent
319 virtual int constituent_index(unsigned int i) const = 0;
320
321 /// returns the ith constituent (as a PseusoJet)
322 const PseudoJet & constituent(unsigned int i) const{
323 return _cs.jets()[constituent_index(i)];
324 }
325
326 // /// returns the scalar pt of this stable cone
327 // virtual double pt_tilde() const = 0;
328
329 /// returns the sm_var2 (signed ordering variable squared) for this stable cone
330 virtual double ordering_var2() const = 0;
331
332protected:
333 const ClusterSequence &_cs; ///< a reference to the CS (for access to the particles)
334};
335
336
337///@ingroup internal
338/// template class derived from UserScaleBase::StryctureType that
339/// works for both SISCone jet classes
340/// implemented below
341template<class Tjet>
343public:
344 UserScaleBaseStructureType(const Tjet &jet, const ClusterSequence &cs)
345 : UserScaleBase::StructureType(cs), _jet(jet){}
346
347 /// empty virtual dtor
349
350 //--------------------------------------------------
351 // additional information relevant for this structure
352
353 /// returns the number of constituents
354 virtual unsigned int size() const{
355 return _jet.n;
356 }
357
358 /// returns the index (in the original particle list) of the ith
359 /// constituent
360 virtual int constituent_index(unsigned int i) const{
361 return _jet.contents[i];
362 }
363
364 // /// returns the scalar pt of this stable cone
365 // virtual double pt_tilde() const{
366 // return _jet.pt_tilde;
367 // }
368
369 /// returns the sm_var2 (signed ordering variable squared) for this stable cone
370 virtual double ordering_var2() const{
371 return _jet.sm_var2;
372 }
373
374protected:
375 const Tjet &_jet; ///< a reference to the internal jet in SISCone
376};
377
378
379FASTJET_END_NAMESPACE // defined in fastjet/internal/base.hh
380
381#endif // __SISCONEBASEPLUGIN_HH__
382
base class to store extra information that plugins may provide
deals with clustering
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.
Definition: PseudoJet.hh:68
int cluster_hist_index() const
return the cluster_hist_index, intended to be used by clustering routines.
Definition: PseudoJet.hh:834
Class that provides extra information about a SISCone clustering.
int pass(const PseudoJet &jet) const
return the # of the pass at which a given jet was found; will return -1 if the pass is invalid
double most_ambiguous_split() const
return the smallest difference in squared distance encountered during splitting between a particle an...
const std::vector< PseudoJet > & stable_cones() const
returns a reference to the vector of stable cones (aka protocones)
std::string description() const
return a brief summary of the contents of the extras object (specifically, the number of protocones.
const std::vector< PseudoJet > & protocones() const
an old name for getting the vector of stable cones (aka protocones)
SISConeBaseExtras(int nparticles)
constructor
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 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)
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.
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