Fork me on GitHub

source: git/external/fastjet/tools/CASubJetTagger.hh@ d7d2da3

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

move branches/ModularDelphes to trunk

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