FastJet 3.4.1
ClusterSequenceAreaBase.hh
1//FJSTARTHEADER
2// $Id$
3//
4// Copyright (c) 2005-2023, 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 __FASTJET_CLUSTERSEQUENCEAREABASE_HH__
32#define __FASTJET_CLUSTERSEQUENCEAREABASE_HH__
33
34#include "fastjet/ClusterSequence.hh"
35#include "fastjet/LimitedWarning.hh"
36#include "fastjet/Selector.hh"
37#include "fastjet/internal/deprecated.hh"
38
39FASTJET_BEGIN_NAMESPACE
40
41/// @ingroup area_classes
42/// \class ClusterSequenceAreaBase
43/// base class that sets interface for extensions of ClusterSequence
44/// that provide information about the area of each jet
45///
46/// the virtual functions here all return 0, since no area determination
47/// is implemented.
49public:
50
51 /// a constructor which just carries out the construction of the
52 /// parent class
53 template<class L> ClusterSequenceAreaBase
54 (const std::vector<L> & pseudojets,
55 const JetDefinition & jet_def_in,
56 const bool & writeout_combinations = false) :
57 ClusterSequence(pseudojets, jet_def_in, writeout_combinations) {}
58
59
60 /// default constructor
62
63
64 /// destructor
66
67
68 /// return the area associated with the given jet; this base class
69 /// returns 0.
70 virtual double area (const PseudoJet & ) const {return 0.0;}
71
72 /// return the error (uncertainty) associated with the determination
73 /// of the area of this jet; this base class returns 0.
74 virtual double area_error (const PseudoJet & ) const {return 0.0;}
75
76 /// return a PseudoJet whose 4-vector is defined by the following integral
77 ///
78 /// \int drap d\phi PseudoJet("rap,phi,pt=one") *
79 /// * Theta("rap,phi inside jet boundary")
80 ///
81 /// where PseudoJet("rap,phi,pt=one") is a 4-vector with the given
82 /// rapidity (rap), azimuth (phi) and pt=1, while Theta("rap,phi
83 /// inside jet boundary") is a function that is 1 when rap,phi
84 /// define a direction inside the jet boundary and 0 otherwise.
85 ///
86 /// This base class returns a null 4-vector.
87 virtual PseudoJet area_4vector(const PseudoJet & ) const {
88 return PseudoJet(0.0,0.0,0.0,0.0);}
89
90 /// true if a jet is made exclusively of ghosts
91 ///
92 /// NB: most area classes do not give any explicit ghost jets, but
93 /// some do, and they should replace this function with their own
94 /// version.
95 virtual bool is_pure_ghost(const PseudoJet & ) const {
96 return false;
97 }
98
99 /// returns true if ghosts are explicitly included within
100 /// jets for this ClusterSequence;
101 ///
102 /// Derived classes that do include explicit ghosts should provide
103 /// an alternative version of this routine and set it properly.
104 virtual bool has_explicit_ghosts() const {
105 return false;
106 }
107
108 /// return the total area, corresponding to the given Selector, that
109 /// is free of jets, in general based on the inclusive jets.
110 ///
111 /// The selector passed as an argument has to have a finite area and
112 /// apply jet-by-jet (see the BackgroundEstimator and Subtractor
113 /// tools for more generic usages)
114 virtual double empty_area(const Selector & selector) const;
115
116 /// return the total area, corresponding to the given Selector, that
117 /// is free of jets, based on the supplied all_jets
118 ///
119 /// The selector passed as an argument has to have a finite area and
120 /// apply jet-by-jet (see the BackgroundEstimator and Subtractor
121 /// tools for more generic usages)
122 double empty_area_from_jets(const std::vector<PseudoJet> & all_jets,
123 const Selector & selector) const;
124
125 /// return something similar to the number of pure ghost jets
126 /// in the given selector's range in an active area case.
127 /// For the local implementation we return empty_area/(0.55 pi R^2),
128 /// based on measured properties of ghost jets with kt and cam
129 /// (cf arXiv:0802.1188).
130 ///
131 /// Note that the number returned is a double.
132 ///
133 /// The selector passed as an argument has to have a finite area and
134 /// apply jet-by-jet (see the BackgroundEstimator and Subtractor
135 /// tools for more generic usages)
136 virtual double n_empty_jets(const Selector & selector) const {
137 double R = jet_def().R();
138 return empty_area(selector)/(0.55*pi*R*R);
139 }
140
141 /// the median of (pt/area) for jets contained within the selector
142 /// range, making use also of the info on n_empty_jets
143 ///
144 /// The selector passed as an argument has to have a finite area and
145 /// apply jet-by-jet (see the BackgroundEstimator and Subtractor
146 /// tools for more generic usages)
147 FASTJET_DEPRECATED_MSG("ClusterSequenceAreaBase::median_pt_per_unit_area(...) is deprecated since FastJet 3.0. Use the BackgroundEstimator series of tools instead",
148 double median_pt_per_unit_area(const Selector & selector) const);
149
150 /// the median of (pt/area_4vector) for jets contained within the
151 /// selector range, making use also of the info on n_empty_jets
152 ///
153 /// The selector passed as an argument has to have a finite area and
154 /// apply jet-by-jet
155 FASTJET_DEPRECATED_MSG("ClusterSequenceAreaBase::median_pt_per_unit_area_4vector(...) is deprecated since FastJet 3.0. Use the BackgroundEstimator series of tools instead",
156 double median_pt_per_unit_area_4vector(const Selector & selector) const);
157
158 /// the function that does the work for median_pt_per_unit_area and
159 /// median_pt_per_unit_area_4vector:
160 /// - something_is_area_4vect = false -> use plain area
161 /// - something_is_area_4vect = true -> use 4-vector area
162 FASTJET_DEPRECATED_MSG("ClusterSequenceAreaBase::median_pt_per_unit_something(...) is deprecated since FastJet 3.0. Use the BackgroundEstimator series of tools instead",
163 double median_pt_per_unit_something(const Selector & selector,
164 bool use_area_4vector) const);
165
166 /// using jets withing the selector range (and with 4-vector areas if
167 /// use_area_4vector), calculate the median pt/area, as well as an
168 /// "error" (uncertainty), which is defined as the 1-sigma
169 /// half-width of the distribution of pt/A, obtained by looking for
170 /// the point below which we have (1-0.6827)/2 of the jets
171 /// (including empty jets).
172 ///
173 /// The subtraction for a jet with uncorrected pt pt^U and area A is
174 ///
175 /// pt^S = pt^U - median*A +- sigma*sqrt(A)
176 ///
177 /// where the error is only that associated with the fluctuations
178 /// in the noise and not that associated with the noise having
179 /// caused changes in the hard-particle content of the jet.
180 ///
181 /// The selector passed as an argument has to have a finite area and
182 /// apply jet-by-jet (see the BackgroundEstimator and Subtractor
183 /// tools for more generic usages)
184 ///
185 /// NB: subtraction may also be done with 4-vector area of course,
186 /// and this is recommended for jets with larger values of R, as
187 /// long as rho has also been determined with a 4-vector area;
188 /// using a scalar area causes one to neglect terms of relative
189 /// order $R^2/8$ in the jet $p_t$.
190 //FASTJET_DEPRECATED_MSG("ClusterSequenceAreaBase::get_median_rho_and_sigma(...) is deprecated since FastJet 3.0. Use the BackgroundEstimator series of tools instead")
191 virtual void get_median_rho_and_sigma(const Selector & selector,
192 bool use_area_4vector,
193 double & median, double & sigma,
194 double & mean_area) const;
195
196 /// a more advanced version of get_median_rho_and_sigma, which allows
197 /// one to use any "view" of the event containing all jets (so that,
198 /// e.g. one might use Cam on a different resolution scale without
199 /// have to rerun the algorithm).
200 ///
201 /// By default it will assume that "all" are not inclusive jets,
202 /// so that in dealing with empty area it has to calculate
203 /// the number of empty jets based on the empty area and the
204 /// the observed <area> of jets rather than a surmised area
205 ///
206 /// Note that for small effective radii, this can cause problems
207 /// because the harder jets get an area >> <ghost-jet-area>
208 /// and so the estimate comes out all wrong. In these situations
209 /// it is highly advisable to use an area with explicit ghosts, since
210 /// then the "empty" jets are actually visible.
211 ///
212 /// The selector passed as an argument has to have a finite area and
213 /// apply jet-by-jet (see the BackgroundEstimator and Subtractor
214 /// tools for more generic usages)
215 //FASTJET_DEPRECATED_MSG("ClusterSequenceAreaBase::get_median_rho_and_sigma(...) is deprecated since FastJet 3.0. Use the BackgroundEstimator series of tools instead")
216 virtual void get_median_rho_and_sigma(const std::vector<PseudoJet> & all_jets,
217 const Selector & selector,
218 bool use_area_4vector,
219 double & median, double & sigma,
220 double & mean_area,
221 bool all_are_inclusive = false) const;
222
223 /// same as the full version of get_median_rho_and_error, but without
224 /// access to the mean_area
225 ///
226 /// The selector passed as an argument has to have a finite area and
227 /// apply jet-by-jet (see the BackgroundEstimator and Subtractor
228 /// tools for more generic usages)
229 //FASTJET_DEPRECATED_MSG("ClusterSequenceAreaBase::get_median_rho_and_sigma(...) is deprecated since FastJet 3.0. Use the BackgroundEstimator series of tools instead")
230 virtual void get_median_rho_and_sigma(const Selector & selector,
231 bool use_area_4vector,
232 double & median, double & sigma) const{
233 return _get_median_rho_and_sigma(selector, use_area_4vector, median, sigma);
234 }
235
236
237 /// fits a form pt_per_unit_area(y) = a + b*y^2 in the selector range.
238 /// exclude_above allows one to exclude large values of pt/area from fit.
239 /// (if negative, the cut is discarded)
240 /// use_area_4vector = true uses the 4vector areas.
241 ///
242 /// The selector passed as an argument has to have a finite area and
243 /// apply jet-by-jet (see the BackgroundEstimator and Subtractor
244 /// tools for more generic usages)
245 //FASTJET_DEPRECATED_MSG("ClusterSequenceAreaBase::parabolic_pt_per_unit_area(...) is deprecated since FastJet 3.0. Use the BackgroundEstimator series of tools instead")
246 virtual void parabolic_pt_per_unit_area(double & a, double & b,
247 const Selector & selector,
248 double exclude_above=-1.0,
249 bool use_area_4vector=false) const;
250
251 /// return a vector of all subtracted jets, using area_4vector, given rho.
252 /// Only inclusive_jets above ptmin are subtracted and returned.
253 /// the ordering is the same as that of sorted_by_pt(cs.inclusive_jets()),
254 /// i.e. not necessarily ordered in pt once subtracted
255 FASTJET_DEPRECATED_MSG("ClusterSequenceAreaBase::subtracted_jets(...) is deprecated since FastJet 3.0. Use the Subtractor tool (with the BackgroundEstimator series of tools) instead",
256 std::vector<PseudoJet> subtracted_jets(const double rho,
257 const double ptmin=0.0) const);
258
259 /// return a vector of subtracted jets, using area_4vector.
260 /// Only inclusive_jets above ptmin are subtracted and returned.
261 /// the ordering is the same as that of sorted_by_pt(cs.inclusive_jets()),
262 /// i.e. not necessarily ordered in pt once subtracted
263 ///
264 /// The selector passed as an argument has to have a finite area and
265 /// apply jet-by-jet (see the BackgroundEstimator and Subtractor
266 /// tools for more generic usages)
267 FASTJET_DEPRECATED_MSG("ClusterSequenceAreaBase::subtracted_jets(...) is deprecated since FastJet 3.0. Use the Subtractor tool (with the BackgroundEstimator series of tools) instead",
268 std::vector<PseudoJet> subtracted_jets(const Selector & selector,
269 const double ptmin=0.0) const);
270
271 /// return a subtracted jet, using area_4vector, given rho
272 FASTJET_DEPRECATED_MSG("ClusterSequenceAreaBase::subtracted_jet(...) is deprecated since FastJet 3.0. Use the Subtractor tool (with the BackgroundEstimator series of tools) instead",
273 PseudoJet subtracted_jet(const PseudoJet & jet,
274 const double rho) const);
275
276 /// return a subtracted jet, using area_4vector; note
277 /// that this is potentially inefficient if repeatedly used for many
278 /// different jets, because rho will be recalculated each time
279 /// around.
280 ///
281 /// The selector passed as an argument has to have a finite area and
282 /// apply jet-by-jet (see the BackgroundEstimator and Subtractor
283 /// tools for more generic usages)
284 FASTJET_DEPRECATED_MSG("ClusterSequenceAreaBase::subtracted_jet(...) is deprecated since FastJet 3.0. Use the Subtractor tool (with the BackgroundEstimator series of tools) instead",
285 PseudoJet subtracted_jet(const PseudoJet & jet,
286 const Selector & selector) const);
287
288 /// return the subtracted pt, given rho
289 FASTJET_DEPRECATED_MSG("ClusterSequenceAreaBase::subtracted_pt(...) is deprecated since FastJet 3.0. Use the Subtractor tool (with the BackgroundEstimator series of tools) instead",
290 double subtracted_pt(const PseudoJet & jet,
291 const double rho,
292 bool use_area_4vector=false) const);
293
294 /// return the subtracted pt; note that this is
295 /// potentially inefficient if repeatedly used for many different
296 /// jets, because rho will be recalculated each time around.
297 ///
298 /// The selector passed as an argument has to have a finite area and
299 /// apply jet-by-jet (see the BackgroundEstimator and Subtractor
300 /// tools for more generic usages)
301 FASTJET_DEPRECATED_MSG("ClusterSequenceAreaBase::subtracted_pt(...) is deprecated since FastJet 3.0. Use the Subtractor tool (with the BackgroundEstimator series of tools) instead",
302 double subtracted_pt(const PseudoJet & jet,
303 const Selector & selector,
304 bool use_area_4vector=false) const);
305
306protected:
307 /// check the selector is suited for the computations i.e. applies jet by jet and has a finite area
308 void _check_selector_good_for_median(const Selector &selector) const;
309
310 // the following set of protected methods are all deprecated. Their
311 // role is simply to hide the corresponding methods (without the
312 // first underscore) from the public interface so that they can be
313 // used internally until all the deprecated methods are removed.
314 // DO NOT USE ANY OF THESE METHODS: THEY ARE DEPRECATED AND WILL BE
315 // REMOVED.
316 virtual void _get_median_rho_and_sigma(const Selector & selector,
317 bool use_area_4vector,
318 double & median, double & sigma,
319 double & mean_area) const;
320 virtual void _get_median_rho_and_sigma(const std::vector<PseudoJet> & all_jets,
321 const Selector & selector,
322 bool use_area_4vector,
323 double & median, double & sigma,
324 double & mean_area,
325 bool all_are_inclusive = false) const;
326 virtual void _get_median_rho_and_sigma(const Selector & selector,
327 bool use_area_4vector,
328 double & median, double & sigma) const {
329 double mean_area;
330 _get_median_rho_and_sigma(selector, use_area_4vector,
331 median, sigma, mean_area);
332 }
333 virtual void _parabolic_pt_per_unit_area(double & a, double & b,
334 const Selector & selector,
335 double exclude_above=-1.0,
336 bool use_area_4vector=false) const;
337
338private:
339 /// handle warning messages
340 static LimitedWarning _warnings;
341 static LimitedWarning _warnings_zero_area;
342 static LimitedWarning _warnings_empty_area;
343
344 /// check the jet algorithm is suitable (and if not issue a warning)
345 void _check_jet_alg_good_for_median() const;
346
347 // the following set of private methods are all deprecated. Their
348 // role is simply to hide the corresponding methods (without the
349 // first underscore) from the public interface so that they can be
350 // used internally until all the deprecated methods are removed.
351 // DO NOT USE ANY OF THESE METHODS: THEY ARE DEPRECATED AND WILL BE
352 // REMOVED.
353 double _median_pt_per_unit_area(const Selector & selector) const;
354 double _median_pt_per_unit_area_4vector(const Selector & selector) const;
355 double _median_pt_per_unit_something(const Selector & selector,
356 bool use_area_4vector) const;
357 std::vector<PseudoJet> _subtracted_jets(const double rho,
358 const double ptmin=0.0) const;
359 PseudoJet _subtracted_jet(const PseudoJet & jet,
360 const double rho) const;
361 PseudoJet _subtracted_jet(const PseudoJet & jet,
362 const Selector & selector) const;
363 double _subtracted_pt(const PseudoJet & jet,
364 const double rho,
365 bool use_area_4vector=false) const;
366};
367
368
369
370FASTJET_END_NAMESPACE
371
372#endif // __FASTJET_CLUSTERSEQUENCEAREABASE_HH__
base class that sets interface for extensions of ClusterSequence that provide information about the a...
virtual void get_median_rho_and_sigma(const Selector &selector, bool use_area_4vector, double &median, double &sigma) const
same as the full version of get_median_rho_and_error, but without access to the mean_area
virtual bool has_explicit_ghosts() const
returns true if ghosts are explicitly included within jets for this ClusterSequence;
virtual void get_median_rho_and_sigma(const std::vector< PseudoJet > &all_jets, const Selector &selector, bool use_area_4vector, double &median, double &sigma, double &mean_area, bool all_are_inclusive=false) const
a more advanced version of get_median_rho_and_sigma, which allows one to use any "view" of the event ...
virtual double n_empty_jets(const Selector &selector) const
return something similar to the number of pure ghost jets in the given selector's range in an active ...
virtual double area_error(const PseudoJet &) const
return the error (uncertainty) associated with the determination of the area of this jet; this base c...
ClusterSequenceAreaBase(const std::vector< L > &pseudojets, const JetDefinition &jet_def_in, const bool &writeout_combinations=false)
a constructor which just carries out the construction of the parent class
virtual double area(const PseudoJet &) const
return the area associated with the given jet; this base class returns 0.
virtual PseudoJet area_4vector(const PseudoJet &) const
return a PseudoJet whose 4-vector is defined by the following integral
virtual bool is_pure_ghost(const PseudoJet &) const
true if a jet is made exclusively of ghosts
deals with clustering
class that is intended to hold a full definition of the jet clusterer
Class to contain pseudojets, including minimal information of use to jet-clustering routines.
Definition: PseudoJet.hh:68
Class that encodes information about cuts and other selection criteria that can be applied to PseudoJ...
Definition: Selector.hh:149