Fork me on GitHub

source: git/external/fastjet/tools/GridMedianBackgroundEstimator.hh@ 51cabe8

ImprovedOutputFile Timing dual_readout llp
Last change on this file since 51cabe8 was 406b698, checked in by Pavel Demin <pavel.demin@…>, 10 years ago

replace map with vector< GridMedianBackgroundEstimator * >

  • Property mode set to 100644
File size: 10.2 KB
RevLine 
[d7d2da3]1#ifndef __GRID_MEDIAN_BACKGROUND_ESTIMATOR_HH__
2#define __GRID_MEDIAN_BACKGROUND_ESTIMATOR_HH__
3
[35cdc46]4//FJSTARTHEADER
5// $Id: GridMedianBackgroundEstimator.hh 3610 2014-08-13 09:49:28Z salam $
[d7d2da3]6//
[35cdc46]7// Copyright (c) 2005-2014, Matteo Cacciari, Gavin P. Salam and Gregory Soyez
[d7d2da3]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
[35cdc46]18// development. They are described in the original FastJet paper,
19// hep-ph/0512210 and in the manual, arXiv:1111.6097. If you use
[d7d2da3]20// FastJet as part of work towards a scientific publication, please
[35cdc46]21// quote the version you use and include a citation to the manual and
22// optionally also to hep-ph/0512210.
[d7d2da3]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//----------------------------------------------------------------------
[35cdc46]32//FJENDHEADER
[d7d2da3]33
34
35#include "fastjet/tools/BackgroundEstimatorBase.hh"
36
[35cdc46]37// if defined then we'll use the RectangularGrid class
38//
39// (For FastJet 3.2, maybe remove the symbol and simply clean up the
40// code below to use exclusively the RectangularGrid)
41#define FASTJET_GMBGE_USEFJGRID
42
43#ifdef FASTJET_GMBGE_USEFJGRID
44#include "fastjet/RectangularGrid.hh"
45#endif
46
47
48
[d7d2da3]49FASTJET_BEGIN_NAMESPACE // defined in fastjet/internal/base.hh
50
51/// @ingroup tools_background
52/// \class GridMedianBackgroundEstimator
53///
54/// Background Estimator based on the median pt/area of a set of grid
55/// cells.
56///
57/// Description of the method:
58/// This background estimator works by projecting the event onto a
59/// grid in rapidity and azimuth. In each grid cell, the scalar pt
60/// sum of the particles in the cell is computed. The background
61/// density is then estimated by the median of (scalar pt sum/cell
62/// area) for all cells.
63///
64/// Parameters:
65/// The class takes 2 arguments: the absolute rapidity extent of the
66/// cells and the size of the grid cells. Note that the size of the cell
67/// will be adjusted in azimuth to satisfy the 2pi periodicity and
68/// in rapidity to match the requested rapidity extent.
69///
70/// Rescaling:
71/// It is possible to use a rescaling profile. In this case, the
72/// profile needs to be set before setting the particles and it will
73/// be applied to each particle (i.e. not to each cell).
74/// Note also that in this case one needs to call rho(jet) instead of
75/// rho() [Without rescaling, they are identical]
76///
[35cdc46]77class GridMedianBackgroundEstimator : public BackgroundEstimatorBase
78#ifdef FASTJET_GMBGE_USEFJGRID
[406b698]79 , public RectangularGrid
[35cdc46]80#endif
81{
82
[d7d2da3]83public:
84 /// @name constructors and destructors
85 //\{
[35cdc46]86#ifdef FASTJET_GMBGE_USEFJGRID
[d7d2da3]87 //----------------------------------------------------------------
88 /// \param ymax maximal absolute rapidity extent of the grid
89 /// \param requested_grid_spacing size of the grid cell. The
90 /// "real" cell size could differ due e.g. to the 2pi
91 /// periodicity in azimuthal angle (size, not area)
[35cdc46]92 GridMedianBackgroundEstimator(double ymax, double requested_grid_spacing) :
93 RectangularGrid(ymax, requested_grid_spacing),
94 _has_particles(false), _enable_rho_m(true) {}
95
[406b698]96 /// ctor with more control over initialisation
97 /// \param rapmin the minimum rapidity extent of the grid
98 /// \param rapmax the maximum rapidity extent of the grid
99 /// \param drap the grid spacing in rapidity
100 /// \param dphi the grid spacing in azimuth
101 /// \param tile_selector optional (geometric) selector to specify
102 /// which tiles are good; a tile is good if
103 /// a massless 4-vector at the center of the tile passes
104 /// the selection
105 GridMedianBackgroundEstimator(double rapmin_in, double rapmax_in, double drap_in, double dphi_in,
106 Selector tile_selector = Selector()) :
107 RectangularGrid(rapmin_in, rapmax_in, drap_in, dphi_in, tile_selector),
108 _has_particles(false), _enable_rho_m(true) {}
109
[35cdc46]110 //----------------------------------------------------------------
111 /// Constructor based on a user's fully specified RectangularGrid
112 GridMedianBackgroundEstimator(const RectangularGrid & grid) :
113 RectangularGrid(grid),
114 _has_particles(false), _enable_rho_m(true) {
115 if (!RectangularGrid::is_initialised())
116 throw Error("attempt to construct GridMedianBackgroundEstimator with uninitialised RectangularGrid");
117 }
118
119#else // alternative in old framework where we didn't have the rectangular grid
[d7d2da3]120 GridMedianBackgroundEstimator(double ymax, double requested_grid_spacing) :
121 _ymin(-ymax), _ymax(ymax),
122 _requested_grid_spacing(requested_grid_spacing),
[35cdc46]123 _has_particles(false), _enable_rho_m(true)
124 {
125 setup_grid();
126 }
127#endif // FASTJET_GMBGE_USEFJGRID
128
[d7d2da3]129 //\}
130
131
132 /// @name setting a new event
133 //\{
134 //----------------------------------------------------------------
135
136 /// tell the background estimator that it has a new event, composed
137 /// of the specified particles.
138 void set_particles(const std::vector<PseudoJet> & particles);
139
[35cdc46]140 /// determine whether the automatic calculation of rho_m and sigma_m
141 /// is enabled (by default true)
142 void set_compute_rho_m(bool enable){ _enable_rho_m = enable;}
143
[d7d2da3]144 //\}
145
146 /// @name retrieving fundamental information
147 //\{
148 //----------------------------------------------------------------
149
150 /// returns rho, the median background density per unit area
151 double rho() const;
152
153 /// returns sigma, the background fluctuations per unit area; must be
154 /// multipled by sqrt(area) to get fluctuations for a region of a
155 /// given area.
156 double sigma() const;
157
158 /// returns rho, the background density per unit area, locally at the
159 /// position of a given jet. Note that this is not const, because a
160 /// user may then wish to query other aspects of the background that
161 /// could depend on the position of the jet last used for a rho(jet)
162 /// determination.
163 double rho(const PseudoJet & jet);
164
165 /// returns sigma, the background fluctuations per unit area, locally at
166 /// the position of a given jet. As for rho(jet), it is non-const.
167 double sigma(const PseudoJet & jet);
168
169 /// returns true if this background estimator has support for
170 /// determination of sigma
171 bool has_sigma() {return true;}
172
[35cdc46]173 //-----------------------------------------------------------------
174 /// Returns rho_m, the purely longitudinal, particle-mass-induced
175 /// component of the background density per unit area
176 double rho_m() const;
177
178 /// returns sigma_m, a measure of the fluctuations in the purely
179 /// longitudinal, particle-mass-induced component of the background
180 /// density per unit area; must be multipled by sqrt(area) to get
181 /// fluctuations for a region of a given area.
182 double sigma_m() const;
183
184 /// Returns rho_m locally at the jet position. As for rho(jet), it is non-const.
185 double rho_m(const PseudoJet & jet);
186
187 /// Returns sigma_m locally at the jet position. As for rho(jet), it is non-const.
188 double sigma_m(const PseudoJet & jet);
189
190 /// Returns true if this background estimator has support for
191 /// determination of rho_m.
192 ///
193 /// Note that support for sigma_m is automatic is one has sigma and
194 /// rho_m support.
195 bool has_rho_m() const {return _enable_rho_m;}
196
197
[d7d2da3]198 /// returns the area of the grid cells (all identical, but
199 /// referred to as "mean" area for uniformity with JetMedianBGE).
[35cdc46]200 double mean_area() const {return mean_tile_area();}
[d7d2da3]201 //\}
202
203 /// @name configuring the behaviour
204 //\{
205 //----------------------------------------------------------------
206
207 /// Set a pointer to a class that calculates the rescaling factor as
208 /// a function of the jet (position). Note that the rescaling factor
209 /// is used both in the determination of the "global" rho (the pt/A
210 /// of each jet is divided by this factor) and when asking for a
211 /// local rho (the result is multiplied by this factor).
212 ///
213 /// The BackgroundRescalingYPolynomial class can be used to get a
214 /// rescaling that depends just on rapidity.
215 ///
216 /// Note that this has to be called BEFORE any attempt to do an
217 /// actual computation
[35cdc46]218 ///
219 /// The same profile will be used for both pt and mt (this is
220 /// probabaly a good approximation since the particle density
221 /// changes is what dominates the rapidity profile)
[d7d2da3]222 virtual void set_rescaling_class(const FunctionOfPseudoJet<double> * rescaling_class);
223
224 //\}
225
226 /// @name description
227 //\{
228 //----------------------------------------------------------------
229
230 /// returns a textual description of the background estimator
231 std::string description() const;
232
233 //\}
234
235
236private:
[35cdc46]237
238#ifndef FASTJET_GMBGE_USEFJGRID
239
[d7d2da3]240 /// configure the grid
241 void setup_grid();
242
243 /// retrieve the grid cell index for a given PseudoJet
[35cdc46]244 int tile_index(const PseudoJet & p) const;
[d7d2da3]245
246 // information about the grid
[35cdc46]247 double _ymin, _ymax, _dy, _dphi, _requested_grid_spacing, _tile_area;
[d7d2da3]248 int _ny, _nphi, _ntotal;
249
[35cdc46]250 int n_tiles() const {return _ntotal;}
251 int n_good_tiles() const {return n_tiles();}
252 int tile_is_good(int /* itile */) const {return true;}
253
254 double mean_tile_area() const {return _tile_area;}
255#endif // FASTJET_GMBGE_USEFJGRID
256
257
258 /// verify that particles have been set and throw an error if not
259 void verify_particles_set() const;
260
[d7d2da3]261 // information abotu the event
[35cdc46]262 //std::vector<double> _scalar_pt;
263 double _rho, _sigma, _rho_m, _sigma_m;
[d7d2da3]264 bool _has_particles;
[35cdc46]265 bool _enable_rho_m;
[d7d2da3]266
[35cdc46]267 // various warnings to inform people of potential dangers
[d7d2da3]268 LimitedWarning _warning_rho_of_jet;
269 LimitedWarning _warning_rescaling;
270};
271
272FASTJET_END_NAMESPACE // defined in fastjet/internal/base.hh
273
274#endif // __GRID_MEDIAN_BACKGROUND_ESTIMATOR_HH__
Note: See TracBrowser for help on using the repository browser.