FastJet 3.4.1
SISConePlugin.hh
1#ifndef __SISCONEPLUGIN_HH__
2#define __SISCONEPLUGIN_HH__
3
4#include "SISConeBasePlugin.hh"
5
6// forward declaration of the siscone classes we'll need
7namespace siscone{
8 class Csiscone;
9 class Cjet;
10}
11
12
13FASTJET_BEGIN_NAMESPACE // defined in fastjet/internal/base.hh
14
15//----------------------------------------------------------------------
16//
17/// @ingroup plugins
18/// \class SISConePlugin
19/// Implementation of the SISCone algorithm (plugin for fastjet v2.1 upwards)
20///
21/// SISConePlugin is a plugin for fastjet (v2.1 upwards) that provides
22/// an interface to the seedless infrared safe cone jet finder by
23/// Gregory Soyez and Gavin Salam.
24///
25/// SISCone uses geometrical techniques to exhaustively consider all
26/// possible distinct cones. It then finds out which ones are stable
27/// and sends the result to the Tevatron Run-II type split-merge
28/// procedure for overlapping cones.
29///
30/// Four parameters govern the "physics" of the algorithm:
31///
32/// - the cone_radius (this should be self-explanatory!)
33///
34/// - the overlap_threshold is the parameter which dictates how much
35/// two jets must overlap (pt_overlap/min(pt1,pt2)) if they are to be
36/// merged
37///
38/// - Not all particles are in stable cones in the first round of
39/// searching for stable cones; one can therefore optionally have the
40/// the jet finder carry out additional passes of searching for
41/// stable cones among particles that were in no stable cone in
42/// previous passes --- the maximum number of passes carried out is
43/// n_pass_max. If this is zero then additional passes are carried
44/// out until no new stable cones are found.
45///
46/// - Protojet ptmin: protojets that are below this ptmin
47/// (default = 0) are discarded before each iteration of the
48/// split-merge loop.
49///
50/// One parameter governs some internal algorithmic shortcuts:
51///
52/// - if "caching" is turned on then the last event clustered by
53/// siscone is stored -- if the current event is identical and the
54/// cone_radius and n_pass_mass are identical, then the only part of
55/// the clustering that needs to be rerun is the split-merge part,
56/// leading to significant speed gains; there is a small (O(N) storage
57/// and speed) penalty for caching, so it should be kept off
58/// (default) if only a single overlap_threshold is used.
59///
60/// The final jets can be accessed by requestion the
61/// inclusive_jets(...) from the ClusterSequence object. Note that
62/// these PseudoJets have their user_index() set to the index of the
63/// pass in which they were found (first pass = 0). NB: This does not
64/// currently work for jets that consist of a single particle.
65///
66/// For further information on the details of the algorithm see the
67/// SISCone paper, arXiv:0704.0292 [JHEP 0705:086,2007].
68///
69/// For documentation about the implementation, see the
70/// siscone/doc/html/index.html file.
71//
73public:
74
75 /// enum for the different split-merge scale choices;
76 /// Note that order _must_ be the same as in siscone
78 SM_pt, ///< transverse momentum (E-scheme), IR unsafe
79 SM_Et, ///< transverse energy (E-scheme), not long. boost invariant
80 ///< original run-II choice [may not be implemented]
81 SM_mt, ///< transverse mass (E-scheme), IR safe except
82 ///< in decays of two identical narrow heavy particles
83 SM_pttilde ///< pt-scheme pt = \sum_{i in jet} |p_{ti}|, should
84 ///< be IR safe in all cases
85 };
86
87
88 /// Main constructor for the SISCone Plugin class.
89 ///
90 /// Note: wrt version prior to 2.4 this constructor differs in that a
91 /// the default value has been removed for overlap_threshold. The
92 /// former has been removed because the old default of 0.5 was found
93 /// to be unsuitable in high-noise environments; so the user should
94 /// now explicitly think about the value for this -- we recommend
95 /// 0.75.
96 ///
97 SISConePlugin (double cone_radius_in,
98 double overlap_threshold_in,
99 int n_pass_max_in = 0,
100 double protojet_ptmin_in = 0.0,
101 bool caching_in = false,
102 SplitMergeScale split_merge_scale_in = SM_pttilde,
103 double split_merge_stopping_scale_in = 0.0){
104 _cone_radius = cone_radius_in;
105 _overlap_threshold = overlap_threshold_in;
106 _n_pass_max = n_pass_max_in;
107 _protojet_ptmin = protojet_ptmin_in;
108 _caching = caching_in;
109 _split_merge_scale = split_merge_scale_in;
110 _split_merge_stopping_scale = split_merge_stopping_scale_in;
111 _ghost_sep_scale = 0.0;
112 _use_pt_weighted_splitting = false;
113 _user_scale = 0;}
114
115
116 /// Backwards compatible constructor for the SISCone Plugin class
117 SISConePlugin (double cone_radius_in,
118 double overlap_threshold_in,
119 int n_pass_max_in,
120 double protojet_ptmin_in,
121 bool caching_in,
122 bool split_merge_on_transverse_mass_in){
123 _cone_radius = cone_radius_in;
124 _overlap_threshold = overlap_threshold_in;
125 _n_pass_max = n_pass_max_in;
126 _protojet_ptmin = protojet_ptmin_in;
127 _caching = caching_in;
128 _split_merge_stopping_scale = 0.0;
129 _split_merge_scale = split_merge_on_transverse_mass_in ? SM_mt : SM_pttilde;
130 _ghost_sep_scale = 0.0;
131 _user_scale = 0;}
132
133 /// backwards compatible constructor for the SISCone Plugin class
134 /// (avoid using this in future).
135 SISConePlugin (double cone_radius_in,
136 double overlap_threshold_in,
137 int n_pass_max_in,
138 bool caching_in) {
139 _cone_radius = cone_radius_in;
140 _overlap_threshold = overlap_threshold_in;
141 _n_pass_max = n_pass_max_in;
142 _protojet_ptmin = 0.0;
143 _caching = caching_in;
144 _split_merge_scale = SM_mt;
145 _split_merge_stopping_scale = 0.0;
146 _ghost_sep_scale = 0.0;
147 _use_pt_weighted_splitting = false;
148 _user_scale = 0;}
149
150 /// minimum pt for a protojet to be considered in the split-merge step
151 /// of the algorithm
152 double protojet_ptmin () const {return _protojet_ptmin ;}
153
154 /// return the scale to be passed to SISCone as the protojet_ptmin
155 /// -- if we have a ghost separation scale that is above the
156 /// protojet_ptmin, then the ghost_separation_scale becomes the
157 /// relevant one to use here
158 double protojet_or_ghost_ptmin () const {return std::max(_protojet_ptmin,
159 _ghost_sep_scale);}
160
161 /// indicates scale used in split-merge
162 SplitMergeScale split_merge_scale() const {return _split_merge_scale;}
163 /// sets scale used in split-merge
164 void set_split_merge_scale(SplitMergeScale sms) {_split_merge_scale = sms;}
165
166 /// indicates whether the split-merge orders on transverse mass or not.
167 /// retained for backwards compatibility with 2.1.0b3
168 bool split_merge_on_transverse_mass() const {return _split_merge_scale == SM_mt ;}
169 void set_split_merge_on_transverse_mass(bool val) {
170 _split_merge_scale = val ? SM_mt : SM_pt;}
171
172 /// indicates whether the split-merge orders on transverse mass or not.
173 /// retained for backwards compatibility with 2.1.0b3
174 bool split_merge_use_pt_weighted_splitting() const {return _use_pt_weighted_splitting;}
175 void set_split_merge_use_pt_weighted_splitting(bool val) {
176 _use_pt_weighted_splitting = val;}
177
178 // the things that are required by base class
179 virtual std::string description () const;
180 virtual void run_clustering(ClusterSequence &) const ;
181
182protected:
183 virtual void reset_stored_plugin() const;
184
185private:
186 double _protojet_ptmin;
187 SplitMergeScale _split_merge_scale;
188
189 bool _use_pt_weighted_splitting;
190
191 // part needed for the cache
192 // variables for caching the results and the input
193 static SharedPtr<SISConePlugin > stored_plugin;
194 static SharedPtr<std::vector<PseudoJet> > stored_particles;
195 static SharedPtr<siscone::Csiscone > stored_siscone;
196};
197
198
199/////\class SISConePlugin::UserScaleBase::StructureType
200///// the structure that allows to store the information contained
201///// into a siscone::Cjet (built internally in SISCone from a stable
202///// cone) into a PseudoJet
203//class SISConePlugin::UserScaleBase::StructureType : public PseudoJetStructureBase {
204//public:
205// StructureType(const siscone::Cjet & jet, const ClusterSequence &cs)
206// : _jet(jet), _cs(cs){}
207//
208// //--------------------------------------------------
209// // members inherited from the base class
210// /// the textual descripotion
211// virtual std::string description() const;
212//
213// /// this structure has constituents
214// virtual bool has_constituents() const {return true;}
215//
216// /// retrieve the constituents
217// ///
218// /// if you simply need to iterate over the constituents, it will be
219// /// faster to access them via constituent(i)
220// virtual std::vector<PseudoJet> constituents(const PseudoJet & /*reference*/) const;
221//
222// //--------------------------------------------------
223// // additional information relevant for this structure
224//
225// /// returns the number of constituents
226// unsigned int size() const;
227//
228// /// returns the index (in the original particle list) of the ith
229// /// constituent
230// int constituent_index(unsigned int i) const;
231//
232// /// returns the ith constituent (as a PseusoJet)
233// const PseudoJet & constituent(unsigned int i) const;
234//
235// /// returns the scalar pt of this stable cone
236// double pt_tilde() const;
237//
238// /// returns the sm_var2 (signed ordering variable squared) for this stable cone
239// double ordering_var2() const;
240//
241//protected:
242// const siscone::Cjet &_jet; ///< a dreference to the internal jet in SISCone
243// const ClusterSequence &_cs; ///< a reference to the CS (for access to the particles)
244//};
245
246//======================================================================
247/// @ingroup extra_info
248/// \class SISConeExtras
249/// Class that provides extra information about a SISCone clustering
251public:
252 /// constructor
253 // it just initialises the pass information
254 SISConeExtras(int nparticles)
255 : SISConeBaseExtras(nparticles){}
256
257 /// access to the siscone jet def plugin (more convenient than
258 /// getting it from the original jet definition, because here it's
259 /// directly of the right type (rather than the base type)
261 return dynamic_cast<const SISConePlugin*>(_jet_def_plugin);
262 }
263
264private:
265 // let us be written to by SISConePlugin
266 friend class SISConePlugin;
267};
268
269FASTJET_END_NAMESPACE // defined in fastjet/internal/base.hh
270
271#endif // __SISCONEPLUGIN_HH__
272
Class that provides extra information about a SISCone clustering.
Class that provides extra information about a SISCone clustering.
SISConeExtras(int nparticles)
constructor
const SISConePlugin * jet_def_plugin() const
access to the siscone jet def plugin (more convenient than getting it from the original jet definitio...
Implementation of the SISCone algorithm (plugin for fastjet v2.1 upwards)
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 ...
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).
bool split_merge_use_pt_weighted_splitting() const
indicates whether the split-merge orders on transverse mass or not.
void set_split_merge_scale(SplitMergeScale sms)
sets scale used in split-merge
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=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.
double protojet_ptmin() const
minimum pt for a protojet to be considered in the split-merge step of the algorithm
SplitMergeScale split_merge_scale() const
indicates scale used in split-merge
bool split_merge_on_transverse_mass() const
indicates whether the split-merge orders on transverse mass or not.
SplitMergeScale
enum for the different split-merge scale choices; Note that order must be the same as in siscone
@ SM_pt
transverse momentum (E-scheme), IR unsafe