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