|
||||
File indexing completed on 2025-01-18 09:57:16
0001 //FJSTARTHEADER 0002 // $Id$ 0003 // 0004 // Copyright (c) 2005-2021, Matteo Cacciari, Gavin P. Salam and Gregory Soyez 0005 // 0006 //---------------------------------------------------------------------- 0007 // This file is part of FastJet. 0008 // 0009 // FastJet is free software; you can redistribute it and/or modify 0010 // it under the terms of the GNU General Public License as published by 0011 // the Free Software Foundation; either version 2 of the License, or 0012 // (at your option) any later version. 0013 // 0014 // The algorithms that underlie FastJet have required considerable 0015 // development. They are described in the original FastJet paper, 0016 // hep-ph/0512210 and in the manual, arXiv:1111.6097. If you use 0017 // FastJet as part of work towards a scientific publication, please 0018 // quote the version you use and include a citation to the manual and 0019 // optionally also to hep-ph/0512210. 0020 // 0021 // FastJet is distributed in the hope that it will be useful, 0022 // but WITHOUT ANY WARRANTY; without even the implied warranty of 0023 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the 0024 // GNU General Public License for more details. 0025 // 0026 // You should have received a copy of the GNU General Public License 0027 // along with FastJet. If not, see <http://www.gnu.org/licenses/>. 0028 //---------------------------------------------------------------------- 0029 //FJENDHEADER 0030 0031 #ifndef __CASUBJET_TAGGER_HH__ 0032 #define __CASUBJET_TAGGER_HH__ 0033 0034 #include "fastjet/PseudoJet.hh" 0035 #include "fastjet/WrappedStructure.hh" 0036 #include "fastjet/tools/Transformer.hh" 0037 #include "fastjet/LimitedWarning.hh" 0038 0039 FASTJET_BEGIN_NAMESPACE 0040 0041 class CASubJetTagger; 0042 class CASubJetTaggerStructure; 0043 0044 //---------------------------------------------------------------------- 0045 /// @ingroup tools_taggers 0046 /// \class CASubJetTagger 0047 /// clean (almost parameter-free) tagger searching for the element in 0048 /// the clustering history that maximises a chosen distance 0049 /// 0050 /// class to help us get a clean (almost parameter-free) handle on 0051 /// substructure inside a C/A jet. It follows the logic described in 0052 /// arXiv:0906.0728 (and is inspired by the original Cambridge 0053 /// algorithm paper in its use of separate angular and dimensionful 0054 /// distances), but provides some extra flexibility. 0055 /// 0056 /// It searches for all splittings that pass a symmetry cut (zcut) and 0057 /// then selects the one with the largest auxiliary scale choice 0058 /// (e.g. jade distance of the splitting, kt distance of the 0059 /// splitting, etc.) 0060 /// 0061 /// By default, the zcut is calculated from the fraction of the child 0062 /// pt carried by the parent jet. If one calls set_absolute_z_cut the 0063 /// fraction of transverse momentum will be computed wrt the original 0064 /// jet. 0065 /// 0066 /// original code copyright (C) 2009 by Gavin Salam, released under 0067 /// the GPL. 0068 /// 0069 /// \section desc Options 0070 /// 0071 /// - the distance choice: options are 0072 /// kt2_distance : usual min(kti^2,ktj^2)DeltaR_{ij}^2 0073 /// jade_distance : kti . ktj DeltaR_{ij}^2 (LI version of jade) 0074 /// jade2_distance : kti . ktj DeltaR_{ij}^4 (LI version of jade * DR^2) 0075 /// plain_distance : DeltaR_{ij}^2 0076 /// mass_drop_distance : m_jet - max(m_parent1,m_parent2) 0077 /// dot_product_distance: parent1.parent2 0078 /// (kt2_distance by default) 0079 /// 0080 /// - the z cut (0 by default) 0081 /// 0082 /// - by calling set_absolute_z_cut(), one can ask that the pt 0083 /// fraction if calculated wrt the original jet 0084 /// 0085 /// - by calling set_dr_min(drmin), one can ask that only the 0086 /// recombinations where the 2 objects are (geometrically) distant 0087 /// by at least drmin are kept in the maximisation. 0088 /// 0089 /// \section input Input conditions 0090 /// 0091 /// - the jet must have been obtained from a Cambridge/Aachen cluster 0092 /// sequence 0093 /// 0094 /// \section output Output/structure 0095 /// 0096 /// - the element of the cluster sequence maximising the requested 0097 /// distance (and satisfying the zcut) is returned. 0098 /// 0099 /// - if the original jet has no parents, it will be returned 0100 /// 0101 /// - the value of the "z" and distance corresponding to that history 0102 /// element are stored and accessible through 0103 /// result.structure_of<CASubJetTagger>().z(); 0104 /// result.structure_of<CASubJetTagger>().distance(); 0105 /// 0106 class CASubJetTagger : public Transformer { 0107 public: 0108 0109 /// the available choices of auxiliary scale with respect to which 0110 /// to order the splittings 0111 enum ScaleChoice { 0112 kt2_distance, //< usual min(kti^2,ktj^2)DeltaR_{ij}^2 0113 jade_distance, //< kti . ktj DeltaR_{ij}^2 (LI version of jade) 0114 jade2_distance, //< kti . ktj DeltaR_{ij}^4 (LI version of jade * DR^2) 0115 plain_distance, //< DeltaR_{ij}^2 0116 mass_drop_distance, //< m_jet - max(m_parent1,m_parent2) 0117 dot_product_distance //< parent1.parent2 0118 }; 0119 0120 /// just constructs 0121 CASubJetTagger(ScaleChoice scale_choice = jade_distance, 0122 double z_threshold = 0.1) 0123 : _scale_choice(scale_choice), _z_threshold(z_threshold), 0124 _dr2_min(0.0), _absolute_z_cut(false){}; 0125 0126 /// sets a minimum delta R below which spliting will be ignored 0127 /// (only relevant if set prior to calling run()) 0128 void set_dr_min(double drmin) {_dr2_min = drmin*drmin;} 0129 0130 /// returns a textual description of the tagger 0131 virtual std::string description() const; 0132 0133 /// If (abs_z_cut) is set to false (the default) then for a 0134 /// splitting to be considered, each subjet must satisfy 0135 /// 0136 /// p_{t,sub} > z_threshold * p_{t,parent} 0137 /// 0138 /// whereas if it is set to true, then each subject must satisfy 0139 /// 0140 /// p_{t,sub} > z_threshold * p_{t,original-jet} 0141 /// 0142 /// where parent is the immediate parent of the splitting, and 0143 /// original jet is the one supplied to the run() function. 0144 /// 0145 /// Only relevant is called prior to run(). 0146 void set_absolute_z_cut(bool abs_z_cut=true) {_absolute_z_cut = abs_z_cut;} 0147 0148 /// runs the tagger on the given jet and 0149 /// returns the tagged PseudoJet if successful, or a PseudoJet==0 otherwise 0150 /// (standard access is through operator()). 0151 virtual PseudoJet result(const PseudoJet & jet) const; 0152 0153 /// the type of Structure returned 0154 typedef CASubJetTaggerStructure StructureType; 0155 0156 protected: 0157 /// class that contains the result internally 0158 class JetAux { 0159 public: 0160 PseudoJet jet; //< the subjet (immediate parent of splitting) 0161 double aux_distance; //< the auxiliary distance between its two subjets 0162 double delta_r; //< the angular distance between its two subjets 0163 double z; //< the transverse momentum fraction 0164 }; 0165 0166 0167 void _recurse_through_jet(const PseudoJet & current_jet, 0168 JetAux &aux_max, 0169 const PseudoJet & original_jet) const; 0170 0171 ScaleChoice _scale_choice; 0172 double _z_threshold; 0173 double _dr2_min; 0174 bool _absolute_z_cut; 0175 0176 static LimitedWarning _non_ca_warnings; 0177 }; 0178 0179 0180 //------------------------------------------------------------------------ 0181 /// @ingroup tools_taggers 0182 /// the structure returned by a CASubJetTagger 0183 /// 0184 /// Since this is directly an element of the ClusterSequence, we keep 0185 /// basically the original ClusterSequenceStructure (wrapped for 0186 /// memory-management reasons) and add information about the pt 0187 /// fraction and distance of the subjet structure 0188 class CASubJetTaggerStructure : public WrappedStructure{ 0189 public: 0190 /// default ctor 0191 /// \param result_jet the jet for which we have to keep the 0192 /// structure 0193 CASubJetTaggerStructure(const PseudoJet & result_jet) 0194 : WrappedStructure(result_jet.structure_shared_ptr()){} 0195 0196 /// returns the scale choice asked for the maximisation 0197 CASubJetTagger::ScaleChoice scale_choice() const {return _scale_choice;} 0198 0199 /// returns the value of the distance measure (corresponding to 0200 /// ScaleChoice) for this jet's splitting 0201 double distance() const {return _distance;} 0202 0203 /// returns the pt fraction contained by the softer of the two component 0204 /// pieces of this jet (normalised relative to this jet) 0205 double z() const {return _z;} 0206 0207 /// returns the pt fraction contained by the softer of the two component 0208 /// pieces of this jet (normalised relative to the original jet) 0209 bool absolute_z() const {return _absolute_z;} 0210 0211 // /// returns the original jet (before tagging) 0212 // const PseudoJet & original() const {return _original_jet;} 0213 0214 protected: 0215 CASubJetTagger::ScaleChoice _scale_choice; ///< the user scale choice 0216 double _distance; ///< the maximal distance associated with the result 0217 bool _absolute_z; ///< whether z is computed wrt to the original jet or not 0218 double _z; ///< the transverse momentum fraction 0219 // PseudoJet _original_jet; ///< the original jet (before tagging) 0220 0221 friend class CASubJetTagger; ///< to allow setting the internal information 0222 }; 0223 0224 FASTJET_END_NAMESPACE 0225 0226 #endif // __CASUBJET_HH__
[ Source navigation ] | [ Diff markup ] | [ Identifier search ] | [ general search ] |
This page was automatically generated by the 2.3.7 LXR engine. The LXR team |