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