FastJet  3.4.0
GridMedianBackgroundEstimator.cc
1 //FJSTARTHEADER
2 // $Id$
3 //
4 // Copyright (c) 2005-2021, 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 
32 #include "fastjet/tools/GridMedianBackgroundEstimator.hh"
33 using namespace std;
34 
35 FASTJET_BEGIN_NAMESPACE // defined in fastjet/internal/base.hh
36 
37 
38 //----------------------------------------------------------------------
39 // setting a new event
40 //----------------------------------------------------------------------
41 // tell the background estimator that it has a new event, composed
42 // of the specified particles.
43 void GridMedianBackgroundEstimator::set_particles(const vector<PseudoJet> & particles) {
44  vector<double> scalar_pt(n_tiles(), 0.0);
45 
46  assert(all_tiles_equal_area());
47  //assert(n_good_tiles() == n_tiles()); // not needed now that we have an implementation
48 
49  _cached_estimate.reset();
50  _cached_estimate.set_has_sigma(true);
51  _cached_estimate.set_mean_area(mean_tile_area());
52 
53  // check if we need to compute only rho or both rho and rho_m
54  if (_enable_rho_m){
55  // both rho and rho_m
56  //
57  // this requires a few other variables
58  vector<double> scalar_dt(n_tiles(), 0.0);
59  double pt, dt;
60  for (unsigned i = 0; i < particles.size(); i++) {
61  int j = tile_index(particles[i]);
62  if (j >= 0){
63  pt = particles[i].pt();
64  dt = particles[i].mt() - pt;
65  if (_rescaling_class == 0){
66  scalar_pt[j] += pt;
67  scalar_dt[j] += dt;
68  } else {
69  double r = (*_rescaling_class)(particles[i]);
70  scalar_pt[j] += pt/r;
71  scalar_dt[j] += dt/r;
72  }
73  }
74  }
75  // sort things for _percentile
76  sort(scalar_dt.begin(), scalar_dt.end());
77 
78  // compute rho_m and sigma_m (see comment below for the
79  // normaliosation of sigma)
80  double p50 = _percentile(scalar_dt, 0.5);
81  _cached_estimate.set_has_rho_m(true);
82  _cached_estimate.set_rho_m(p50 / mean_tile_area());
83  _cached_estimate.set_sigma_m((p50-_percentile(scalar_dt, (1.0-0.6827)/2.0))/sqrt(mean_tile_area()));
84  } else {
85  // only rho
86  //fill(_scalar_pt.begin(), _scalar_pt.end(), 0.0);
87  for (unsigned i = 0; i < particles.size(); i++) {
88  int j = tile_index(particles[i]);
89  if (j >= 0){
90  if (_rescaling_class == 0){
91  scalar_pt[j] += particles[i].pt();
92  } else {
93  scalar_pt[j] += particles[i].pt()/(*_rescaling_class)(particles[i]);
94  }
95  }
96  }
97  }
98 
99  // if there are some "bad" tiles, then we need to exclude them from
100  // the calculation of the median. We'll do this by condensing the
101  // scalar_pt vector down to just the values for the tiles that are
102  // good.
103  //
104  // tested answers look right in "issue" 2014-08-08-testing-rect-grid
105  if (n_good_tiles() != n_tiles()) {
106  int newn = 0;
107  for (unsigned i = 0; i < scalar_pt.size(); i++) {
108  if (tile_is_good(i)) {
109  // clang gets confused with the SharedPtr swap if we don't
110  // have std:: here
111  std::swap(scalar_pt[i],scalar_pt[newn]);
112  newn++;
113  }
114  }
115  scalar_pt.resize(newn);
116  }
117 
118  // in all cases, carry on with the computation of rho
119  //
120  // first sort
121  sort(scalar_pt.begin(), scalar_pt.end());
122 
123  // then compute rho
124  //
125  // watch out: by definition, our sigma is the standard deviation of
126  // the pt density multiplied by the square root of the cell area
127  double p50 = _percentile(scalar_pt, 0.5);
128  _cached_estimate.set_rho(p50 / mean_tile_area());
129  _cached_estimate.set_sigma((p50-_percentile(scalar_pt, (1.0-0.6827)/2.0))/sqrt(mean_tile_area()));
130 
131  _cache_available = true;
132 }
133 
134 
135 //----------------------------------------------------------------------
136 // retrieving fundamental information
137 //----------------------------------------------------------------------
138 
139 // get the full set of background properties
140 BackgroundEstimate GridMedianBackgroundEstimator::estimate() const{
141  verify_particles_set();
142  return _cached_estimate;
143 }
144 
145 // get the full set of background properties for a given reference jet
146 BackgroundEstimate GridMedianBackgroundEstimator::estimate(const PseudoJet &jet) const{
147  verify_particles_set();
148  if (_rescaling_class == 0)
149  return _cached_estimate;
150 
151  BackgroundEstimate local_estimate = _cached_estimate;
152  local_estimate.apply_rescaling_factor((*_rescaling_class)(jet));
153  return local_estimate;
154 }
155 
156 
157 // get rho, the median background density per unit area
158 double GridMedianBackgroundEstimator::rho() const {
159  verify_particles_set();
160  return _cached_estimate.rho();
161 }
162 
163 
164 //----------------------------------------------------------------------
165 // get sigma, the background fluctuations per unit area; must be
166 // multipled by sqrt(area) to get fluctuations for a region of a
167 // given area.
168 double GridMedianBackgroundEstimator::sigma() const{
169  verify_particles_set();
170  return _cached_estimate.sigma();
171 }
172 
173 //----------------------------------------------------------------------
174 // get rho, the background density per unit area, locally at the
175 // position of a given jet. Note that this is not const, because a
176 // user may then wish to query other aspects of the background that
177 // could depend on the position of the jet last used for a rho(jet)
178 // determination.
179 double GridMedianBackgroundEstimator::rho(const PseudoJet & jet) {
180  //verify_particles_set();
181  double rescaling = (_rescaling_class == 0) ? 1.0 : (*_rescaling_class)(jet);
182  return rescaling*rho();
183 }
184 
185 
186 //----------------------------------------------------------------------
187 // get sigma, the background fluctuations per unit area, locally at
188 // the position of a given jet. As for rho(jet), it is non-const.
189 double GridMedianBackgroundEstimator::sigma(const PseudoJet & jet){
190  //verify_particles_set();
191  double rescaling = (_rescaling_class == 0) ? 1.0 : (*_rescaling_class)(jet);
192  return rescaling*sigma();
193 }
194 
195 //----------------------------------------------------------------------
196 // returns rho_m (particle-masses contribution to the 4-vector density)
197 double GridMedianBackgroundEstimator::rho_m() const {
198  if (! _enable_rho_m){
199  throw Error("GridMediamBackgroundEstimator: rho_m requested but rho_m calculation has been disabled.");
200  }
201  verify_particles_set();
202  return _cached_estimate.rho_m();
203 }
204 
205 
206 //----------------------------------------------------------------------
207 // returns sigma_m (particle-masses contribution to the 4-vector
208 // density); must be multipled by sqrt(area) to get fluctuations
209 // for a region of a given area.
210 double GridMedianBackgroundEstimator::sigma_m() const{
211  if (! _enable_rho_m){
212  throw Error("GridMediamBackgroundEstimator: sigma_m requested but rho_m/sigma_m calculation has been disabled.");
213  }
214  verify_particles_set();
215  return _cached_estimate.sigma_m();
216 }
217 
218 //----------------------------------------------------------------------
219 // returns rho_m locally at the position of a given jet. As for
220 // rho(jet), it is non-const.
221 double GridMedianBackgroundEstimator::rho_m(const PseudoJet & jet) {
222  //verify_particles_set();
223  double rescaling = (_rescaling_class == 0) ? 1.0 : (*_rescaling_class)(jet);
224  return rescaling*rho_m();
225 }
226 
227 
228 //----------------------------------------------------------------------
229 // returns sigma_m locally at the position of a given jet. As for
230 // rho(jet), it is non-const.
231 double GridMedianBackgroundEstimator::sigma_m(const PseudoJet & jet){
232  //verify_particles_set();
233  double rescaling = (_rescaling_class == 0) ? 1.0 : (*_rescaling_class)(jet);
234  return rescaling*sigma_m();
235 }
236 
237 //----------------------------------------------------------------------
238 // verify that particles have been set and throw an error if not
239 void GridMedianBackgroundEstimator::verify_particles_set() const {
240  if (!_cache_available) throw Error("GridMedianBackgroundEstimator::rho() or sigma() called without particles having been set");
241 }
242 
243 
244 //----------------------------------------------------------------------
245 // description
246 //----------------------------------------------------------------------
247 string GridMedianBackgroundEstimator::description() const {
248  ostringstream desc;
249  desc << "GridMedianBackgroundEstimator, with " << RectangularGrid::description();
250  return desc.str();
251 }
252 
253 
254 //----------------------------------------------------------------------
255 // configuring the behaviour
256 //----------------------------------------------------------------------
257 // Set a pointer to a class that calculates the rescaling factor as
258 // a function of the jet (position). Note that the rescaling factor
259 // is used both in the determination of the "global" rho (the pt/A
260 // of each jet is divided by this factor) and when asking for a
261 // local rho (the result is multiplied by this factor).
262 //
263 // The BackgroundRescalingYPolynomial class can be used to get a
264 // rescaling that depends just on rapidity.
265 //
266 // Note that this has to be called BEFORE any attempt to do an
267 // actual computation
268 void GridMedianBackgroundEstimator::set_rescaling_class(const FunctionOfPseudoJet<double> * rescaling_class_in) {
269  // The rescaling is taken into account when particles are set. So
270  // you need to call set_particles again if you set the rescaling
271  // class. We thus warn if there are already some available
272  // particles
273  if (_cache_available)
274  _warning_rescaling.warn("GridMedianBackgroundEstimator::set_rescaling_class(): trying to set the rescaling class when there are already particles that have been set is dangerous: the rescaling will not affect the already existing particles resulting in mis-estimation of rho. You need to call set_particles() again before proceeding with any background estimation.");
275 
276  BackgroundEstimatorBase::set_rescaling_class(rescaling_class_in);
277 }
278 
279 
280 
281 FASTJET_END_NAMESPACE // defined in fastjet/internal/base.hh
/// a class that holds the result of the calculation
void apply_rescaling_factor(double rescaling_factor)
apply a rescaling factor (to rho, rho_m, sigma, sigma_m)
base class corresponding to errors that can be thrown by FastJet
Definition: Error.hh:52
Class to contain pseudojets, including minimal information of use to jet-clustering routines.
Definition: PseudoJet.hh:68