1#ifndef __FASTJET_BACKGROUND_ESTIMATOR_HH__
2#define __FASTJET_BACKGROUND_ESTIMATOR_HH__
35#include "fastjet/config.h"
36#include "fastjet/ClusterSequenceAreaBase.hh"
37#include "fastjet/AreaDefinition.hh"
38#include "fastjet/FunctionOfPseudoJet.hh"
39#include "fastjet/Selector.hh"
40#include "fastjet/tools/BackgroundEstimatorBase.hh"
43#ifdef FASTJET_HAVE_THREAD_SAFETY
47FASTJET_BEGIN_NAMESPACE
132 _enable_rho_m(true){ reset(); }
147 virtual void set_particles(
const std::vector<PseudoJet> & particles) FASTJET_OVERRIDE;
152 virtual void set_particles_with_seed(
const std::vector<PseudoJet> & particles,
const std::vector<int> & seed) FASTJET_OVERRIDE;
176 void set_jets(
const std::vector<PseudoJet> &jets);
180 _rho_range = rho_range_selector;
181 _set_cache_unavailable();
187 _enable_rho_m = enable;
188 _set_cache_unavailable();
215 double rho() const FASTJET_OVERRIDE;
218 double sigma() const FASTJET_OVERRIDE;
226 double rho(const
PseudoJet & jet) FASTJET_OVERRIDE;
234 double sigma(const
PseudoJet &jet) FASTJET_OVERRIDE;
238 virtual
bool has_sigma() const FASTJET_OVERRIDE {
return true;}
245 virtual double rho_m() const FASTJET_OVERRIDE;
251 virtual
double sigma_m() const FASTJET_OVERRIDE;
254 virtual
double rho_m(const
PseudoJet & ) FASTJET_OVERRIDE;
257 virtual
double sigma_m(const
PseudoJet & ) FASTJET_OVERRIDE;
267 virtual
bool has_rho_m() const FASTJET_OVERRIDE {
return _enable_rho_m && (_jet_density_class == 0);}
276 double mean_area()
const;
281 unsigned int n_jets_used()
const;
285 std::vector<PseudoJet> jets_used()
const;
301 double empty_area()
const;
315 double n_empty_jets()
const;
340 _use_area_4vector = use_it;
341 _set_cache_unavailable();
352 _provide_fj2_sigma = provide_fj2_sigma;
353 _set_cache_unavailable();
367 return _jet_density_class;
379 BackgroundEstimatorBase::set_rescaling_class(rescaling_class_in);
380 _set_cache_unavailable();
390 std::string description() const FASTJET_OVERRIDE;
399 : _reference_jet(
PseudoJet()), _n_jets_used(0),
400 _n_empty_jets(0.0), _empty_area(0.0) {}
414 void set_reference_jet(
const PseudoJet &reference_jet_in){
415 _reference_jet = reference_jet_in;
417 void set_n_jets_used(
int n_jets_used_in){ _n_jets_used=n_jets_used_in;}
418 void set_n_empty_jets(
double n_empty_jets_in){ _n_empty_jets=n_empty_jets_in;}
419 void set_empty_area(
double empty_area_in){ _empty_area=empty_area_in;}
447 void _compute_and_cache_no_overwrite()
const;
465 void _check_csa_alive()
const;
470 void _check_jet_alg_good_for_median()
const;
476 std::vector<PseudoJet> _included_jets;
479 bool _use_area_4vector;
480 bool _provide_fj2_sigma;
507 virtual std::string
description()
const {
return "BackgroundJetPtDensity";}
528 virtual double result(
const PseudoJet & jet)
const;
530 virtual std::string description()
const;
549 std::vector<PseudoJet> constituents = jet.
constituents();
550 double scalar_ptm = 0;
551 for (
unsigned i = 0; i < constituents.size(); i++) {
552 scalar_ptm += constituents[i].mperp() - constituents[i].perp();
554 return scalar_ptm / jet.
area();
557 virtual std::string
description()
const {
return "BackgroundPtMDensity";}
class that holds a generic area definition
/// a class that holds the result of the calculation
Abstract base class that provides the basic interface for classes that estimate levels of background ...
Class that implements pt/area_4vector.perp() for background estimation (this is a preliminary class).
virtual std::string description() const
returns a description of the function (an empty string by default)
virtual double result(const PseudoJet &jet) const
the action of the function this has to be overloaded in derived classes
Class that implements for background estimation (this is a preliminary class).
virtual std::string description() const
returns a description of the function (an empty string by default)
virtual double result(const PseudoJet &jet) const
the action of the function this has to be overloaded in derived classes
Class that implements (scalar pt sum of jet)/(scalar area of jet) for background estimation (this is ...
BackgroundJetScalarPtDensity(double n)
Constructor to provide background estimation based on .
BackgroundJetScalarPtDensity()
Default constructor provides background estimation with scalar pt sum.
base class that sets interface for extensions of ClusterSequence that provide information about the a...
class that is intended to hold a full definition of the jet clusterer
class to provide facilities for giving warnings up to some maximum number of times and to provide glo...
Class to contain pseudojets, including minimal information of use to jet-clustering routines.
virtual std::vector< PseudoJet > constituents() const
retrieve the constituents.
double perp() const
returns the scalar transverse momentum
virtual double area() const
return the jet (scalar) area.
virtual PseudoJet area_4vector() const
return the jet 4-vector area.
Class that encodes information about cuts and other selection criteria that can be applied to PseudoJ...
An implementation of shared pointers that is broadly similar to C++11 shared_ptr (https://en....