FastJet 3.4.3
Loading...
Searching...
No Matches
ClusterSequenceAreaBase.cc
1
2//FJSTARTHEADER
3// $Id$
4//
5// Copyright (c) 2005-2024, Matteo Cacciari, Gavin P. Salam and Gregory Soyez
6//
7//----------------------------------------------------------------------
8// This file is part of FastJet.
9//
10// FastJet is free software; you can redistribute it and/or modify
11// it under the terms of the GNU General Public License as published by
12// the Free Software Foundation; either version 2 of the License, or
13// (at your option) any later version.
14//
15// The algorithms that underlie FastJet have required considerable
16// development. They are described in the original FastJet paper,
17// hep-ph/0512210 and in the manual, arXiv:1111.6097. If you use
18// FastJet as part of work towards a scientific publication, please
19// quote the version you use and include a citation to the manual and
20// optionally also to hep-ph/0512210.
21//
22// FastJet is distributed in the hope that it will be useful,
23// but WITHOUT ANY WARRANTY; without even the implied warranty of
24// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
25// GNU General Public License for more details.
26//
27// You should have received a copy of the GNU General Public License
28// along with FastJet. If not, see <http://www.gnu.org/licenses/>.
29//----------------------------------------------------------------------
30//FJENDHEADER
31
32
33
34
35#include "fastjet/ClusterSequenceAreaBase.hh"
36#include <algorithm>
37
38FASTJET_BEGIN_NAMESPACE
39
40using namespace std;
41
42
43/// allow for warnings
44LimitedWarning ClusterSequenceAreaBase::_warnings;
45LimitedWarning ClusterSequenceAreaBase::_warnings_zero_area;
46LimitedWarning ClusterSequenceAreaBase::_warnings_empty_area;
47
48//----------------------------------------------------------------------
49/// return the total area, within the selector's range, that is free
50/// of jets.
51///
52/// Calculate this as (range area) - \sum_{i in range} A_i
53///
54/// for ClusterSequences with explicit ghosts, assume that there will
55/// never be any empty area, i.e. it is always filled in by pure
56/// ghosts jets. This holds for seq.rec. algorithms
57double ClusterSequenceAreaBase::empty_area(const Selector & selector) const {
58
59 if (has_explicit_ghosts()) {return 0.0;}
60 else { return empty_area_from_jets(inclusive_jets(0.0), selector);}
61
62}
63
64//----------------------------------------------------------------------
65/// return the total area, within range, that is free of jets.
66///
67/// Calculate this as (range area) - \sum_{i in range} A_i
68///
69double ClusterSequenceAreaBase::empty_area_from_jets(
70 const std::vector<PseudoJet> & all_jets,
71 const Selector & selector) const {
72 _check_selector_good_for_median(selector);
73
74 double empty = selector.area();
75 for (unsigned i = 0; i < all_jets.size(); i++) {
76 if (selector.pass(all_jets[i])) empty -= area(all_jets[i]);
77 }
78 return empty;
79}
80
81// this is deprecated but used by other deprecated methods. So we hide
82// the implementation in a protected method so that (i) it can still
83// be used internally (without generating a compile-time warning when
84// building FastJet) and the interface can be marked as deprecated.
85// This can disappear once all the public interfaces have disappeared.
86double ClusterSequenceAreaBase::median_pt_per_unit_area(const Selector & selector) const {
87 return _median_pt_per_unit_area(selector);
88}
89
90// the hidden implementation
91double ClusterSequenceAreaBase::_median_pt_per_unit_area(const Selector & selector) const {
92 return _median_pt_per_unit_something(selector,false);
93}
94
95
96
97// this is deprecated but used by other deprecated methods. So we hide
98// the implementation in a protected method so that (i) it can still
99// be used internally (without generating a compile-time warning when
100// building FastJet) and the interface can be marked as deprecated.
101// This can disappear once all the public interfaces have disappeared.
102double ClusterSequenceAreaBase::median_pt_per_unit_area_4vector(const Selector & selector) const {
103 return _median_pt_per_unit_area_4vector(selector);
104}
105
106// the deprecated interface
107double ClusterSequenceAreaBase::_median_pt_per_unit_area_4vector(const Selector & selector) const {
108 return _median_pt_per_unit_something(selector,true);
109}
110
111
112//----------------------------------------------------------------------
113// this is deprecated but used by other deprecated methods. So we hide
114// the implementation in a protected method so that (i) it can still
115// be used internally (without generating a compile-time warning when
116// building FastJet) and the interface can be marked as deprecated.
117// This can disappear once all the public interfaces have disappeared.
118double ClusterSequenceAreaBase::median_pt_per_unit_something(
119 const Selector & selector, bool use_area_4vector) const {
120 return _median_pt_per_unit_something(selector, use_area_4vector);
121}
122
123// the median of (pt/area) for jets contained within range, counting
124// the empty area as if it were made up of a collection of empty
125// jets each of area (0.55 * pi R^2).
126double ClusterSequenceAreaBase::_median_pt_per_unit_something(
127 const Selector & selector, bool use_area_4vector) const {
128 double median, sigma, mean_area;
129 _get_median_rho_and_sigma(selector, use_area_4vector, median, sigma, mean_area);
130 return median;
131}
132
133
134//----------------------------------------------------------------------
135/// fits a form pt_per_unit_area(y) = a + b*y^2 for jets in range.
136/// exclude_above allows one to exclude large values of pt/area from
137/// fit. use_area_4vector = true uses the 4vector areas.
138void ClusterSequenceAreaBase::parabolic_pt_per_unit_area(
139 double & a, double & b, const Selector & selector,
140 double exclude_above, bool use_area_4vector) const {
141 return _parabolic_pt_per_unit_area(a, b, selector, exclude_above, use_area_4vector);
142}
143
144void ClusterSequenceAreaBase::_parabolic_pt_per_unit_area(
145 double & a, double & b, const Selector & selector,
146 double exclude_above, bool use_area_4vector) const {
147 // sanity check on the selector: we require a finite area and that
148 // it applies jet by jet (see BackgroundEstimator for more advanced
149 // usage)
150 _check_selector_good_for_median(selector);
151
152 int n=0;
153 double mean_f=0, mean_x2=0, mean_x4=0, mean_fx2=0;
154
155 vector<PseudoJet> incl_jets = inclusive_jets();
156
157 for (unsigned i = 0; i < incl_jets.size(); i++) {
158 if (selector.pass(incl_jets[i])) {
159 double this_area;
160 if ( use_area_4vector ) {
161 this_area = area_4vector(incl_jets[i]).perp();
162 } else {
163 this_area = area(incl_jets[i]);
164 }
165 double f = incl_jets[i].perp()/this_area;
166 if (exclude_above <= 0.0 || f < exclude_above) {
167 double x = incl_jets[i].rap(); double x2 = x*x;
168 mean_f += f;
169 mean_x2 += x2;
170 mean_x4 += x2*x2;
171 mean_fx2 += f*x2;
172 n++;
173 }
174 }
175 }
176
177 if (n <= 1) {
178 // meaningful results require at least two jets inside the
179 // area -- mind you if there are empty jets we should be in
180 // any case doing something special...
181 a = 0.0;
182 b = 0.0;
183 } else {
184 mean_f /= n;
185 mean_x2 /= n;
186 mean_x4 /= n;
187 mean_fx2 /= n;
188
189 b = (mean_f*mean_x2 - mean_fx2)/(mean_x2*mean_x2 - mean_x4);
190 a = mean_f - b*mean_x2;
191 }
192}
193
194
195//----------------------------------------------------------------------
196void ClusterSequenceAreaBase::get_median_rho_and_sigma(
197 const Selector & selector, bool use_area_4vector,
198 double & median, double & sigma, double & mean_area) const {
199 _get_median_rho_and_sigma(selector, use_area_4vector, median, sigma, mean_area);
200}
201
202void ClusterSequenceAreaBase::_get_median_rho_and_sigma(
203 const Selector & selector, bool use_area_4vector,
204 double & median, double & sigma, double & mean_area) const {
205
206 vector<PseudoJet> incl_jets = inclusive_jets();
207 _get_median_rho_and_sigma(incl_jets, selector, use_area_4vector,
208 median, sigma, mean_area, true);
209}
210
211void ClusterSequenceAreaBase::get_median_rho_and_sigma(
212 const vector<PseudoJet> & all_jets,
213 const Selector & selector, bool use_area_4vector,
214 double & median, double & sigma, double & mean_area,
215 bool all_are_incl) const {
216 _get_median_rho_and_sigma(all_jets, selector, use_area_4vector,
217 median, sigma, mean_area, all_are_incl);
218}
219
220void ClusterSequenceAreaBase::_get_median_rho_and_sigma(
221 const vector<PseudoJet> & all_jets,
222 const Selector & selector, bool use_area_4vector,
223 double & median, double & sigma, double & mean_area,
224 bool all_are_incl) const {
225
226 _check_jet_alg_good_for_median();
227
228 // sanity check on the selector: we require a finite area and that
229 // it applies jet by jet (see BackgroundEstimator for more advanced
230 // usage)
231 _check_selector_good_for_median(selector);
232
233 vector<double> pt_over_areas;
234 double total_area = 0.0;
235 double total_njets = 0;
236
237 for (unsigned i = 0; i < all_jets.size(); i++) {
238 if (selector.pass(all_jets[i])) {
239 double this_area;
240 if (use_area_4vector) {
241 this_area = area_4vector(all_jets[i]).perp();
242 } else {
243 this_area = area(all_jets[i]);
244 }
245
246 if (this_area>0) {
247 pt_over_areas.push_back(all_jets[i].perp()/this_area);
248 } else {
249 _warnings_zero_area.warn("ClusterSequenceAreaBase::get_median_rho_and_sigma(...): discarded jet with zero area. Zero-area jets may be due to (i) too large a ghost area (ii) a jet being outside the ghost range (iii) the computation not being done using an appropriate algorithm (kt;C/A).");
250 }
251
252 total_area += this_area;
253 total_njets += 1.0;
254 }
255 }
256
257 // there is nothing inside our region, so answer will always be zero
258 if (pt_over_areas.size() == 0) {
259 median = 0.0;
260 sigma = 0.0;
261 mean_area = 0.0;
262 return;
263 }
264
265 // get median (pt/area) [this is the "old" median definition. It considers
266 // only the "real" jets in calculating the median, i.e. excluding the
267 // only-ghost ones; it will be supplemented with more info below]
268 sort(pt_over_areas.begin(), pt_over_areas.end());
269
270 // now get the median & error, accounting for empty jets
271 // define the fractions of distribution at median, median-1sigma
272 double posn[2] = {0.5, (1.0-0.6827)/2.0};
273 double res[2];
274
275 double n_empty, empty_a;
276 if (has_explicit_ghosts()) {
277 // NB: the following lines of code are potentially incorrect in cases
278 // where there are unclustered particles (empty_area would do a better job,
279 // at least for active areas). This is not an issue with kt or C/A, or other
280 // algorithms that cluster all particles (and the median estimation should in
281 // any case only be done with kt or C/A!)
282 empty_a = 0.0;
283 n_empty = 0;
284 } else if (all_are_incl) {
285 // the default case
286 empty_a = empty_area(selector);
287 n_empty = n_empty_jets(selector);
288 } else {
289 // this one is intended to be used when e.g. one runs C/A, then looks at its
290 // exclusive jets in order to get an effective smaller R value, and passes those
291 // to this routine.
292 empty_a = empty_area_from_jets(all_jets, selector);
293 mean_area = total_area / total_njets; // temporary value
294 n_empty = empty_a / mean_area;
295 }
296 //cout << "*** tot_area = " << total_area << ", empty_a = " << empty_a << endl;
297 //cout << "*** n_empty = " << n_empty << ", ntotal = " << total_njets << endl;
298 total_njets += n_empty;
299 total_area += empty_a;
300
301 // we need an int (rather than an unsigned int) with the size of the
302 // pt_over_areas array, because we'll often be doing subtraction of
303 // -1, negating it, etc. All of these operations go crazy with unsigned ints.
304 int pt_over_areas_size = pt_over_areas.size();
305 if (n_empty < -pt_over_areas_size/4.0) {
306 _warnings_empty_area.warn("ClusterSequenceAreaBase::get_median_rho_and_sigma(...): the estimated empty area is suspiciously large and negative and may lead to an over-estimation of rho. This may be due to (i) a rare statistical fluctuation or (ii) too small a range used to estimate the background properties.");
307 }
308
309 for (int i = 0; i < 2; i++) {
310 double nj_median_pos =
311 (pt_over_areas_size-1.0 + n_empty)*posn[i] - n_empty;
312 double nj_median_ratio;
313 if (nj_median_pos >= 0 && pt_over_areas_size > 1) {
314 int int_nj_median = int(nj_median_pos);
315
316 // avoid potential overflow issues
317 if (int_nj_median+1 > pt_over_areas_size-1){
318 int_nj_median = pt_over_areas_size-2;
319 nj_median_pos = pt_over_areas_size-1;
320 }
321
322 nj_median_ratio =
323 pt_over_areas[int_nj_median] * (int_nj_median+1-nj_median_pos)
324 + pt_over_areas[int_nj_median+1] * (nj_median_pos - int_nj_median);
325 } else {
326 nj_median_ratio = 0.0;
327 }
328 res[i] = nj_median_ratio;
329 }
330 median = res[0];
331 double error = res[0] - res[1];
332 mean_area = total_area / total_njets;
333 sigma = error * sqrt( max(0.0, mean_area) );
334}
335
336
337// return a vector of all subtracted jets, using area_4vector, given rho.
338// Only inclusive_jets above ptmin are subtracted and returned.
339// the ordering is the same as that of sorted_by_pt(cs.inclusive_jets()),
340// i.e. not necessarily ordered in pt once subtracted
341vector<PseudoJet> ClusterSequenceAreaBase::subtracted_jets(const double rho,
342 const double ptmin)
343 const {
344 return _subtracted_jets(rho,ptmin);
345}
346
347vector<PseudoJet> ClusterSequenceAreaBase::_subtracted_jets(const double rho,
348 const double ptmin)
349 const {
350 vector<PseudoJet> sub_jets;
351 vector<PseudoJet> jets_local = sorted_by_pt(inclusive_jets(ptmin));
352 for (unsigned i=0; i<jets_local.size(); i++) {
353 PseudoJet sub_jet = _subtracted_jet(jets_local[i],rho);
354 sub_jets.push_back(sub_jet);
355 }
356 return sub_jets;
357}
358
359// return a vector of subtracted jets, using area_4vector.
360// Only inclusive_jets above ptmin are subtracted and returned.
361// the ordering is the same as that of sorted_by_pt(cs.inclusive_jets()),
362// i.e. not necessarily ordered in pt once subtracted
363vector<PseudoJet> ClusterSequenceAreaBase::subtracted_jets(
364 const Selector & selector,
365 const double ptmin)
366 const {
367 double rho = _median_pt_per_unit_area_4vector(selector);
368 return _subtracted_jets(rho,ptmin);
369}
370
371
372/// return a subtracted jet, using area_4vector, given rho
373PseudoJet ClusterSequenceAreaBase::subtracted_jet(const PseudoJet & jet,
374 const double rho) const {
375 return _subtracted_jet(jet, rho);
376}
377
378PseudoJet ClusterSequenceAreaBase::_subtracted_jet(const PseudoJet & jet,
379 const double rho) const {
380 PseudoJet area4vect = area_4vector(jet);
381 PseudoJet sub_jet;
382 // sanity check
383 if (rho*area4vect.perp() < jet.perp() ) {
384 sub_jet = jet - rho*area4vect;
385 } else { sub_jet = PseudoJet(0.0,0.0,0.0,0.0); }
386
387 // make sure the subtracted jet has the same index (cluster, user, csw)
388 // (i.e. "looks like") the original jet
390 sub_jet.set_user_index(jet.user_index());
391 // do not use CS::_set_structure_shared_ptr here, which should
392 // only be called to maintain the tally during construction
394 return sub_jet;
395}
396
397
398/// return a subtracted jet, using area_4vector; note that this is
399/// potentially inefficient if repeatedly used for many different
400/// jets, because rho will be recalculated each time around.
401PseudoJet ClusterSequenceAreaBase::subtracted_jet(const PseudoJet & jet,
402 const Selector & selector) const {
403 return _subtracted_jet(jet, selector);
404}
405
406PseudoJet ClusterSequenceAreaBase::_subtracted_jet(const PseudoJet & jet,
407 const Selector & selector) const {
408 double rho = _median_pt_per_unit_area_4vector(selector);
409 PseudoJet sub_jet = _subtracted_jet(jet, rho);
410 return sub_jet;
411}
412
413
414/// return the subtracted pt, given rho
415double ClusterSequenceAreaBase::subtracted_pt(const PseudoJet & jet,
416 const double rho,
417 bool use_area_4vector) const {
418 return _subtracted_pt(jet, rho, use_area_4vector);
419}
420
421double ClusterSequenceAreaBase::_subtracted_pt(const PseudoJet & jet,
422 const double rho,
423 bool use_area_4vector) const {
424 if ( use_area_4vector ) {
425 PseudoJet sub_jet = _subtracted_jet(jet,rho);
426 return sub_jet.perp();
427 } else {
428 return jet.perp() - rho*area(jet);
429 }
430}
431
432
433/// return the subtracted pt; note that this is
434/// potentially inefficient if repeatedly used for many different
435/// jets, because rho will be recalculated each time around.
436double ClusterSequenceAreaBase::subtracted_pt(const PseudoJet & jet,
437 const Selector & selector,
438 bool use_area_4vector) const {
439 if ( use_area_4vector ) {
440 PseudoJet sub_jet = _subtracted_jet(jet,selector);
441 return sub_jet.perp();
442 } else {
443 double rho = _median_pt_per_unit_area(selector);
444 return _subtracted_pt(jet,rho,false);
445 }
446}
447
448// check the selector is suited for the computations i.e. applies jet
449// by jet and has a finite area
450void ClusterSequenceAreaBase::_check_selector_good_for_median(const Selector &selector) const{
451 // make sure the selector has a finite area
452 if ((! has_explicit_ghosts()) && (! selector.has_finite_area())){
453 throw Error("ClusterSequenceAreaBase: empty area can only be computed from selectors with a finite area");
454 }
455
456 // make sure the selector applies jet by jet
457 if (! selector.applies_jet_by_jet()){
458 throw Error("ClusterSequenceAreaBase: empty area can only be computed from selectors that apply jet by jet");
459 }
460}
461
462
463/// check the jet algorithm is suitable (and if not issue a warning)
464void ClusterSequenceAreaBase::_check_jet_alg_good_for_median() const {
465 if (jet_def().jet_algorithm() != kt_algorithm
466 && jet_def().jet_algorithm() != cambridge_algorithm
467 && jet_def().jet_algorithm() != cambridge_for_passive_algorithm) {
468 _warnings.warn("ClusterSequenceAreaBase: jet_def being used may not be suitable for estimating diffuse backgrounds (good options are kt, cam)");
469 }
470}
471
472
473
474FASTJET_END_NAMESPACE
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
void set_user_index(const int index)
set the user_index, intended to allow the user to add simple identifying information to a particle/je...
Definition PseudoJet.hh:392
double perp() const
returns the scalar transverse momentum
Definition PseudoJet.hh:158
void set_cluster_hist_index(const int index)
set the cluster_hist_index, intended to be used by clustering routines.
Definition PseudoJet.hh:836
int cluster_hist_index() const
return the cluster_hist_index, intended to be used by clustering routines.
Definition PseudoJet.hh:834
void set_structure_shared_ptr(const SharedPtr< PseudoJetStructureBase > &structure_in)
set the associated structure
Definition PseudoJet.cc:561
int user_index() const
return the user_index,
Definition PseudoJet.hh:389
const SharedPtr< PseudoJetStructureBase > & structure_shared_ptr() const
return a reference to the shared pointer to the PseudoJetStructureBase associated with this PseudoJet
Definition PseudoJet.cc:614
Class that encodes information about cuts and other selection criteria that can be applied to PseudoJ...
Definition Selector.hh:149
bool has_finite_area() const
returns true if it has a meaningful and finite area (i.e.
Definition Selector.hh:250
bool pass(const PseudoJet &jet) const
return true if the jet passes the selection
Definition Selector.hh:178
bool applies_jet_by_jet() const
returns true if this can be applied jet by jet
Definition Selector.hh:215
double area() const
returns the rapidity-phi area associated with the Selector (throws InvalidArea if the area does not m...
Definition Selector.cc:189
vector< PseudoJet > sorted_by_pt(const vector< PseudoJet > &jets)
return a vector of jets sorted into decreasing kt2
Definition PseudoJet.cc:871