FastJet  3.4.0-beta.1
JetMedianBackgroundEstimator.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 #include "fastjet/tools/JetMedianBackgroundEstimator.hh"
32 #include "fastjet/ClusterSequenceArea.hh"
33 #include "fastjet/ClusterSequenceStructure.hh"
34 #include <iostream>
35 #include <sstream>
36 
37 using namespace std;
38 
39 FASTJET_BEGIN_NAMESPACE // defined in fastjet/internal/base.hh
40 
41 double BackgroundJetScalarPtDensity::result(const PseudoJet & jet) const {
42  // do not include the ghosts in the list of constituents to have a
43  // correct behaviour when _pt_power is <= 0
44  std::vector<PseudoJet> constituents = (!SelectorIsPureGhost())(jet.constituents());
45  double scalar_pt = 0;
46  for (unsigned i = 0; i < constituents.size(); i++) {
47  scalar_pt += pow(constituents[i].perp(), _pt_power);
48  }
49  return scalar_pt / jet.area();
50 }
51 
52 
53 std::string BackgroundJetScalarPtDensity::description() const {
54  ostringstream oss;
55  oss << "BackgroundScalarJetPtDensity";
56  if (_pt_power != 1.0) oss << " with pt_power = " << _pt_power;
57  return oss.str();
58 }
59 
60 
61 //----------------------------------------------------------------------
62 double BackgroundRescalingYPolynomial::result(const PseudoJet & jet) const {
63  double y = jet.rap();
64  double y2 = y*y;
65  double rescaling = _a0 + _a1*y + _a2*y2 + _a3*y2*y + _a4*y2*y2;
66  return rescaling;
67 }
68 
69 /// allow for warnings
70 LimitedWarning JetMedianBackgroundEstimator::_warnings;
71 LimitedWarning JetMedianBackgroundEstimator::_warnings_zero_area;
72 LimitedWarning JetMedianBackgroundEstimator::_warnings_preliminary;
73 
74 
75 //---------------------------------------------------------------------
76 // class JetMedianBackgroundEstimator
77 // Class to estimate the density of the background per unit area
78 //---------------------------------------------------------------------
79 
80 //----------------------------------------------------------------------
81 // ctors and dtors
82 //----------------------------------------------------------------------
83 // ctor that allows to set only the particles later on
84 JetMedianBackgroundEstimator::JetMedianBackgroundEstimator(const Selector &rho_range,
85  const JetDefinition &jet_def,
86  const AreaDefinition &area_def)
87  : _rho_range(rho_range), _jet_def(jet_def), _area_def(area_def){
88 
89  // initialise things decently
90  reset();
91 
92  // make a few checks
93  _check_jet_alg_good_for_median();
94 }
95 
96 
97 //----------------------------------------------------------------------
98 // ctor from a cluster sequence
99 // - csa the ClusterSequenceArea to use
100 // - rho_range the Selector specifying which jets will be considered
102  : _rho_range(rho_range), _jet_def(JetDefinition()){
103 
104  // initialise things properly
105  reset();
106 
107  // tell the BGE about the cluster sequence
109 }
110 
111 
112 //----------------------------------------------------------------------
113 // setting a new event
114 //----------------------------------------------------------------------
115 // tell the background estimator that it has a new event, composed
116 // of the specified particles.
117 void JetMedianBackgroundEstimator::set_particles(const vector<PseudoJet> & particles) {
118  // pass an empty seed vector to the full set_particles method to tell it to use
119  // default seeds rather than fixed seeds
120  vector<int> seed;
121  set_particles_with_seed(particles, seed);
122 }
123 
124 
125 // tell the background estimator that it has a new event, composed
126 // of the specified particles and use the supplied seed for the
127 // generation of ghosts. If the seed is empty, it is ignored.
128 void JetMedianBackgroundEstimator::set_particles_with_seed(const vector<PseudoJet> & particles, const vector<int> & seed) {
129  // make sure that we have been provided a genuine jet definition
130  if (_jet_def.jet_algorithm() == undefined_jet_algorithm)
131  throw Error("JetMedianBackgroundEstimator::set_particles can only be called if you set the jet (and area) definition explicitly through the class constructor");
132 
133  // initialise things decently (including setting uptodate to false!)
134  //reset();
135 
136  // cluster the particles
137  //
138  // One may argue that it is better to cache the particles and only
139  // do the clustering later but clustering the particles now has 2
140  // practical advantages:
141  // - it allows us to use only '_included_jets' in all that follows
142  // - it avoids adding another flag to ensure particles are
143  // clustered only once
144  ClusterSequenceArea *csa;
145  if (seed.size() == 0) {
146  csa = new ClusterSequenceArea(particles, _jet_def, _area_def);
147  } else {
148  csa = new ClusterSequenceArea(particles, _jet_def, _area_def.with_fixed_seed(seed));
149  }
150 //THREAD-SAFETY-QUESTION: #ifdef FASTJET_HAVE_LIMITED_THREAD_SAFETY
151 //THREAD-SAFETY-QUESTION: // before caching thing, lock things down to avoid concurrency issues
152 //THREAD-SAFETY-QUESTION: std::lock_guard<std::mutex> guard(_jets_caching_mutex);
153 //THREAD-SAFETY-QUESTION: #endif
154  _included_jets = csa->inclusive_jets();
155 
156  // store the CS for later on
157  _csi = csa->structure_shared_ptr();
159 
160  _set_cache_unavailable();
161 
162 //THREAD-SAFETY-QUESTION: // in thread-safe mode, the lock will automatically be released here
163 }
164 
165 //----------------------------------------------------------------------
166 // (re)set the cluster sequence (with area support) to be used by
167 // future calls to rho() etc.
168 //
169 // \param csa the cluster sequence area
170 //
171 // Pre-conditions:
172 // - one should be able to estimate the "empty area" (i.e. the area
173 // not occupied by jets). This is feasible if at least one of the following
174 // conditions is satisfied:
175 // ( i) the ClusterSequence has explicit ghosts
176 // (ii) the range selected has a computable area.
177 // - the jet algorithm must be suited for median computation
178 // (otherwise a warning will be issues)
179 //
180 // Note that selectors with e.g. hardest-jets exclusion do not have
181 // a well-defined area. For this reasons, it is STRONGLY advised to
182 // use an area with explicit ghosts.
184  // sanity checks
185  //---------------
186  // (i) check that, if there are no explicit ghosts, the selector has a finite area
187  if ((!csa.has_explicit_ghosts()) && (!_rho_range.has_finite_area())){
188  throw Error("JetMedianBackgroundEstimator: either an area with explicit ghosts (recommended) or a Selector with finite area is needed (to allow for the computation of the empty area)");
189  }
190 
191 //THREAD-SAFETY-QUESTION: #ifdef FASTJET_HAVE_LIMITED_THREAD_SAFETY
192 //THREAD-SAFETY-QUESTION: // before caching thing, lock things down to avoid concurrency issues
193 //THREAD-SAFETY-QUESTION: std::lock_guard<std::mutex> guard(_jets_caching_mutex);
194 //THREAD-SAFETY-QUESTION: #endif
195  _csi = csa.structure_shared_ptr();
196 
197  // (ii) check the alg is appropriate
198  _check_jet_alg_good_for_median();
199 
200 
201  // get the initial list of jets
202  _included_jets = csa.inclusive_jets();
203 
204  _set_cache_unavailable();
205 
206 //THREAD-SAFETY-QUESTION: // in thread-safe mode, the lock will automatically be released here
207 }
208 
209 
210 //----------------------------------------------------------------------
211 // (re)set the jets (which must have area support) to be used by future
212 // calls to rho() etc.; for the conditions that must be satisfied
213 // by the jets, see the Constructor that takes jets.
214 void JetMedianBackgroundEstimator::set_jets(const vector<PseudoJet> &jets) {
215 
216  if (! jets.size())
217  throw Error("JetMedianBackgroundEstimator::JetMedianBackgroundEstimator: At least one jet is needed to compute the background properties");
218 
219  // sanity checks
220  //---------------
221  // (o) check that there is an underlying CS shared by all the jets
222  if (! (jets[0].has_associated_cluster_sequence()) && (jets[0].has_area()))
223  throw Error("JetMedianBackgroundEstimator::JetMedianBackgroundEstimator: the jets used to estimate the background properties must be associated with a valid ClusterSequenceAreaBase");
224 
225  SharedPtr<PseudoJetStructureBase> csi_shared = jets[0].structure_shared_ptr();
226  ClusterSequenceStructure * csi = dynamic_cast<ClusterSequenceStructure*>(csi_shared.get());
227  const ClusterSequenceAreaBase * csab = csi->validated_csab();
228 
229  for (unsigned int i=1;i<jets.size(); i++){
230  if (! jets[i].has_associated_cluster_sequence()) // area automatic if the next test succeeds
231  throw Error("JetMedianBackgroundEstimator::set_jets(...): the jets used to estimate the background properties must be associated with a valid ClusterSequenceAreaBase");
232 
233  if (jets[i].structure_shared_ptr().get() != csi_shared.get())
234  throw Error("JetMedianBackgroundEstimator::set_jets(...): all the jets used to estimate the background properties must share the same ClusterSequence");
235  }
236 
237  // (i) check that, if there are no explicit ghosts, the selector has a finite area
238  if ((!csab->has_explicit_ghosts()) && (!_rho_range.has_finite_area())){
239  throw Error("JetMedianBackgroundEstimator: either an area with explicit ghosts (recommended) or a Selector with finite area is needed (to allow for the computation of the empty area)");
240  }
241 
242 //THREAD-SAFETY-QUESTION: #ifdef FASTJET_HAVE_LIMITED_THREAD_SAFETY
243 //THREAD-SAFETY-QUESTION: // before caching thing, lock things down to avoid concurrency issues
244 //THREAD-SAFETY-QUESTION: std::lock_guard<std::mutex> guard(_jets_caching_mutex);
245 //THREAD-SAFETY-QUESTION: #endif
246  _csi = csi_shared;
247 
248  // (ii) check the alg is appropriate
249  _check_jet_alg_good_for_median();
250 
251  // get the initial list of jets
252  _included_jets = jets;
253 
254  // ensure recalculation of quantities that need it
255  _set_cache_unavailable();
256 
257 //THREAD-SAFETY-QUESTION: // in thread-safe mode, the lock will automatically be released here
258 }
259 
260 //----------------------------------------------------------------------
261 // retrieving fundamental information
262 //----------------------------------------------------------------------
263 
264 // get the full set of background properties
265 //
266 // For background estimators using a local ranges, this throws an
267 // error (use estimate(jet) instead)
268 // In the presence of a rescaling, the rescaling factor is not taken
269 // into account
271  if (_rho_range.takes_reference())
272  throw Error("The background estimation is obtained from a selector that takes a reference jet. estimate(PseudoJet) should be used in that case");
273 
274  if (!_cache_available) _compute_and_cache_no_overwrite();
275  return _cached_estimate;
276 }
277 
278 // get the full set of background properties for a given reference jet
279 // This does not affect the cache
281  // first compute an optional rescaling factor
282  double rescaling_factor = (_rescaling_class != 0)
283  ? (*_rescaling_class)(jet) : 1.0;
284  BackgroundEstimate local_estimate;
285 
286  // adopt a different strategy for ranges taking a reference and others
287  if (_rho_range.takes_reference()){
288  // we compute the background and rescale it (no caching)
289  local_estimate = _compute(jet);
290  } else {
291  // otherwise, we're in a situation where things can be cached once
292  // and for all and then the cache can be used frely
293  if (!_cache_available) _compute_and_cache_no_overwrite();
294  local_estimate = _cached_estimate;
295  }
296  local_estimate.apply_rescaling_factor(rescaling_factor);
297  return local_estimate;
298 }
299 
300 
301 //------
302 // get rho, the median background density per unit area
304  if (_rho_range.takes_reference())
305  throw Error("The background estimation is obtained from a selector that takes a reference jet. rho(PseudoJet) should be used in that case");
306 
307  // we are in a situation where the cache only needs to be computed
308  // once, but once it has been computed, we can use it freely.
309  if (!_cache_available) _compute_and_cache_no_overwrite();
310  return _cached_estimate.rho();
311 }
312 
313 // get sigma, the background fluctuations per unit area
315  if (_rho_range.takes_reference())
316  throw Error("The background estimation is obtained from a selector that takes a reference jet. sigma(PseudoJet) should be used in that case");
317  if (!_cache_available) _compute_and_cache_no_overwrite();
318  return _cached_estimate.sigma();
319 }
320 
321 // get rho, the median background density per unit area, locally at
322 // the position of a given jet.
323 //
324 // If the Selector associated with the range takes a reference jet
325 // (i.e. is relocatable), then for subsequent operations the
326 // Selector has that jet set as its reference.
328  // first compute an optional rescaling factor
329  double rescaling_factor = (_rescaling_class != 0)
330  ? (*_rescaling_class)(jet) : 1.0;
331 
332  // adopt a different strategy for ranges taking a reference and others
333  if (_rho_range.takes_reference()){
334  // we compute the background and use it
335  BackgroundEstimate estimate = _compute_and_cache_if_needed(jet);
336  return rescaling_factor * estimate.rho();
337  }
338 
339  // otherwise, we're in a situation where things can be cached once
340  // and for all and then the cache can be used frely
341  if (!_cache_available) _compute_and_cache_no_overwrite();
342  return rescaling_factor * _cached_estimate.rho();
343 }
344 
345 // get sigma, the background fluctuations per unit area,
346 // locally at the position of a given jet.
347 //
348 // If the Selector associated with the range takes a reference jet
349 // (i.e. is relocatable), then for subsequent operations the
350 // Selector has that jet set as its reference.
352  // first compute an optional rescaling factor
353  double rescaling_factor = (_rescaling_class != 0)
354  ? (*_rescaling_class)(jet) : 1.0;
355 
356  // see "rho(jet)" for a descrtiption of the strategy
357  if (_rho_range.takes_reference()){
358  BackgroundEstimate estimate = _compute_and_cache_if_needed(jet);
359  return rescaling_factor * estimate.sigma();
360  }
361  // otherwise, cache things once and for all
362  if (!_cache_available) _compute_and_cache_no_overwrite();
363  return rescaling_factor * _cached_estimate.sigma();
364 }
365 
366 
367 //----------------------------------------------------------------------
368 // returns rho_m (particle-masses contribution to the 4-vector density)
370  if (! has_rho_m()){
371  throw Error("JetMediamBackgroundEstimator: rho_m requested but rho_m calculation is disabled (either eplicitly or due to the presence of a jet density class).");
372  }
373  if (_rho_range.takes_reference()){
374  throw Error("The background estimation is obtained from a selector that takes a reference jet. rho_m(PseudoJet) should be used in that case");
375  }
376  if (!_cache_available) _compute_and_cache_no_overwrite();
377  return _cached_estimate.rho_m();
378 }
379 
380 
381 //----------------------------------------------------------------------
382 // returns sigma_m (particle-masses contribution to the 4-vector
383 // density); must be multipled by sqrt(area) to get fluctuations
384 // for a region of a given area.
386  if (! has_rho_m()){
387  throw Error("JetMediamBackgroundEstimator: sigma_m requested but rho_m/sigma_m calculation is disabled (either explicitly or due to the presence of a jet density class).");
388  }
389  if (_rho_range.takes_reference())
390  throw Error("The background estimation is obtained from a selector that takes a reference jet. sigma_m(PseudoJet) should be used in that case");
391  if (!_cache_available) _compute_and_cache_no_overwrite();
392  return _cached_estimate.sigma_m();
393 }
394 
395 //----------------------------------------------------------------------
396 // returns rho_m locally at the position of a given jet. As for
397 // rho(jet), it is non-const.
399  // first compute an optional rescaling factor
400  double rescaling_factor = (_rescaling_class != 0)
401  ? (*_rescaling_class)(jet) : 1.0;
402 
403  // see "rho(jet)" for a descrtiption of the strategy
404  if (_rho_range.takes_reference()){
405  BackgroundEstimate estimate = _compute_and_cache_if_needed(jet);
406  return rescaling_factor * estimate.rho_m();
407  }
408  // otherwise, cache things once and for all
409  if (!_cache_available) _compute_and_cache_no_overwrite();
410  return rescaling_factor * _cached_estimate.rho_m();
411 }
412 
413 
414 //----------------------------------------------------------------------
415 // returns sigma_m locally at the position of a given jet. As for
416 // rho(jet), it is non-const.
418  // first compute an optional rescaling factor
419  double rescaling_factor = (_rescaling_class != 0)
420  ? (*_rescaling_class)(jet) : 1.0;
421 
422  // see "rho(jet)" for a descrtiption of the strategy
423  if (_rho_range.takes_reference()){
424  BackgroundEstimate estimate = _compute_and_cache_if_needed(jet);
425  return rescaling_factor * estimate.sigma_m();
426  }
427  // otherwise, cache things once and for all
428  if (!_cache_available) _compute_and_cache_no_overwrite();
429  return rescaling_factor * _cached_estimate.sigma_m();
430 }
431 
432 
433 //----------------------------------------------------------------
434 /// Returns the mean area of the jets used to actually compute the
435 /// background properties in the last call of rho() or sigma()
436 /// If the configuration has changed in the meantime, throw an error.
438  // if the selector takes a reference, we need to use the cache
439  if (_rho_range.takes_reference()){
440  // lock to make sure no other thread interferes w the caching
441  _lock_if_needed();
442  if (!_cache_available){
443  _unlock_if_needed();
444  throw Error("Calls to JetMedianBackgroundEstimator::mean_area() in cases where the background estimation uses a selector that takes a reference jet need to call a method that fills the cached estimate (rho(jet), sigma(jet), ...).");
445  }
446  double return_value = _cached_estimate.mean_area();
447  _unlock_if_needed();
448  return return_value;
449  }
450  if (!_cache_available) _compute_and_cache_no_overwrite();
451  return _cached_estimate.mean_area();
452 }
453 
454 /// returns the number of jets used to actually compute the
455 /// background properties in the last call of rho() or sigma()
456 /// If the configuration has changed in the meantime, throw an error.
458  // if the selector takes a reference, we need to use the cache
459  if (_rho_range.takes_reference()){
460  // lock to make sure no other thread interferes w the caching
461  _lock_if_needed();
462  if (!_cache_available){
463  _unlock_if_needed();
464  throw Error("Calls to JetMedianBackgroundEstimator::n_jets_used() in cases where the background estimation uses a selector that takes a reference jet need to call a method that fills the cached estimate (rho(jet), sigma(jet), ...).");
465  }
466  unsigned int return_value = _cached_estimate.extras<JetMedianBackgroundEstimator>().n_jets_used();
467  _unlock_if_needed();
468  return return_value;
469  }
470  if (!_cache_available) _compute_and_cache_no_overwrite();
472 }
473 
474 /// returns the jets used to actually compute the background
475 /// properties
476 std::vector<PseudoJet> JetMedianBackgroundEstimator::jets_used() const{
477  vector<PseudoJet> tmp_jets;
478 
479  // if the selector takes a reference, we need to use the cache
480  if (_rho_range.takes_reference()){
481  // lock to make sure no other thread interferes w the caching
482  _lock_if_needed();
483  if (!_cache_available){
484  _unlock_if_needed();
485  throw Error("Calls to JetMedianBackgroundEstimator::jets_used() in cases where the background estimation uses a selector that takes a reference jet need to call a method that fills the cached estimate (rho(jet), sigma(jet), ...).");
486  }
487  PseudoJet reference_jet = _cached_estimate.extras<JetMedianBackgroundEstimator>().reference_jet();
488  _unlock_if_needed();
489  Selector local_rho_range = _rho_range;
490  local_rho_range.set_reference(reference_jet);
491  tmp_jets = local_rho_range(_included_jets);
492  } else {
493  if (!_cache_available) _compute_and_cache_no_overwrite();
494  tmp_jets = _rho_range(_included_jets);
495  }
496 
497  std::vector<PseudoJet> used_jets;
498  for (unsigned int i=0; i<tmp_jets.size(); i++){
499  if (tmp_jets[i].area()>0) used_jets.push_back(tmp_jets[i]);
500  }
501  return used_jets;
502 }
503 
504 /// Returns the estimate of the area (within the range defined by
505 /// the selector) that is not occupied by jets. The value is that
506 /// for the last call of rho() or sigma()
507 /// If the configuration has changed in the meantime, throw an error.
508 ///
509 /// The answer is defined to be zero if the area calculation
510 /// involved explicit ghosts; if the area calculation was an active
511 /// area, then use is made of the active area's internal list of
512 /// pure ghost jets (taking those that pass the selector); otherwise
513 /// it is based on the difference between the selector's total area
514 /// and the area of the jets that pass the selector.
515 ///
516 /// The result here is just the cached result of the corresponding
517 /// call to the ClusterSequenceAreaBase function.
519  // if the selector takes a reference, we need to use the cache
520  if (_rho_range.takes_reference()){
521  // lock to make sure no other thread interferes w the caching
522  _lock_if_needed();
523  if (!_cache_available){
524  _unlock_if_needed();
525  throw Error("Calls to JetMedianBackgroundEstimator::empty_area() in cases where the background estimation uses a selector that takes a reference jet need to call a method that fills the cached estimate (rho(jet), sigma(jet), ...).");
526  }
528  _unlock_if_needed();
529  return return_value;
530  }
531  if (!_cache_available) _compute_and_cache_no_overwrite();
533 }
534 
535 /// Returns the number of empty jets used when computing the
536 /// background properties. The value is that for the last call of
537 /// rho() or sigma().
538 /// If the configuration has changed in the meantime, throw an error.
539 ///
540 /// If the area has explicit ghosts the result is zero; for active
541 /// areas it is the number of internal pure ghost jets that pass the
542 /// selector; otherwise it is deduced from the empty area, divided by
543 /// \f$ 0.55 \pi R^2 \f$ (the average pure-ghost-jet area).
544 ///
545 /// The result here is just the cached result of the corresponding
546 /// call to the ClusterSequenceAreaBase function.
548  // if the selector takes a reference, we need to use the cache
549  if (_rho_range.takes_reference()){
550  // lock to make sure no other thread interferes w the caching
551  _lock_if_needed();
552  if (!_cache_available){
553  _unlock_if_needed();
554  throw Error("Calls to JetMedianBackgroundEstimator::n_empty_jets() in cases where the background estimation uses a selector that takes a reference jet need to call a method that fills the cached estimate (rho(jet), sigma(jet), ...).");
555  }
557  _unlock_if_needed();
558  return return_value;
559  }
560  if (!_cache_available) _compute_and_cache_no_overwrite();
562 }
563 
564 
565 //----------------------------------------------------------------------
566 // configuring behaviour
567 //----------------------------------------------------------------------
568 // reset to default values
569 //
570 // set the variou options to their default values
572  // set the remaining default parameters
573  set_use_area_4vector(); // true by default
574  set_provide_fj2_sigma(false);
575 
576  _enable_rho_m = true;
577 
578  // reset the computed values
579  _included_jets.clear();
580 
581  _jet_density_class = 0; // null pointer
582  _rescaling_class = 0; // null pointer
583 
584  _set_cache_unavailable();
585 }
586 
587 
588 // Set a pointer to a class that calculates the quantity whose
589 // median will be calculated; if the pointer is null then pt/area
590 // is used (as occurs also if this function is not called).
592  _warnings_preliminary.warn("JetMedianBackgroundEstimator::set_jet_density_class: density classes are still preliminary in FastJet 3.1. Their interface may differ in future releases (without guaranteeing backward compatibility). Note that since FastJet 3.1, rho_m and sigma_m are accessible direclty in JetMedianBackgroundEstimator and GridMedianBackgroundEstimator(with no need for a density class).");
593  _jet_density_class = jet_density_class_in;
594  _set_cache_unavailable();
595 }
596 
597 
598 
599 //----------------------------------------------------------------------
600 // description
601 //----------------------------------------------------------------------
603  ostringstream desc;
604  desc << "JetMedianBackgroundEstimator, using " << _jet_def.description()
605  << " with " << _area_def.description()
606  << " and selecting jets with " << _rho_range.description();
607  return desc.str();
608 }
609 
610 
611 
612 //----------------------------------------------------------------------
613 // computation of the background properties
614 //----------------------------------------------------------------------
615 // do the actual job
616 BackgroundEstimate JetMedianBackgroundEstimator::_compute(const PseudoJet &jet) const {
617  // prepare a local structure to hold temporarily the results
618  // (by design, this comes with default values of 0 for each property)
619  BackgroundEstimate local_estimate;
620 
621  // check if the clustersequence is still valid
622  _check_csa_alive();
623 
624  local_estimate.set_has_sigma(has_sigma());
625  local_estimate.set_has_rho_m(has_rho_m());
626 
627  // structure to hold the extra info associated w this BGE (the call
628  // below initialises everything to 0)
629  Extras * extras = new Extras;
630  local_estimate.set_extras(extras);
631  extras->set_reference_jet(jet);
632 
633  // fill the vector of pt/area (or the quantity from the jet density class)
634  // - in the range
635  vector<double> vector_for_median_pt;
636  vector<double> vector_for_median_dt;
637  double total_area = 0.0;
638 
639  // apply the selector to the included jets
640  vector<PseudoJet> selected_jets;
641  if (_rho_range.takes_reference()){
642  Selector local_rho_range = _rho_range;
643  selected_jets = local_rho_range.set_reference(jet)(_included_jets);
644  } else {
645  selected_jets = _rho_range(_included_jets);
646  }
647 
648  // compute the pt/area for the selected jets
649  double median_input_pt, median_input_dt=0.0;
650  BackgroundJetPtMDensity m_density;
651  bool do_rho_m = has_rho_m();
652  unsigned int njets_used = 0;
653  for (unsigned i = 0; i < selected_jets.size(); i++) {
654  const PseudoJet & current_jet = selected_jets[i];
655 
656  double this_area = (_use_area_4vector) ? current_jet.area_4vector().perp() : current_jet.area();
657  if (this_area>0){
658  // for the pt component, we either use pt or the user-provided
659  // density class
660  if (_jet_density_class == 0) {
661  median_input_pt = current_jet.perp()/this_area;
662  } else {
663  median_input_pt = (*_jet_density_class)(current_jet);
664  }
665 
666  // handle the rho_m part if requested
667  // note that we're using the scalar area as a normalisation inside the
668  // density class!
669  if (do_rho_m)
670  median_input_dt = m_density(current_jet);
671 
672  // perform rescaling if needed
673  if (_rescaling_class != 0) {
674  double resc = (*_rescaling_class)(current_jet);
675  median_input_pt /= resc;
676  median_input_dt /= resc;
677  }
678 
679  // store the result for future computation of the median
680  vector_for_median_pt.push_back(median_input_pt);
681  if (do_rho_m)
682  vector_for_median_dt.push_back(median_input_dt);
683 
684  total_area += this_area;
685  njets_used++;
686  } else {
687  _warnings_zero_area.warn("JetMedianBackgroundEstimator::_compute(...): 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).");
688  }
689  }
690 
691  // there is nothing inside our region, so answer will always be zero
692  if (vector_for_median_pt.size() == 0) {
693  // record that the computation has been performed
694  return local_estimate;
695  }
696 
697  // determine the number of empty jets
698  // If we have explicit ghosts, this is 0 (i.e. the default)
699  const ClusterSequenceAreaBase * csab = (dynamic_cast<ClusterSequenceStructure*>(_csi.get()))->validated_csab();
700  if (! (csab->has_explicit_ghosts())) {
701  if (_rho_range.takes_reference()){
702  Selector local_rho_range = _rho_range;
703  local_rho_range.set_reference(jet);
704  extras->set_empty_area (csab->empty_area(local_rho_range));
705  extras->set_n_empty_jets(csab->n_empty_jets(local_rho_range));
706  } else {
707  extras->set_empty_area (csab->empty_area(_rho_range));
708  extras->set_n_empty_jets(csab->n_empty_jets(_rho_range));
709  }
710  }
711 
712  extras->set_n_jets_used(njets_used);
713  double total_njets = extras->n_jets_used() + extras->n_empty_jets();
714  total_area += extras->empty_area();
715 
716  double rho_tmp, stand_dev;
717  _median_and_stddev(vector_for_median_pt, extras->n_empty_jets(),
718  rho_tmp, stand_dev, _provide_fj2_sigma);
719  local_estimate.set_rho(rho_tmp);
720 
721  // process and store the results (_rho was already stored above)
722  local_estimate.set_mean_area(total_area / total_njets);
723  local_estimate.set_sigma(stand_dev * sqrt(local_estimate.mean_area()));
724 
725  // compute the rho_m part now
726  if (do_rho_m){
727  _median_and_stddev(vector_for_median_dt, extras->n_empty_jets(),
728  rho_tmp, stand_dev,
729  _provide_fj2_sigma);
730  local_estimate.set_rho_m(rho_tmp);
731  local_estimate.set_sigma_m(stand_dev * sqrt(local_estimate.mean_area()));
732  }
733 
734  return local_estimate;
735 }
736 
737 
738 void JetMedianBackgroundEstimator::_cache_no_overwrite(const BackgroundEstimate &estimate) const {
739  /// this is meant to be called if the selector is not local
740  assert(!(_rho_range.takes_reference()));
741 
742  // we need to write to the cache, so set a lock if needed
743  _lock_if_needed();
744 
745  // we only need to write if someone else did not do it earlier
746  //
747  // Doing this avoids potential "read" problems when two threads try
748  // to compute the cache at the same time. The first one may return
749  // and try to read the result when the second actually writes. The
750  // lines below guarantee that the first thread would have set
751  // _cache_available to true before releasing the lock and therefore
752  // the second thread will not attempt to write
753  if (!_cache_available){
755  _cache_available = true;
756  }
757 
758  // release the lock
759  _unlock_if_needed();
760 }
761 
762 void JetMedianBackgroundEstimator::_compute_and_cache_no_overwrite() const {
763  /// this is meant to be called if the selector is not local
764  assert(!(_rho_range.takes_reference()));
765 
766  // get the result (for a dummy PseudoJet which will anyway not be used)
767  // and cache it
768  _cache_no_overwrite(_compute(PseudoJet()));
769 }
770 
771 void JetMedianBackgroundEstimator::_cache(const BackgroundEstimate &estimate) const {
772  /// this is meant to be called if the selector is local
773  assert(_rho_range.takes_reference());
774 
775  // we need to write to the cache, so set a lock if needed
776  _lock_if_needed();
777 
778  // if we overwrite the cache, we need to make sure that other
779  // parts of the code use r/w accesses that respect the lock that
780  // we have acquired.
781  //
782  // In practice, this will only be used in the case wheere we have
783  // a local selector, in queries that require access to the cache.
785  _cache_available = true;
786 
787  // release the lock
788  _unlock_if_needed();
789 }
790 
791 BackgroundEstimate JetMedianBackgroundEstimator::_compute_and_cache_if_needed(const PseudoJet &jet) const {
792  /// this is meant to be called if the selector is local
793  assert(_rho_range.takes_reference());
794 
795  BackgroundEstimate local_estimate;
796 
797  _lock_if_needed();
798  if ((_cache_available) && (_cached_estimate.extras<JetMedianBackgroundEstimator>().reference_jet() == jet)){
799  local_estimate = _cached_estimate;
800  _unlock_if_needed();
801  return local_estimate;
802  }
803  _unlock_if_needed();
804 
805  local_estimate = _compute(jet);
806  _cache(local_estimate);
807  return local_estimate;
808 }
809 
810 // check that the underlying structure is still alive;
811 // throw an error otherwise
812 void JetMedianBackgroundEstimator::_check_csa_alive() const{
813  ClusterSequenceStructure* csa = dynamic_cast<ClusterSequenceStructure*>(_csi.get());
814  if (csa == 0) {
815  throw Error("JetMedianBackgroundEstimator: there is no cluster sequence associated with the JetMedianBackgroundEstimator");
816  }
817  if (! csa->has_associated_cluster_sequence())
818  throw Error("JetMedianBackgroundEstimator: modifications are no longer possible as the underlying ClusterSequence has gone out of scope");
819 }
820 
821 
822 // check that the algorithm used for the clustering is suitable for
823 // background estimation (i.e. either kt or C/A).
824 // Issue a warning otherwise
825 void JetMedianBackgroundEstimator::_check_jet_alg_good_for_median() const{
826  const JetDefinition * jet_def = &_jet_def;
827 
828  // if no explicit jet def has been provided, fall back on the
829  // cluster sequence
830  if (_jet_def.jet_algorithm() == undefined_jet_algorithm){
831  const ClusterSequence * cs = dynamic_cast<ClusterSequenceStructure*>(_csi.get())->validated_cs();
832  jet_def = &(cs->jet_def());
833  }
834 
835  if (jet_def->jet_algorithm() != kt_algorithm
836  && jet_def->jet_algorithm() != cambridge_algorithm
837  && jet_def->jet_algorithm() != cambridge_for_passive_algorithm) {
838  _warnings.warn("JetMedianBackgroundEstimator: jet_def being used may not be suitable for estimating diffuse backgrounds (good alternatives are kt, cam)");
839  }
840 }
841 
842 
843 FASTJET_END_NAMESPACE
844 
845 
fastjet::AreaDefinition::description
std::string description() const
return a description of the current area definition
Definition: AreaDefinition.cc:48
fastjet::PseudoJet::area
virtual double area() const
return the jet (scalar) area.
Definition: PseudoJet.cc:829
fastjet::JetMedianBackgroundEstimator::mean_area
double mean_area() const
Returns the mean area of the jets used to actually compute the background properties in the last call...
Definition: JetMedianBackgroundEstimator.cc:437
fastjet::BackgroundEstimate::apply_rescaling_factor
void apply_rescaling_factor(double rescaling_factor)
apply a rescaling factor (to rho, rho_m, sigma, sigma_m)
Definition: BackgroundEstimatorBase.hh:144
fastjet::ClusterSequenceAreaBase::has_explicit_ghosts
virtual bool has_explicit_ghosts() const
returns true if ghosts are explicitly included within jets for this ClusterSequence;
Definition: ClusterSequenceAreaBase.hh:104
fastjet::JetMedianBackgroundEstimator::sigma_m
virtual double sigma_m() const override
returns sigma_m, a measure of the fluctuations in the purely longitudinal, particle-mass-induced comp...
Definition: JetMedianBackgroundEstimator.cc:385
fastjet::ClusterSequence::inclusive_jets
std::vector< PseudoJet > inclusive_jets(const double ptmin=0.0) const
return a vector of all jets (in the sense of the inclusive algorithm) with pt >= ptmin.
Definition: ClusterSequence.cc:916
fastjet::Selector::takes_reference
bool takes_reference() const
returns true if this can be applied jet by jet
Definition: Selector.hh:287
fastjet::ClusterSequence::delete_self_when_unused
void delete_self_when_unused()
by calling this routine you tell the ClusterSequence to delete itself when all the Pseudojets associa...
Definition: ClusterSequence.cc:1778
fastjet::JetMedianBackgroundEstimator::set_use_area_4vector
void set_use_area_4vector(bool use_it=true)
By default when calculating pt/Area for a jet, it is the transverse component of the 4-vector area th...
Definition: JetMedianBackgroundEstimator.hh:339
fastjet::JetMedianBackgroundEstimator::set_jets
void set_jets(const std::vector< PseudoJet > &jets)
(re)set the jets (which must have area support) to be used by future calls to rho() etc.
Definition: JetMedianBackgroundEstimator.cc:214
fastjet::JetDefinition
class that is intended to hold a full definition of the jet clusterer
Definition: JetDefinition.hh:250
fastjet::AreaDefinition::with_fixed_seed
AreaDefinition with_fixed_seed(const std::vector< int > &iseed) const
return a copy of this AreaDefinition with a user-defined set of seeds
Definition: AreaDefinition.hh:148
fastjet::BackgroundEstimate::mean_area
double mean_area() const
mean area of the patches used to compute the background properties
Definition: BackgroundEstimatorBase.hh:89
fastjet::SharedPtr::get
T * get() const
get the stored pointer
Definition: SharedPtr.hh:473
fastjet::Selector::description
std::string description() const
returns a textual description of the selector
Definition: Selector.hh:237
fastjet::PseudoJet::constituents
virtual std::vector< PseudoJet > constituents() const
retrieve the constituents.
Definition: PseudoJet.cc:692
fastjet::JetMedianBackgroundEstimator
Class to estimate the pt density of the background per unit area, using the median of the distributio...
Definition: JetMedianBackgroundEstimator.hh:86
fastjet::JetMedianBackgroundEstimator::has_sigma
virtual bool has_sigma() const override
returns true if this background estimator has support for determination of sigma
Definition: JetMedianBackgroundEstimator.hh:238
fastjet::JetMedianBackgroundEstimator::has_rho_m
virtual bool has_rho_m() const override
Returns true if this background estimator has support for determination of rho_m.
Definition: JetMedianBackgroundEstimator.hh:267
fastjet::LimitedWarning::warn
void warn(const char *warning)
outputs a warning to standard error (or the user's default warning stream if set)
Definition: LimitedWarning.hh:71
fastjet::JetMedianBackgroundEstimator::estimate
BackgroundEstimate estimate() const override
get the full set of background properties
Definition: JetMedianBackgroundEstimator.cc:270
fastjet::BackgroundEstimatorBase::_cached_estimate
BackgroundEstimate _cached_estimate
all the info about what is computed
Definition: BackgroundEstimatorBase.hh:374
fastjet::JetMedianBackgroundEstimator::empty_area
double empty_area() const
Returns the estimate of the area (within the range defined by the selector) that is not occupied by j...
Definition: JetMedianBackgroundEstimator.cc:518
fastjet::BackgroundEstimate
/// a class that holds the result of the calculation
Definition: BackgroundEstimatorBase.hh:53
fastjet::cambridge_algorithm
@ cambridge_algorithm
the longitudinally invariant variant of the cambridge algorithm (aka Aachen algoithm).
Definition: JetDefinition.hh:142
fastjet::BackgroundEstimate::extras
const BackgroundEstimator::Extras & extras() const
returns a reference to the extra information associated with a given BackgroundEstimator.
Definition: BackgroundEstimatorBase.hh:119
fastjet::SharedPtr
An implementation of shared pointers that is broadly similar to C++11 shared_ptr (https://en....
Definition: SharedPtr.hh:341
fastjet::JetDefinition::description
std::string description() const
return a textual description of the current jet definition
Definition: JetDefinition.cc:106
fastjet::BackgroundEstimatorBase::_median_and_stddev
void _median_and_stddev(const std::vector< double > &quantity_vector, double n_empty_jets, double &median, double &stand_dev_if_gaussian, bool do_fj2_calculation=false) const
given a quantity in a vector (e.g.
Definition: BackgroundEstimatorBase.cc:57
fastjet::AreaDefinition
class that holds a generic area definition
Definition: AreaDefinition.hh:82
fastjet::ClusterSequence::structure_shared_ptr
const SharedPtr< PseudoJetStructureBase > & structure_shared_ptr() const
retrieve a shared pointer to the wrapper to this ClusterSequence
Definition: ClusterSequence.hh:542
fastjet::Selector::set_reference
const Selector & set_reference(const PseudoJet &reference)
set the reference jet for this Selector
Definition: Selector.hh:292
fastjet::JetMedianBackgroundEstimator::set_particles
virtual void set_particles(const std::vector< PseudoJet > &particles) override
tell the background estimator that it has a new event, composed of the specified particles.
Definition: JetMedianBackgroundEstimator.cc:117
fastjet::ClusterSequenceArea
General class for user to obtain ClusterSequence with additional area information.
Definition: ClusterSequenceArea.hh:51
fastjet::ClusterSequenceStructure::validated_csab
virtual const ClusterSequenceAreaBase * validated_csab() const override
if the jet has valid area information then return a pointer to the associated ClusterSequenceAreaBase...
Definition: ClusterSequenceStructure.cc:269
fastjet::FunctionOfPseudoJet< double >
fastjet::BackgroundEstimate::sigma_m
double sigma_m() const
fluctuations in the purely longitudinal (particle-mass-induced) component of the background density p...
Definition: BackgroundEstimatorBase.hh:82
fastjet::JetDefinition::jet_algorithm
JetAlgorithm jet_algorithm() const
return information about the definition...
Definition: JetDefinition.hh:432
fastjet::JetMedianBackgroundEstimator::n_jets_used
unsigned int n_jets_used() const
returns the number of jets used to actually compute the background properties in the last call of rho...
Definition: JetMedianBackgroundEstimator.cc:457
fastjet::JetMedianBackgroundEstimator::n_empty_jets
double n_empty_jets() const
Returns the number of empty jets used when computing the background properties.
Definition: JetMedianBackgroundEstimator.cc:547
fastjet::JetMedianBackgroundEstimator::rho
double rho() const override
get rho, the median background density per unit area
Definition: JetMedianBackgroundEstimator.cc:303
fastjet::JetMedianBackgroundEstimator::set_jet_density_class
void set_jet_density_class(const FunctionOfPseudoJet< double > *jet_density_class)
Set a pointer to a class that calculates the quantity whose median will be calculated; if the pointer...
Definition: JetMedianBackgroundEstimator.cc:591
fastjet::JetMedianBackgroundEstimator::reset
void reset()
Resets the class to its default state, including the choice to use 4-vector areas.
Definition: JetMedianBackgroundEstimator.cc:571
fastjet::BackgroundEstimate::rho_m
double rho_m() const
purely longitudinal (particle-mass-induced) component of the background density per unit area
Definition: BackgroundEstimatorBase.hh:78
fastjet::PseudoJet
Class to contain pseudojets, including minimal information of use to jet-clustering routines.
Definition: PseudoJet.hh:68
fastjet::ClusterSequenceAreaBase
base class that sets interface for extensions of ClusterSequence that provide information about the a...
Definition: ClusterSequenceAreaBase.hh:48
fastjet::JetMedianBackgroundEstimator::description
std::string description() const override
returns a textual description of the background estimator
Definition: JetMedianBackgroundEstimator.cc:602
fastjet::PseudoJet::rap
double rap() const
returns the rapidity or some large value when the rapidity is infinite
Definition: PseudoJet.hh:138
fastjet::JetMedianBackgroundEstimator::jets_used
std::vector< PseudoJet > jets_used() const
returns the jets used to actually compute the background properties
Definition: JetMedianBackgroundEstimator.cc:476
fastjet::Selector
Class that encodes information about cuts and other selection criteria that can be applied to PseudoJ...
Definition: Selector.hh:149
fastjet::LimitedWarning
class to provide facilities for giving warnings up to some maximum number of times and to provide glo...
Definition: LimitedWarning.hh:54
fastjet::BackgroundEstimate::set_extras
void set_extras(Extras *extras_in)
sets the extra info based on the provided pointer
Definition: BackgroundEstimatorBase.hh:155
fastjet::undefined_jet_algorithm
@ undefined_jet_algorithm
the value for the jet algorithm in a JetDefinition for which no algorithm has yet been defined
Definition: JetDefinition.hh:172
fastjet::Selector::has_finite_area
bool has_finite_area() const
returns true if it has a meaningful and finite area (i.e.
Definition: Selector.hh:250
fastjet::BackgroundEstimate::rho
double rho() const
background density per unit area
Definition: BackgroundEstimatorBase.hh:66
fastjet::SelectorIsPureGhost
Selector SelectorIsPureGhost()
select objects that are (or are only made of) ghosts.
Definition: Selector.cc:1420
fastjet::JetMedianBackgroundEstimator::JetMedianBackgroundEstimator
JetMedianBackgroundEstimator(const Selector &rho_range, const JetDefinition &jet_def, const AreaDefinition &area_def)
Constructor that sets the rho range as well as the jet definition and area definition to be used to c...
Definition: JetMedianBackgroundEstimator.cc:84
fastjet::JetMedianBackgroundEstimator::set_cluster_sequence
void set_cluster_sequence(const ClusterSequenceAreaBase &csa)
(re)set the cluster sequence (with area support) to be used by future calls to rho() etc.
Definition: JetMedianBackgroundEstimator.cc:183
fastjet::JetMedianBackgroundEstimator::sigma
double sigma() const override
get sigma, the background fluctuations per unit area
Definition: JetMedianBackgroundEstimator.cc:314
fastjet::JetMedianBackgroundEstimator::rho_m
virtual double rho_m() const override
returns rho_m, the purely longitudinal, particle-mass-induced component of the background density per...
Definition: JetMedianBackgroundEstimator.cc:369
fastjet::ClusterSequenceStructure
Contains any information related to the clustering that should be directly accessible to PseudoJet.
Definition: ClusterSequenceStructure.hh:61
fastjet::JetMedianBackgroundEstimator::set_provide_fj2_sigma
void set_provide_fj2_sigma(bool provide_fj2_sigma=true)
The FastJet v2.X sigma calculation had a small spurious offset in the limit of a small number of jets...
Definition: JetMedianBackgroundEstimator.hh:351
fastjet::Error
base class corresponding to errors that can be thrown by FastJet
Definition: Error.hh:51
fastjet::JetMedianBackgroundEstimator::set_particles_with_seed
virtual void set_particles_with_seed(const std::vector< PseudoJet > &particles, const std::vector< int > &seed) override
an alternative call that takes a random number generator seed (typically a vector of length 2) to ens...
Definition: JetMedianBackgroundEstimator.cc:128
fastjet::cambridge_for_passive_algorithm
@ cambridge_for_passive_algorithm
a version of cambridge with a special distance measure for particles whose pt is < extra_param(); thi...
Definition: JetDefinition.hh:156
fastjet::BackgroundEstimate::sigma
double sigma() const
background fluctuations per unit square-root area must be multipled by sqrt(area) to get fluctuations...
Definition: BackgroundEstimatorBase.hh:71
fastjet::kt_algorithm
@ kt_algorithm
the longitudinally invariant kt algorithm
Definition: JetDefinition.hh:139