source: HiSusy/trunk/Delphes-3.0.0/external/fastjet/plugins/SISCone/fastjet/SISConePlugin.hh @ 1

Last change on this file since 1 was 1, checked in by zerwas, 11 years ago

first import of structure, PYTHIA8 and DELPHES

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