|
||||
File indexing completed on 2025-01-18 09:57:16
0001 #ifndef __FASTJET_TOOLS_PRUNER_HH__ 0002 #define __FASTJET_TOOLS_PRUNER_HH__ 0003 0004 //FJSTARTHEADER 0005 // $Id$ 0006 // 0007 // Copyright (c) 2005-2021, Matteo Cacciari, Gavin P. Salam and Gregory Soyez 0008 // 0009 //---------------------------------------------------------------------- 0010 // This file is part of FastJet. 0011 // 0012 // FastJet is free software; you can redistribute it and/or modify 0013 // it under the terms of the GNU General Public License as published by 0014 // the Free Software Foundation; either version 2 of the License, or 0015 // (at your option) any later version. 0016 // 0017 // The algorithms that underlie FastJet have required considerable 0018 // development. They are described in the original FastJet paper, 0019 // hep-ph/0512210 and in the manual, arXiv:1111.6097. If you use 0020 // FastJet as part of work towards a scientific publication, please 0021 // quote the version you use and include a citation to the manual and 0022 // optionally also to hep-ph/0512210. 0023 // 0024 // FastJet is distributed in the hope that it will be useful, 0025 // but WITHOUT ANY WARRANTY; without even the implied warranty of 0026 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the 0027 // GNU General Public License for more details. 0028 // 0029 // You should have received a copy of the GNU General Public License 0030 // along with FastJet. If not, see <http://www.gnu.org/licenses/>. 0031 //---------------------------------------------------------------------- 0032 //FJENDHEADER 0033 0034 #include "fastjet/ClusterSequence.hh" 0035 #include "fastjet/WrappedStructure.hh" 0036 #include "fastjet/tools/Transformer.hh" 0037 #include <iostream> 0038 #include <string> 0039 0040 FASTJET_BEGIN_NAMESPACE // defined in fastjet/internal/base.hh 0041 0042 // fwd declarations 0043 class Pruner; 0044 class PrunerStructure; 0045 class PruningRecombiner; 0046 class PruningPlugin; 0047 0048 // This tells third-party code that the pruner structure 0049 // stores Rcut info; the alternative is for the user to 0050 // get the information from the version number 0051 #define FASTJET_PRUNER_STRUCTURE_STORES_RCUT 0052 0053 //---------------------------------------------------------------------- 0054 /// @ingroup tools_generic 0055 /// \class Pruner 0056 /// Transformer that prunes a jet 0057 /// 0058 /// This transformer prunes a jet according to the ideas presented in 0059 /// arXiv:0903.5081 (S.D. Ellis, C.K. Vermilion and J.R. Walsh). 0060 /// 0061 /// The jet's constituents are reclustered with a user-specified jet 0062 /// definition, with the modification that objects i and j are only 0063 /// recombined if at least one of the following two criteria is 0064 /// satisfied: 0065 /// 0066 /// - the geometric distance between i and j is smaller than 'Rcut' 0067 /// with Rcut = Rcut_factor*2m/pt (Rcut_factor is a parameter of 0068 /// the Pruner and m and pt obtained from the jet being pruned) 0069 /// - the transverse momenta of i and j are at least 'zcut' p_t(i+j) 0070 /// 0071 /// If both these criteria fail, i and j are not recombined, the 0072 /// harder of i and j is kept, and the softer is rejected. 0073 /// 0074 /// Usage: 0075 /// \code 0076 /// Pruner pruner(jet_def, zcut, Rcut_factor); 0077 /// PseudoJet pruned_jet = pruner(jet); 0078 /// \endcode 0079 /// 0080 /// The pruned_jet has a valid associated cluster sequence. In addition 0081 /// the subjets of the original jet that have been vetoed by pruning 0082 /// (i.e. have been 'pruned away') can be accessed using 0083 /// 0084 /// \code 0085 /// vector<PseudoJet> rejected_subjets = pruned_jet.structure_of<Pruner>().rejected(); 0086 /// \endcode 0087 /// 0088 /// If the re-clustering happens to find more than a single inclusive 0089 /// jet (this should normally not happen if the radius of the jet 0090 /// definition used for the reclustering was set large enough), 0091 /// the hardest of these jets is retured as the result of the 0092 /// Pruner. The other jets can be accessed through 0093 /// 0094 /// \code 0095 /// vector<PseudoJet> extra_jets = pruned_jet.structure_of<Pruner>().extra_jets(); 0096 /// \endcode 0097 /// 0098 /// Instead of using Rcut_factor and zcut, one can alternatively 0099 /// construct a Pruner by passing two (pointers to) functions of 0100 /// PseudoJet that dynamically compute the Rcut and zcut to 0101 /// be used for the jet being pruned. 0102 /// 0103 /// When the jet being pruned has area support and explicit ghosts, 0104 /// the resulting pruned jet will likewise have area. 0105 /// 0106 //---------------------------------------------------------------------- 0107 class Pruner : public Transformer{ 0108 public: 0109 /// minimal constructor, which takes a jet algorithm, sets the radius 0110 /// to JetDefinition::max_allowable_R (practically equivalent to 0111 /// infinity) and also tries to use a recombiner based on the one in 0112 /// the jet definition of the particular jet being pruned. 0113 /// 0114 /// \param jet_alg the jet algorithm for the internal clustering 0115 /// \param zcut pt-fraction cut in the pruning 0116 /// \param Rcut_factor the angular distance cut in the pruning will be 0117 /// Rcut_factor * 2m/pt 0118 Pruner(const JetAlgorithm jet_alg, double zcut, double Rcut_factor) 0119 : _jet_def(jet_alg, JetDefinition::max_allowable_R), 0120 _zcut(zcut), _Rcut_factor(Rcut_factor), 0121 _zcut_dyn(0), _Rcut_dyn(0), _get_recombiner_from_jet(true) {} 0122 0123 0124 /// alternative ctor in which the full reclustering jet definition can 0125 /// be specified. 0126 /// 0127 /// \param jet_def the jet definition for the internal clustering 0128 /// \param zcut pt-fraction cut in the pruning 0129 /// \param Rcut_factor the angular distance cut in the pruning will be 0130 /// Rcut_factor * 2m/pt 0131 Pruner(const JetDefinition &jet_def, double zcut, double Rcut_factor) 0132 : _jet_def(jet_def), 0133 _zcut(zcut), _Rcut_factor(Rcut_factor), 0134 _zcut_dyn(0), _Rcut_dyn(0), _get_recombiner_from_jet(false) {} 0135 0136 0137 /// alternative ctor in which the pt-fraction cut and angular distance 0138 /// cut are functions of the jet being pruned. 0139 /// 0140 /// \param jet_def the jet definition for the internal clustering 0141 /// \param zcut_dyn dynamic pt-fraction cut in the pruning 0142 /// \param Rcut_dyn dynamic angular distance cut in the pruning 0143 Pruner(const JetDefinition &jet_def, 0144 const FunctionOfPseudoJet<double> *zcut_dyn, 0145 const FunctionOfPseudoJet<double> *Rcut_dyn); 0146 0147 /// action on a single jet 0148 virtual PseudoJet result(const PseudoJet &jet) const; 0149 0150 /// description 0151 virtual std::string description() const; 0152 0153 // the type of the associated structure 0154 typedef PrunerStructure StructureType; 0155 0156 private: 0157 /// check if the jet has explicit_ghosts (knowing that there is an 0158 /// area support) 0159 bool _check_explicit_ghosts(const PseudoJet &jet) const; 0160 0161 /// see if there is a common recombiner among the pieces; if there 0162 /// is return true and set jet_def_for_recombiner so that the 0163 /// recombiner can be taken from that JetDefinition. Otherwise, 0164 /// return false. 'assigned' is initially false; when true, each 0165 /// time we meet a new jet definition, we'll check it shares the 0166 /// same recombiner as jet_def_for_recombiner. 0167 bool _check_common_recombiner(const PseudoJet &jet, 0168 JetDefinition &jet_def_for_recombiner, 0169 bool assigned=false) const; 0170 0171 JetDefinition _jet_def; ///< the internal jet definition 0172 double _zcut; ///< the pt-fraction cut 0173 double _Rcut_factor; ///< the angular separation cut factor 0174 const FunctionOfPseudoJet<double> *_zcut_dyn; ///< dynamic zcut 0175 const FunctionOfPseudoJet<double> *_Rcut_dyn; ///< dynamic Rcut 0176 bool _get_recombiner_from_jet; ///< true for minimal constructor, 0177 ///< causes recombiner to be set equal 0178 ///< to that already used in the jet 0179 ///< (if it can be deduced) 0180 }; 0181 0182 0183 //---------------------------------------------------------------------- 0184 /// @ingroup tools_generic 0185 /// \class PrunerStructure 0186 /// The structure associated with a PseudoJet thas has gone through a 0187 /// Pruner transformer 0188 //---------------------------------------------------------------------- 0189 class PrunerStructure : public WrappedStructure{ 0190 public: 0191 /// default ctor 0192 /// \param result_jet the jet for which we have to keep the structure 0193 PrunerStructure(const PseudoJet & result_jet) 0194 : WrappedStructure(result_jet.structure_shared_ptr()){} 0195 0196 /// description 0197 virtual std::string description() const{ return "Pruned PseudoJet";} 0198 0199 /// return the constituents that have been rejected 0200 std::vector<PseudoJet> rejected() const{ 0201 return validated_cs()->childless_pseudojets(); 0202 } 0203 0204 /// return the other jets that may have been found along with the 0205 /// result of the pruning 0206 /// The resulting vector is sorted in pt 0207 std::vector<PseudoJet> extra_jets() const; 0208 0209 /// return the value of Rcut that was used for this specific pruning. 0210 double Rcut() const {return _Rcut;} 0211 0212 /// return the value of Rcut that was used for this specific pruning. 0213 double zcut() const {return _zcut;} 0214 0215 protected: 0216 friend class Pruner; ///< to allow setting the internal information 0217 0218 private: 0219 double _Rcut, _zcut; 0220 }; 0221 0222 //---------------------------------------------------------------------- 0223 /// \if internal_doc 0224 /// @ingroup internal 0225 /// \class PruningRecombiner 0226 /// recombines the objects that are not vetoed by pruning 0227 /// 0228 /// This recombiner only recombines, using the provided 'recombiner', 0229 /// objects (i and j) that pass at least one of the following two criteria: 0230 /// 0231 /// - the geometric distance between i and j is smaller than 'Rcut' 0232 /// - the transverse momenta of i and j are at least 'zcut' p_t(i+j) 0233 /// 0234 /// If both these criteria fail, the hardest of i and j is kept and 0235 /// the softest is rejected. 0236 /// 0237 /// Note that this in not meant for standalone use [in particular 0238 /// because it could lead to memory (and thread-safety) issues due to 0239 /// the rejected indices stored internally (as a mutable variable)]. 0240 /// 0241 /// \endif 0242 class PruningRecombiner : public JetDefinition::Recombiner{ 0243 public: 0244 /// ctor 0245 /// \param zcut transverse momentum fraction cut 0246 /// \param Rcut angular separation cut 0247 /// \param recomb pointer to a recombiner to use to cluster pairs 0248 PruningRecombiner(double zcut, double Rcut, 0249 const JetDefinition::Recombiner *recombiner) 0250 : _zcut2(zcut*zcut), _Rcut2(Rcut*Rcut), 0251 _recombiner(recombiner){} 0252 0253 /// perform a recombination taking into account the pruning 0254 /// conditions 0255 virtual void recombine(const PseudoJet &pa, 0256 const PseudoJet &pb, 0257 PseudoJet &pab) const; 0258 0259 /// returns the description of the recombiner 0260 virtual std::string description() const; 0261 0262 /// return the history indices that have been pruned away 0263 const std::vector<unsigned int> & rejected() const{ return _rejected;} 0264 0265 /// clears the list of rejected indices 0266 /// 0267 /// If one decides to use this recombiner standalone, one has to 0268 /// call this after each clustering in order for the rejected() vector 0269 /// to remain sensible and not grow to infinite size. 0270 void clear_rejected(){ _rejected.clear();} 0271 0272 private: 0273 double _zcut2; ///< transverse momentum fraction cut (squared) 0274 double _Rcut2; ///< angular separation cut (squared) 0275 const JetDefinition::Recombiner *_recombiner; ///< the underlying recombiner to use 0276 mutable std::vector<unsigned int> _rejected; ///< list of rejected history indices 0277 }; 0278 0279 0280 //---------------------------------------------------------------------- 0281 /// \if internal_doc 0282 /// @ingroup internal 0283 /// \class PruningPlugin 0284 /// FastJet internal plugin that clusters the particles using the 0285 /// PruningRecombiner. 0286 /// 0287 /// See PruningRecombiner for a description of what pruning does. 0288 /// 0289 /// Note that this is an internal FastJet class used by the Pruner 0290 /// transformer and it is not meant to be used as a standalone clustering 0291 /// tool. 0292 /// 0293 /// \endif 0294 //---------------------------------------------------------------------- 0295 class PruningPlugin : public JetDefinition::Plugin{ 0296 public: 0297 /// ctor 0298 /// \param jet_def the jet definition to be used for the 0299 /// internal clustering 0300 /// \param zcut transverse momentum fraction cut 0301 /// \param Rcut angular separation cut 0302 PruningPlugin(const JetDefinition &jet_def, double zcut, double Rcut) 0303 : _jet_def(jet_def), _zcut(zcut), _Rcut(Rcut){} 0304 0305 /// the actual clustering work for the plugin 0306 virtual void run_clustering(ClusterSequence &input_cs) const; 0307 0308 /// description of the plugin 0309 virtual std::string description() const; 0310 0311 /// returns the radius 0312 virtual double R() const {return _jet_def.R();} 0313 0314 private: 0315 JetDefinition _jet_def; ///< the internal jet definition 0316 double _zcut; ///< transverse momentum fraction cut 0317 double _Rcut; ///< angular separation cut 0318 }; 0319 0320 0321 0322 FASTJET_END_NAMESPACE // defined in fastjet/internal/base.hh 0323 0324 #endif // __FASTJET_TOOLS_PRUNER_HH__
[ Source navigation ] | [ Diff markup ] | [ Identifier search ] | [ general search ] |
This page was automatically generated by the 2.3.7 LXR engine. The LXR team |