FastJet  3.4.0
GridMedianBackgroundEstimator.hh
1 #ifndef __GRID_MEDIAN_BACKGROUND_ESTIMATOR_HH__
2 #define __GRID_MEDIAN_BACKGROUND_ESTIMATOR_HH__
3 
4 //FJSTARTHEADER
5 // $Id$
6 //
7 // Copyright (c) 2005-2021, Matteo Cacciari, Gavin P. Salam and Gregory Soyez
8 //
9 //----------------------------------------------------------------------
10 // This file is part of FastJet.
11 //
12 // FastJet is free software; you can redistribute it and/or modify
13 // it under the terms of the GNU General Public License as published by
14 // the Free Software Foundation; either version 2 of the License, or
15 // (at your option) any later version.
16 //
17 // The algorithms that underlie FastJet have required considerable
18 // development. They are described in the original FastJet paper,
19 // hep-ph/0512210 and in the manual, arXiv:1111.6097. If you use
20 // FastJet as part of work towards a scientific publication, please
21 // quote the version you use and include a citation to the manual and
22 // optionally also to hep-ph/0512210.
23 //
24 // FastJet is distributed in the hope that it will be useful,
25 // but WITHOUT ANY WARRANTY; without even the implied warranty of
26 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
27 // GNU General Public License for more details.
28 //
29 // You should have received a copy of the GNU General Public License
30 // along with FastJet. If not, see <http://www.gnu.org/licenses/>.
31 //----------------------------------------------------------------------
32 //FJENDHEADER
33 
34 
35 #include "fastjet/tools/BackgroundEstimatorBase.hh"
36 #include "fastjet/RectangularGrid.hh"
37 
38 
39 FASTJET_BEGIN_NAMESPACE // defined in fastjet/internal/base.hh
40 
41 /// @ingroup tools_background
42 /// \class GridMedianBackgroundEstimator
43 ///
44 /// Background Estimator based on the median pt/area of a set of grid
45 /// cells.
46 ///
47 /// Description of the method:
48 /// This background estimator works by projecting the event onto a
49 /// grid in rapidity and azimuth. In each grid cell, the scalar pt
50 /// sum of the particles in the cell is computed. The background
51 /// density is then estimated by the median of (scalar pt sum/cell
52 /// area) for all cells.
53 ///
54 /// Parameters:
55 /// The class takes 2 arguments: the absolute rapidity extent of the
56 /// cells and the size of the grid cells. Note that the size of the cell
57 /// will be adjusted in azimuth to satisfy the 2pi periodicity and
58 /// in rapidity to match the requested rapidity extent.
59 ///
60 /// Rescaling:
61 /// It is possible to use a rescaling profile. In this case, the
62 /// profile needs to be set before setting the particles and it will
63 /// be applied to each particle (i.e. not to each cell).
64 /// Note also that in this case one needs to call rho(jet) instead of
65 /// rho() [Without rescaling, they are identical]
66 ///
68  , public RectangularGrid
69 {
70 
71 public:
72  /// @name constructors and destructors
73  //\{
74  //----------------------------------------------------------------
75  /// \param ymax maximal absolute rapidity extent of the grid
76  /// \param requested_grid_spacing size of the grid cell. The
77  /// "real" cell size could differ due e.g. to the 2pi
78  /// periodicity in azimuthal angle (size, not area)
79  GridMedianBackgroundEstimator(double ymax, double requested_grid_spacing) :
80  RectangularGrid(ymax, requested_grid_spacing),
81  _enable_rho_m(true) {}
82 
83  //----------------------------------------------------------------
84  /// Constructor based on a user's fully specified RectangularGrid
86  RectangularGrid(grid), _enable_rho_m(true) {
87  if (!RectangularGrid::is_initialised())
88  throw Error("attempt to construct GridMedianBackgroundEstimator with uninitialised RectangularGrid");
89  }
90 
91  //----------------------------------------------------------------
92  /// Constructor with the explicit parameters for the underlying
93  /// RectangularGrid
94  ///
95  /// \param rapmin the minimum rapidity extent of the grid
96  /// \param rapmax the maximum rapidity extent of the grid
97  /// \param drap the grid spacing in rapidity
98  /// \param dphi the grid spacing in azimuth
99  /// \param tile_selector optional (geometric) selector to specify
100  /// which tiles are good; a tile is good if
101  /// a massless 4-vector at the center of the tile passes
102  /// the selection
103  GridMedianBackgroundEstimator(double rapmin_in, double rapmax_in, double drap_in, double dphi_in,
104  Selector tile_selector = Selector()) :
105  RectangularGrid(rapmin_in, rapmax_in, drap_in, dphi_in, tile_selector), _enable_rho_m(true) {}
106 
107  //\}
108 
109 
110  /// @name setting a new event
111  //\{
112  //----------------------------------------------------------------
113 
114  /// tell the background estimator that it has a new event, composed
115  /// of the specified particles.
116  void set_particles(const std::vector<PseudoJet> & particles) FASTJET_OVERRIDE;
117 
118  /// determine whether the automatic calculation of rho_m and sigma_m
119  /// is enabled (by default true)
120  void set_compute_rho_m(bool enable){ _enable_rho_m = enable; }
121 
122  //\}
123 
124  /// return a pointer to a copy of this BGE; the user is responsible
125  /// for eventually deleting the resulting object.
126  BackgroundEstimatorBase * copy() const FASTJET_OVERRIDE {
127  return new GridMedianBackgroundEstimator(*this);
128  };
129 
130 
131 
132  /// @name retrieving fundamental information
133  //\{
134  //----------------------------------------------------------------
135  /// get the full set of background properties
136  BackgroundEstimate estimate() const FASTJET_OVERRIDE;
137 
138  /// get the full set of background properties for a given reference jet
139  BackgroundEstimate estimate(const PseudoJet &jet) const FASTJET_OVERRIDE;
140 
141  /// returns rho, the median background density per unit area
142  double rho() const FASTJET_OVERRIDE;
143 
144  /// returns sigma, the background fluctuations per unit area; must be
145  /// multipled by sqrt(area) to get fluctuations for a region of a
146  /// given area.
147  double sigma() const FASTJET_OVERRIDE;
148 
149  /// returns rho, the background density per unit area, locally at the
150  /// position of a given jet. Note that this is not const, because a
151  /// user may then wish to query other aspects of the background that
152  /// could depend on the position of the jet last used for a rho(jet)
153  /// determination.
154  double rho(const PseudoJet & jet) FASTJET_OVERRIDE;
155 
156  /// returns sigma, the background fluctuations per unit area, locally at
157  /// the position of a given jet. As for rho(jet), it is non-const.
158  double sigma(const PseudoJet & jet) FASTJET_OVERRIDE;
159 
160  /// returns true if this background estimator has support for
161  /// determination of sigma
162  bool has_sigma() const FASTJET_OVERRIDE {return true;}
163 
164  //-----------------------------------------------------------------
165  /// Returns rho_m, the purely longitudinal, particle-mass-induced
166  /// component of the background density per unit area
167  double rho_m() const FASTJET_OVERRIDE;
168 
169  /// returns sigma_m, a measure of the fluctuations in the purely
170  /// longitudinal, particle-mass-induced component of the background
171  /// density per unit area; must be multipled by sqrt(area) to get
172  /// fluctuations for a region of a given area.
173  double sigma_m() const FASTJET_OVERRIDE;
174 
175  /// Returns rho_m locally at the jet position. As for rho(jet), it is non-const.
176  double rho_m(const PseudoJet & jet) FASTJET_OVERRIDE;
177 
178  /// Returns sigma_m locally at the jet position. As for rho(jet), it is non-const.
179  double sigma_m(const PseudoJet & jet) FASTJET_OVERRIDE;
180 
181  /// Returns true if this background estimator has support for
182  /// determination of rho_m.
183  ///
184  /// Note that support for sigma_m is automatic if one has sigma and
185  /// rho_m support.
186  bool has_rho_m() const FASTJET_OVERRIDE {return _enable_rho_m;}
187 
188 
189  /// returns the area of the grid cells (all identical, but
190  /// referred to as "mean" area for uniformity with JetMedianBGE).
191  double mean_area() const {return mean_tile_area();}
192  //\}
193 
194  /// @name configuring the behaviour
195  //\{
196  //----------------------------------------------------------------
197 
198  /// Set a pointer to a class that calculates the rescaling factor as
199  /// a function of the jet (position). Note that the rescaling factor
200  /// is used both in the determination of the "global" rho (the pt/A
201  /// of each jet is divided by this factor) and when asking for a
202  /// local rho (the result is multiplied by this factor).
203  ///
204  /// The BackgroundRescalingYPolynomial class can be used to get a
205  /// rescaling that depends just on rapidity.
206  ///
207  /// Note that this has to be called BEFORE any attempt to do an
208  /// actual computation
209  ///
210  /// The same profile will be used for both pt and mt (this is
211  /// probabaly a good approximation since the particle density
212  /// changes is what dominates the rapidity profile)
213  virtual void set_rescaling_class(const FunctionOfPseudoJet<double> * rescaling_class) FASTJET_OVERRIDE;
214 
215  //\}
216 
217  /// @name description
218  //\{
219  //----------------------------------------------------------------
220 
221  /// returns a textual description of the background estimator
222  std::string description() const FASTJET_OVERRIDE;
223 
224  //\}
225 
226 
227 private:
228 
229  /// verify that particles have been set and throw an error if not
230  void verify_particles_set() const;
231 
232  // information about the event
233  //std::vector<double> _scalar_pt;
234  //double _rho, _sigma, _rho_m, _sigma_m;
235  BackgroundEstimate _cached_estimate;
236  //bool _has_particles;
237  bool _enable_rho_m;
238 
239  // various warnings to inform people of potential dangers
240  LimitedWarning _warning_rho_of_jet;
241  LimitedWarning _warning_rescaling;
242 };
243 
244 FASTJET_END_NAMESPACE // defined in fastjet/internal/base.hh
245 
246 #endif // __GRID_MEDIAN_BACKGROUND_ESTIMATOR_HH__
/// a class that holds the result of the calculation
Abstract base class that provides the basic interface for classes that estimate levels of background ...
base class corresponding to errors that can be thrown by FastJet
Definition: Error.hh:52
Background Estimator based on the median pt/area of a set of grid cells.
void set_compute_rho_m(bool enable)
determine whether the automatic calculation of rho_m and sigma_m is enabled (by default true)
GridMedianBackgroundEstimator(const RectangularGrid &grid)
Constructor based on a user's fully specified RectangularGrid.
GridMedianBackgroundEstimator(double rapmin_in, double rapmax_in, double drap_in, double dphi_in, Selector tile_selector=Selector())
Constructor with the explicit parameters for the underlying RectangularGrid.
double mean_area() const
returns the area of the grid cells (all identical, but referred to as "mean" area for uniformity with...
BackgroundEstimatorBase * copy() const override
return a pointer to a copy of this BGE; the user is responsible for eventually deleting the resulting...
GridMedianBackgroundEstimator(double ymax, double requested_grid_spacing)
class to provide facilities for giving warnings up to some maximum number of times and to provide glo...
Class to contain pseudojets, including minimal information of use to jet-clustering routines.
Definition: PseudoJet.hh:68
Class that holds a generic rectangular tiling.
Class that encodes information about cuts and other selection criteria that can be applied to PseudoJ...
Definition: Selector.hh:149