FastJet  3.2.2
Pruner.hh
1 #ifndef __FASTJET_TOOLS_PRUNER_HH__
2 #define __FASTJET_TOOLS_PRUNER_HH__
3 
4 //FJSTARTHEADER
5 // $Id: Pruner.hh 3481 2014-07-29 17:24:12Z soyez $
6 //
7 // Copyright (c) 2005-2014, Matteo Cacciari, Gavin P. Salam and Gregory Soyez
8 //
9 //----------------------------------------------------------------------
10 // This file is part of FastJet.
11 //
12 // FastJet is free software; you can redistribute it and/or modify
13 // it under the terms of the GNU General Public License as published by
14 // the Free Software Foundation; either version 2 of the License, or
15 // (at your option) any later version.
16 //
17 // The algorithms that underlie FastJet have required considerable
18 // development. They are described in the original FastJet paper,
19 // hep-ph/0512210 and in the manual, arXiv:1111.6097. If you use
20 // FastJet as part of work towards a scientific publication, please
21 // quote the version you use and include a citation to the manual and
22 // optionally also to hep-ph/0512210.
23 //
24 // FastJet is distributed in the hope that it will be useful,
25 // but WITHOUT ANY WARRANTY; without even the implied warranty of
26 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
27 // GNU General Public License for more details.
28 //
29 // You should have received a copy of the GNU General Public License
30 // along with FastJet. If not, see <http://www.gnu.org/licenses/>.
31 //----------------------------------------------------------------------
32 //FJENDHEADER
33 
34 #include "fastjet/ClusterSequence.hh"
35 #include "fastjet/WrappedStructure.hh"
36 #include "fastjet/tools/Transformer.hh"
37 #include <iostream>
38 #include <string>
39 
40 FASTJET_BEGIN_NAMESPACE // defined in fastjet/internal/base.hh
41 
42 // fwd declarations
43 class Pruner;
44 class PrunerStructure;
45 class PruningRecombiner;
46 class PruningPlugin;
47 
48 // This tells third-party code that the pruner structure
49 // stores Rcut info; the alternative is for the user to
50 // get the information from the version number
51 #define FASTJET_PRUNER_STRUCTURE_STORES_RCUT
52 
53 //----------------------------------------------------------------------
54 /// @ingroup tools_generic
55 /// \class Pruner
56 /// Transformer that prunes a jet
57 ///
58 /// This transformer prunes a jet according to the ideas presented in
59 /// arXiv:0903.5081 (S.D. Ellis, C.K. Vermilion and J.R. Walsh).
60 ///
61 /// The jet's constituents are reclustered with a user-specified jet
62 /// definition, with the modification that objects i and j are only
63 /// recombined if at least one of the following two criteria is
64 /// satisfied:
65 ///
66 /// - the geometric distance between i and j is smaller than 'Rcut'
67 /// with Rcut = Rcut_factor*2m/pt (Rcut_factor is a parameter of
68 /// the Pruner and m and pt obtained from the jet being pruned)
69 /// - the transverse momenta of i and j are at least 'zcut' p_t(i+j)
70 ///
71 /// If both these criteria fail, i and j are not recombined, the
72 /// harder of i and j is kept, and the softer is rejected.
73 ///
74 /// Usage:
75 /// \code
76 /// Pruner pruner(jet_def, zcut, Rcut_factor);
77 /// PseudoJet pruned_jet = pruner(jet);
78 /// \endcode
79 ///
80 /// The pruned_jet has a valid associated cluster sequence. In addition
81 /// the subjets of the original jet that have been vetoed by pruning
82 /// (i.e. have been 'pruned away') can be accessed using
83 ///
84 /// \code
85 /// vector<PseudoJet> rejected_subjets = pruned_jet.structure_of<Pruner>().rejected();
86 /// \endcode
87 ///
88 /// If the re-clustering happens to find more than a single inclusive
89 /// jet (this should normally not happen if the radius of the jet
90 /// definition used for the reclustering was set large enough),
91 /// the hardest of these jets is retured as the result of the
92 /// Pruner. The other jets can be accessed through
93 ///
94 /// \code
95 /// vector<PseudoJet> extra_jets = pruned_jet.structure_of<Pruner>().extra_jets();
96 /// \endcode
97 ///
98 /// Instead of using Rcut_factor and zcut, one can alternatively
99 /// construct a Pruner by passing two (pointers to) functions of
100 /// PseudoJet that dynamically compute the Rcut and zcut to
101 /// be used for the jet being pruned.
102 ///
103 /// When the jet being pruned has area support and explicit ghosts,
104 /// the resulting pruned jet will likewise have area.
105 ///
106 //----------------------------------------------------------------------
107 class Pruner : public Transformer{
108 public:
109  /// minimal constructor, which takes a jet algorithm, sets the radius
110  /// to JetDefinition::max_allowable_R (practically equivalent to
111  /// infinity) and also tries to use a recombiner based on the one in
112  /// the jet definition of the particular jet being pruned.
113  ///
114  /// \param jet_alg the jet algorithm for the internal clustering
115  /// \param zcut pt-fraction cut in the pruning
116  /// \param Rcut_factor the angular distance cut in the pruning will be
117  /// Rcut_factor * 2m/pt
118  Pruner(const JetAlgorithm jet_alg, double zcut, double Rcut_factor)
119  : _jet_def(jet_alg, JetDefinition::max_allowable_R),
120  _zcut(zcut), _Rcut_factor(Rcut_factor),
121  _zcut_dyn(0), _Rcut_dyn(0), _get_recombiner_from_jet(true) {}
122 
123 
124  /// alternative ctor in which the full reclustering jet definition can
125  /// be specified.
126  ///
127  /// \param jet_def the jet definition for the internal clustering
128  /// \param zcut pt-fraction cut in the pruning
129  /// \param Rcut_factor the angular distance cut in the pruning will be
130  /// Rcut_factor * 2m/pt
131  Pruner(const JetDefinition &jet_def, double zcut, double Rcut_factor)
132  : _jet_def(jet_def),
133  _zcut(zcut), _Rcut_factor(Rcut_factor),
134  _zcut_dyn(0), _Rcut_dyn(0), _get_recombiner_from_jet(false) {}
135 
136 
137  /// alternative ctor in which the pt-fraction cut and angular distance
138  /// cut are functions of the jet being pruned.
139  ///
140  /// \param jet_def the jet definition for the internal clustering
141  /// \param zcut_dyn dynamic pt-fraction cut in the pruning
142  /// \param Rcut_dyn dynamic angular distance cut in the pruning
143  Pruner(const JetDefinition &jet_def,
144  const FunctionOfPseudoJet<double> *zcut_dyn,
145  const FunctionOfPseudoJet<double> *Rcut_dyn);
146 
147  /// action on a single jet
148  virtual PseudoJet result(const PseudoJet &jet) const;
149 
150  /// description
151  virtual std::string description() const;
152 
153  // the type of the associated structure
155 
156 private:
157  /// check if the jet has explicit_ghosts (knowing that there is an
158  /// area support)
159  bool _check_explicit_ghosts(const PseudoJet &jet) const;
160 
161  /// see if there is a common recombiner among the pieces; if there
162  /// is return true and set jet_def_for_recombiner so that the
163  /// recombiner can be taken from that JetDefinition. Otherwise,
164  /// return false. 'assigned' is initially false; when true, each
165  /// time we meet a new jet definition, we'll check it shares the
166  /// same recombiner as jet_def_for_recombiner.
167  bool _check_common_recombiner(const PseudoJet &jet,
168  JetDefinition &jet_def_for_recombiner,
169  bool assigned=false) const;
170 
171  JetDefinition _jet_def; ///< the internal jet definition
172  double _zcut; ///< the pt-fraction cut
173  double _Rcut_factor; ///< the angular separation cut factor
174  const FunctionOfPseudoJet<double> *_zcut_dyn; ///< dynamic zcut
175  const FunctionOfPseudoJet<double> *_Rcut_dyn; ///< dynamic Rcut
176  bool _get_recombiner_from_jet; ///< true for minimal constructor,
177  ///< causes recombiner to be set equal
178  ///< to that already used in the jet
179  ///< (if it can be deduced)
180 };
181 
182 
183 //----------------------------------------------------------------------
184 /// @ingroup tools_generic
185 /// \class PrunerStructure
186 /// The structure associated with a PseudoJet thas has gone through a
187 /// Pruner transformer
188 //----------------------------------------------------------------------
190 public:
191  /// default ctor
192  /// \param result_jet the jet for which we have to keep the structure
193  PrunerStructure(const PseudoJet & result_jet)
194  : WrappedStructure(result_jet.structure_shared_ptr()){}
195 
196  /// description
197  virtual std::string description() const{ return "Pruned PseudoJet";}
198 
199  /// return the constituents that have been rejected
200  std::vector<PseudoJet> rejected() const{
201  return validated_cs()->childless_pseudojets();
202  }
203 
204  /// return the other jets that may have been found along with the
205  /// result of the pruning
206  /// The resulting vector is sorted in pt
207  std::vector<PseudoJet> extra_jets() const;
208 
209  /// return the value of Rcut that was used for this specific pruning.
210  double Rcut() const {return _Rcut;}
211 
212  /// return the value of Rcut that was used for this specific pruning.
213  double zcut() const {return _zcut;}
214 
215 protected:
216  friend class Pruner; ///< to allow setting the internal information
217 
218 private:
219  double _Rcut, _zcut;
220 };
221 
222 //----------------------------------------------------------------------
223 /// \if internal_doc
224 /// @ingroup internal
225 /// \class PruningRecombiner
226 /// recombines the objects that are not vetoed by pruning
227 ///
228 /// This recombiner only recombines, using the provided 'recombiner',
229 /// objects (i and j) that pass at least one of the following two criteria:
230 ///
231 /// - the geometric distance between i and j is smaller than 'Rcut'
232 /// - the transverse momenta of i and j are at least 'zcut' p_t(i+j)
233 ///
234 /// If both these criteria fail, the hardest of i and j is kept and
235 /// the softest is rejected.
236 ///
237 /// Note that this in not meant for standalone use [in particular
238 /// because it could lead to memory issues due to the rejected indices
239 /// stored internally].
240 ///
241 /// \endif
242 class PruningRecombiner : public JetDefinition::Recombiner{
243 public:
244  /// ctor
245  /// \param zcut transverse momentum fraction cut
246  /// \param Rcut angular separation cut
247  /// \param recomb pointer to a recombiner to use to cluster pairs
248  PruningRecombiner(double zcut, double Rcut,
249  const JetDefinition::Recombiner *recombiner)
250  : _zcut2(zcut*zcut), _Rcut2(Rcut*Rcut),
251  _recombiner(recombiner){}
252 
253  /// perform a recombination taking into account the pruning
254  /// conditions
255  virtual void recombine(const PseudoJet &pa,
256  const PseudoJet &pb,
257  PseudoJet &pab) const;
258 
259  /// returns the description of the recombiner
260  virtual std::string description() const;
261 
262  /// return the history indices that have been pruned away
263  const std::vector<unsigned int> & rejected() const{ return _rejected;}
264 
265  /// clears the list of rejected indices
266  ///
267  /// If one decides to use this recombiner standalone, one has to
268  /// call this after each clustering in order for the rejected() vector
269  /// to remain sensible and not grow to infinite size.
270  void clear_rejected(){ _rejected.clear();}
271 
272 private:
273  double _zcut2; ///< transverse momentum fraction cut (squared)
274  double _Rcut2; ///< angular separation cut (squared)
275  const JetDefinition::Recombiner *_recombiner; ///< the underlying recombiner to use
276  mutable std::vector<unsigned int> _rejected; ///< list of rejected history indices
277 };
278 
279 
280 //----------------------------------------------------------------------
281 /// \if internal_doc
282 /// @ingroup internal
283 /// \class PruningPlugin
284 /// FastJet internal plugin that clusters the particles using the
285 /// PruningRecombiner.
286 ///
287 /// See PruningRecombiner for a description of what pruning does.
288 ///
289 /// Note that this is an internal FastJet class used by the Pruner
290 /// transformer and it is not meant to be used as a standalone clustering
291 /// tool.
292 ///
293 /// \endif
294 //----------------------------------------------------------------------
295 class PruningPlugin : public JetDefinition::Plugin{
296 public:
297  /// ctor
298  /// \param jet_def the jet definition to be used for the
299  /// internal clustering
300  /// \param zcut transverse momentum fraction cut
301  /// \param Rcut angular separation cut
302  PruningPlugin(const JetDefinition &jet_def, double zcut, double Rcut)
303  : _jet_def(jet_def), _zcut(zcut), _Rcut(Rcut){}
304 
305  /// the actual clustering work for the plugin
306  virtual void run_clustering(ClusterSequence &input_cs) const;
307 
308  /// description of the plugin
309  virtual std::string description() const;
310 
311  /// returns the radius
312  virtual double R() const {return _jet_def.R();}
313 
314 private:
315  JetDefinition _jet_def; ///< the internal jet definition
316  double _zcut; ///< transverse momentum fraction cut
317  double _Rcut; ///< angular separation cut
318 };
319 
320 
321 
322 FASTJET_END_NAMESPACE // defined in fastjet/internal/base.hh
323 
324 #endif // __FASTJET_TOOLS_PRUNER_HH__
Pruner(const JetDefinition &jet_def, double zcut, double Rcut_factor)
alternative ctor in which the full reclustering jet definition can be specified.
Definition: Pruner.hh:131
deals with clustering
double zcut() const
return the value of Rcut that was used for this specific pruning.
Definition: Pruner.hh:213
This wraps a (shared) pointer to an underlying structure.
Base (abstract) class for a jet transformer.
Definition: Transformer.hh:71
double Rcut() const
return the value of Rcut that was used for this specific pruning.
Definition: Pruner.hh:210
Transformer that prunes a jet.
Definition: Pruner.hh:107
PrunerStructure(const PseudoJet &result_jet)
default ctor
Definition: Pruner.hh:193
The structure associated with a PseudoJet thas has gone through a Pruner transformer.
Definition: Pruner.hh:189
virtual std::string description() const
description
Definition: Pruner.hh:197
Pruner(const JetAlgorithm jet_alg, double zcut, double Rcut_factor)
minimal constructor, which takes a jet algorithm, sets the radius to JetDefinition::max_allowable_R (...
Definition: Pruner.hh:118
std::vector< PseudoJet > rejected() const
return the constituents that have been rejected
Definition: Pruner.hh:200
a class that allows a user to introduce their own "plugin" jet finder
An abstract base class that will provide the recombination scheme facilities and/or allow a user to e...
JetAlgorithm
the various families of jet-clustering algorithm
Class to contain pseudojets, including minimal information of use to jet-clustering routines...
Definition: PseudoJet.hh:67
class that is intended to hold a full definition of the jet clusterer