1 | //FJSTARTHEADER
|
---|
2 | // $Id: CASubJetTagger.hh 4442 2020-05-05 07:50:11Z soyez $
|
---|
3 | //
|
---|
4 | // Copyright (c) 2005-2020, Matteo Cacciari, Gavin P. Salam and Gregory Soyez
|
---|
5 | //
|
---|
6 | //----------------------------------------------------------------------
|
---|
7 | // This file is part of FastJet.
|
---|
8 | //
|
---|
9 | // FastJet is free software; you can redistribute it and/or modify
|
---|
10 | // it under the terms of the GNU General Public License as published by
|
---|
11 | // the Free Software Foundation; either version 2 of the License, or
|
---|
12 | // (at your option) any later version.
|
---|
13 | //
|
---|
14 | // The algorithms that underlie FastJet have required considerable
|
---|
15 | // development. They are described in the original FastJet paper,
|
---|
16 | // hep-ph/0512210 and in the manual, arXiv:1111.6097. If you use
|
---|
17 | // FastJet as part of work towards a scientific publication, please
|
---|
18 | // quote the version you use and include a citation to the manual and
|
---|
19 | // optionally also to hep-ph/0512210.
|
---|
20 | //
|
---|
21 | // FastJet is distributed in the hope that it will be useful,
|
---|
22 | // but WITHOUT ANY WARRANTY; without even the implied warranty of
|
---|
23 | // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
|
---|
24 | // GNU General Public License for more details.
|
---|
25 | //
|
---|
26 | // You should have received a copy of the GNU General Public License
|
---|
27 | // along with FastJet. If not, see <http://www.gnu.org/licenses/>.
|
---|
28 | //----------------------------------------------------------------------
|
---|
29 | //FJENDHEADER
|
---|
30 |
|
---|
31 | #ifndef __CASUBJET_TAGGER_HH__
|
---|
32 | #define __CASUBJET_TAGGER_HH__
|
---|
33 |
|
---|
34 | #include <fastjet/PseudoJet.hh>
|
---|
35 | #include <fastjet/WrappedStructure.hh>
|
---|
36 | #include <fastjet/tools/Transformer.hh>
|
---|
37 | #include "fastjet/LimitedWarning.hh"
|
---|
38 |
|
---|
39 | FASTJET_BEGIN_NAMESPACE
|
---|
40 |
|
---|
41 | class CASubJetTagger;
|
---|
42 | class CASubJetTaggerStructure;
|
---|
43 |
|
---|
44 | //----------------------------------------------------------------------
|
---|
45 | /// @ingroup tools_taggers
|
---|
46 | /// \class CASubJetTagger
|
---|
47 | /// clean (almost parameter-free) tagger searching for the element in
|
---|
48 | /// the clustering history that maximises a chosen distance
|
---|
49 | ///
|
---|
50 | /// class to help us get a clean (almost parameter-free) handle on
|
---|
51 | /// substructure inside a C/A jet. It follows the logic described in
|
---|
52 | /// arXiv:0906.0728 (and is inspired by the original Cambridge
|
---|
53 | /// algorithm paper in its use of separate angular and dimensionful
|
---|
54 | /// distances), but provides some extra flexibility.
|
---|
55 | ///
|
---|
56 | /// It searches for all splittings that pass a symmetry cut (zcut) and
|
---|
57 | /// then selects the one with the largest auxiliary scale choice
|
---|
58 | /// (e.g. jade distance of the splitting, kt distance of the
|
---|
59 | /// splitting, etc.)
|
---|
60 | ///
|
---|
61 | /// By default, the zcut is calculated from the fraction of the child
|
---|
62 | /// pt carried by the parent jet. If one calls set_absolute_z_cut the
|
---|
63 | /// fraction of transverse momentum will be computed wrt the original
|
---|
64 | /// jet.
|
---|
65 | ///
|
---|
66 | /// original code copyright (C) 2009 by Gavin Salam, released under
|
---|
67 | /// the GPL.
|
---|
68 | ///
|
---|
69 | /// \section desc Options
|
---|
70 | ///
|
---|
71 | /// - the distance choice: options are
|
---|
72 | /// kt2_distance : usual min(kti^2,ktj^2)DeltaR_{ij}^2
|
---|
73 | /// jade_distance : kti . ktj DeltaR_{ij}^2 (LI version of jade)
|
---|
74 | /// jade2_distance : kti . ktj DeltaR_{ij}^4 (LI version of jade * DR^2)
|
---|
75 | /// plain_distance : DeltaR_{ij}^2
|
---|
76 | /// mass_drop_distance : m_jet - max(m_parent1,m_parent2)
|
---|
77 | /// dot_product_distance: parent1.parent2
|
---|
78 | /// (kt2_distance by default)
|
---|
79 | ///
|
---|
80 | /// - the z cut (0 by default)
|
---|
81 | ///
|
---|
82 | /// - by calling set_absolute_z_cut(), one can ask that the pt
|
---|
83 | /// fraction if calculated wrt the original jet
|
---|
84 | ///
|
---|
85 | /// - by calling set_dr_min(drmin), one can ask that only the
|
---|
86 | /// recombinations where the 2 objects are (geometrically) distant
|
---|
87 | /// by at least drmin are kept in the maximisation.
|
---|
88 | ///
|
---|
89 | /// \section input Input conditions
|
---|
90 | ///
|
---|
91 | /// - the jet must have been obtained from a Cambridge/Aachen cluster
|
---|
92 | /// sequence
|
---|
93 | ///
|
---|
94 | /// \section output Output/structure
|
---|
95 | ///
|
---|
96 | /// - the element of the cluster sequence maximising the requested
|
---|
97 | /// distance (and satisfying the zcut) is returned.
|
---|
98 | ///
|
---|
99 | /// - if the original jet has no parents, it will be returned
|
---|
100 | ///
|
---|
101 | /// - the value of the "z" and distance corresponding to that history
|
---|
102 | /// element are stored and accessible through
|
---|
103 | /// result.structure_of<CASubJetTagger>().z();
|
---|
104 | /// result.structure_of<CASubJetTagger>().distance();
|
---|
105 | ///
|
---|
106 | class CASubJetTagger : public Transformer {
|
---|
107 | public:
|
---|
108 |
|
---|
109 | /// the available choices of auxiliary scale with respect to which
|
---|
110 | /// to order the splittings
|
---|
111 | enum ScaleChoice {
|
---|
112 | kt2_distance, //< usual min(kti^2,ktj^2)DeltaR_{ij}^2
|
---|
113 | jade_distance, //< kti . ktj DeltaR_{ij}^2 (LI version of jade)
|
---|
114 | jade2_distance, //< kti . ktj DeltaR_{ij}^4 (LI version of jade * DR^2)
|
---|
115 | plain_distance, //< DeltaR_{ij}^2
|
---|
116 | mass_drop_distance, //< m_jet - max(m_parent1,m_parent2)
|
---|
117 | dot_product_distance //< parent1.parent2
|
---|
118 | };
|
---|
119 |
|
---|
120 | /// just constructs
|
---|
121 | CASubJetTagger(ScaleChoice scale_choice = jade_distance,
|
---|
122 | double z_threshold = 0.1)
|
---|
123 | : _scale_choice(scale_choice), _z_threshold(z_threshold),
|
---|
124 | _dr2_min(0.0), _absolute_z_cut(false){};
|
---|
125 |
|
---|
126 | /// sets a minimum delta R below which spliting will be ignored
|
---|
127 | /// (only relevant if set prior to calling run())
|
---|
128 | void set_dr_min(double drmin) {_dr2_min = drmin*drmin;}
|
---|
129 |
|
---|
130 | /// returns a textual description of the tagger
|
---|
131 | virtual std::string description() const;
|
---|
132 |
|
---|
133 | /// If (abs_z_cut) is set to false (the default) then for a
|
---|
134 | /// splitting to be considered, each subjet must satisfy
|
---|
135 | ///
|
---|
136 | /// p_{t,sub} > z_threshold * p_{t,parent}
|
---|
137 | ///
|
---|
138 | /// whereas if it is set to true, then each subject must satisfy
|
---|
139 | ///
|
---|
140 | /// p_{t,sub} > z_threshold * p_{t,original-jet}
|
---|
141 | ///
|
---|
142 | /// where parent is the immediate parent of the splitting, and
|
---|
143 | /// original jet is the one supplied to the run() function.
|
---|
144 | ///
|
---|
145 | /// Only relevant is called prior to run().
|
---|
146 | void set_absolute_z_cut(bool abs_z_cut=true) {_absolute_z_cut = abs_z_cut;}
|
---|
147 |
|
---|
148 | /// runs the tagger on the given jet and
|
---|
149 | /// returns the tagged PseudoJet if successful, or a PseudoJet==0 otherwise
|
---|
150 | /// (standard access is through operator()).
|
---|
151 | virtual PseudoJet result(const PseudoJet & jet) const;
|
---|
152 |
|
---|
153 | /// the type of Structure returned
|
---|
154 | typedef CASubJetTaggerStructure StructureType;
|
---|
155 |
|
---|
156 | protected:
|
---|
157 | /// class that contains the result internally
|
---|
158 | class JetAux {
|
---|
159 | public:
|
---|
160 | PseudoJet jet; //< the subjet (immediate parent of splitting)
|
---|
161 | double aux_distance; //< the auxiliary distance between its two subjets
|
---|
162 | double delta_r; //< the angular distance between its two subjets
|
---|
163 | double z; //< the transverse momentum fraction
|
---|
164 | };
|
---|
165 |
|
---|
166 |
|
---|
167 | void _recurse_through_jet(const PseudoJet & current_jet,
|
---|
168 | JetAux &aux_max,
|
---|
169 | const PseudoJet & original_jet) const;
|
---|
170 |
|
---|
171 | ScaleChoice _scale_choice;
|
---|
172 | double _z_threshold;
|
---|
173 | double _dr2_min;
|
---|
174 | bool _absolute_z_cut;
|
---|
175 |
|
---|
176 | static LimitedWarning _non_ca_warnings;
|
---|
177 | };
|
---|
178 |
|
---|
179 |
|
---|
180 | //------------------------------------------------------------------------
|
---|
181 | /// @ingroup tools_taggers
|
---|
182 | /// the structure returned by a CASubJetTagger
|
---|
183 | ///
|
---|
184 | /// Since this is directly an element of the ClusterSequence, we keep
|
---|
185 | /// basically the original ClusterSequenceStructure (wrapped for
|
---|
186 | /// memory-management reasons) and add information about the pt
|
---|
187 | /// fraction and distance of the subjet structure
|
---|
188 | class CASubJetTaggerStructure : public WrappedStructure{
|
---|
189 | public:
|
---|
190 | /// default ctor
|
---|
191 | /// \param result_jet the jet for which we have to keep the
|
---|
192 | /// structure
|
---|
193 | CASubJetTaggerStructure(const PseudoJet & result_jet)
|
---|
194 | : WrappedStructure(result_jet.structure_shared_ptr()){}
|
---|
195 |
|
---|
196 | /// returns the scale choice asked for the maximisation
|
---|
197 | CASubJetTagger::ScaleChoice scale_choice() const {return _scale_choice;}
|
---|
198 |
|
---|
199 | /// returns the value of the distance measure (corresponding to
|
---|
200 | /// ScaleChoice) for this jet's splitting
|
---|
201 | double distance() const {return _distance;}
|
---|
202 |
|
---|
203 | /// returns the pt fraction contained by the softer of the two component
|
---|
204 | /// pieces of this jet (normalised relative to this jet)
|
---|
205 | double z() const {return _z;}
|
---|
206 |
|
---|
207 | /// returns the pt fraction contained by the softer of the two component
|
---|
208 | /// pieces of this jet (normalised relative to the original jet)
|
---|
209 | bool absolute_z() const {return _absolute_z;}
|
---|
210 |
|
---|
211 | // /// returns the original jet (before tagging)
|
---|
212 | // const PseudoJet & original() const {return _original_jet;}
|
---|
213 |
|
---|
214 | protected:
|
---|
215 | CASubJetTagger::ScaleChoice _scale_choice; ///< the user scale choice
|
---|
216 | double _distance; ///< the maximal distance associated with the result
|
---|
217 | bool _absolute_z; ///< whether z is computed wrt to the original jet or not
|
---|
218 | double _z; ///< the transverse momentum fraction
|
---|
219 | // PseudoJet _original_jet; ///< the original jet (before tagging)
|
---|
220 |
|
---|
221 | friend class CASubJetTagger; ///< to allow setting the internal information
|
---|
222 | };
|
---|
223 |
|
---|
224 | FASTJET_END_NAMESPACE
|
---|
225 |
|
---|
226 | #endif // __CASUBJET_HH__
|
---|