Fork me on GitHub

source: git/external/fastjet/tools/Pruner.hh@ e2a76ae

ImprovedOutputFile Timing dual_readout llp
Last change on this file since e2a76ae was d7d2da3, checked in by pavel <pavel@…>, 12 years ago

move branches/ModularDelphes to trunk

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