FastJet 3.4.1
List of all members
fastjet::BackgroundEstimatorBase Class Referenceabstract

Abstract base class that provides the basic interface for classes that estimate levels of background radiation in hadron and heavy-ion collider events. More...

#include <fastjet/tools/BackgroundEstimatorBase.hh>

Inheritance diagram for fastjet::BackgroundEstimatorBase:
Inheritance graph
[legend]
Collaboration diagram for fastjet::BackgroundEstimatorBase:
Collaboration graph
[legend]

Public Member Functions

constructors and destructors
 BackgroundEstimatorBase ()
 
virtual ~BackgroundEstimatorBase ()
 a default virtual destructor that does nothing More...
 
setting a new event
virtual void set_particles (const std::vector< PseudoJet > &particles)=0
 tell the background estimator that it has a new event, composed of the specified particles. More...
 
virtual void set_particles_with_seed (const std::vector< PseudoJet > &particles, const std::vector< int > &)
 an alternative call that takes a random number generator seed (typically a vector of length 2) to ensure reproducibility of background estimators that rely on random numbers (specifically JetMedianBackgroundEstimator with ghosted areas) More...
 
virtual BackgroundEstimatorBasecopy () const =0
 return a pointer to a copy of this BGE; the user is responsible for eventually deleting the resulting object. More...
 
retrieving fundamental information
virtual BackgroundEstimate estimate () const =0
 get the full set of background properties More...
 
virtual BackgroundEstimate estimate (const PseudoJet &jet) const =0
 get the full set of background properties for a given reference jet More...
 
virtual double rho () const =0
 get rho, the background density per unit area More...
 
virtual double sigma () const
 get sigma, the background fluctuations per unit square-root area; must be multipled by sqrt(area) to get fluctuations for a region of a given area. More...
 
virtual double rho (const PseudoJet &jet)=0
 get rho, the background density per unit area, locally at the position of a given jet. More...
 
virtual double sigma (const PseudoJet &)
 get sigma, the background fluctuations per unit area, locally at the position of a given jet. More...
 
virtual bool has_sigma () const
 returns true if this background estimator has support for determination of sigma More...
 
virtual double rho_m () const
 returns rho_m, the purely longitudinal, particle-mass-induced component of the background density per unit area More...
 
virtual double sigma_m () const
 returns sigma_m, a measure of the fluctuations in the purely longitudinal, particle-mass-induced component of the background density per unit square-root area; must be multipled by sqrt(area) to get fluctuations for a region of a given area. More...
 
virtual double rho_m (const PseudoJet &)
 Returns rho_m locally at the jet position. As for rho(jet), it is non-const. More...
 
virtual double sigma_m (const PseudoJet &)
 Returns sigma_m locally at the jet position. As for rho(jet), it is non-const. More...
 
virtual bool has_rho_m () const
 Returns true if this background estimator has support for determination of rho_m. More...
 
configuring the behaviour
virtual void set_rescaling_class (const FunctionOfPseudoJet< double > *rescaling_class_in)
 Set a pointer to a class that calculates the rescaling factor as a function of the jet (position). More...
 
const FunctionOfPseudoJet< double > * rescaling_class () const
 return the pointer to the jet density class More...
 
description
virtual std::string description () const =0
 returns a textual description of the background estimator More...
 

Protected Member Functions

helpers for thread safety

Note that these helpers are related to median-based estimation of the background, so there is no guarantee that they will remain in this base class in the long term

void _set_cache_unavailable ()
 
void _lock_if_needed () const
 
void _unlock_if_needed () const
 

helpers for derived classes

Note that these helpers are related to median-based estimation of the background, so there is no guarantee that they will remain in this base class in the long term

const FunctionOfPseudoJet< double > * _rescaling_class
 
bool _cache_available
 
BackgroundEstimate _cached_estimate
 all the info about what is computed More...
 
static LimitedWarning _warnings_empty_area
 
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. More...
 
double _percentile (const std::vector< double > &sorted_quantity_vector, const double percentile, const double nempty=0.0, const bool do_fj2_calculation=false) const
 computes a percentile of a given sorted vector More...
 

Detailed Description

Abstract base class that provides the basic interface for classes that estimate levels of background radiation in hadron and heavy-ion collider events.

Definition at line 183 of file BackgroundEstimatorBase.hh.

Constructor & Destructor Documentation

◆ BackgroundEstimatorBase()

fastjet::BackgroundEstimatorBase::BackgroundEstimatorBase ( )
inline

Definition at line 188 of file BackgroundEstimatorBase.hh.

◆ ~BackgroundEstimatorBase()

virtual fastjet::BackgroundEstimatorBase::~BackgroundEstimatorBase ( )
inlinevirtual

a default virtual destructor that does nothing

Definition at line 199 of file BackgroundEstimatorBase.hh.

Member Function Documentation

◆ set_particles()

virtual void fastjet::BackgroundEstimatorBase::set_particles ( const std::vector< PseudoJet > &  particles)
pure virtual

tell the background estimator that it has a new event, composed of the specified particles.

Implemented in fastjet::GridMedianBackgroundEstimator, and fastjet::JetMedianBackgroundEstimator.

◆ set_particles_with_seed()

virtual void fastjet::BackgroundEstimatorBase::set_particles_with_seed ( const std::vector< PseudoJet > &  particles,
const std::vector< int > &   
)
inlinevirtual

an alternative call that takes a random number generator seed (typically a vector of length 2) to ensure reproducibility of background estimators that rely on random numbers (specifically JetMedianBackgroundEstimator with ghosted areas)

Reimplemented in fastjet::JetMedianBackgroundEstimator.

Definition at line 214 of file BackgroundEstimatorBase.hh.

◆ copy()

virtual BackgroundEstimatorBase * fastjet::BackgroundEstimatorBase::copy ( ) const
pure virtual

return a pointer to a copy of this BGE; the user is responsible for eventually deleting the resulting object.

Implemented in fastjet::GridMedianBackgroundEstimator, and fastjet::JetMedianBackgroundEstimator.

◆ estimate() [1/2]

virtual BackgroundEstimate fastjet::BackgroundEstimatorBase::estimate ( ) const
pure virtual

get the full set of background properties

Implemented in fastjet::GridMedianBackgroundEstimator, and fastjet::JetMedianBackgroundEstimator.

◆ estimate() [2/2]

virtual BackgroundEstimate fastjet::BackgroundEstimatorBase::estimate ( const PseudoJet jet) const
pure virtual

get the full set of background properties for a given reference jet

Implemented in fastjet::GridMedianBackgroundEstimator, and fastjet::JetMedianBackgroundEstimator.

◆ rho() [1/2]

virtual double fastjet::BackgroundEstimatorBase::rho ( ) const
pure virtual

get rho, the background density per unit area

Implemented in fastjet::GridMedianBackgroundEstimator, and fastjet::JetMedianBackgroundEstimator.

◆ sigma() [1/2]

virtual double fastjet::BackgroundEstimatorBase::sigma ( ) const
inlinevirtual

get sigma, the background fluctuations per unit square-root area; must be multipled by sqrt(area) to get fluctuations for a region of a given area.

Reimplemented in fastjet::GridMedianBackgroundEstimator, and fastjet::JetMedianBackgroundEstimator.

Definition at line 239 of file BackgroundEstimatorBase.hh.

◆ rho() [2/2]

virtual double fastjet::BackgroundEstimatorBase::rho ( const PseudoJet jet)
pure virtual

get rho, the background density per unit area, locally at the position of a given jet.

Note that this is not const, because a user may then wish to query other aspects of the background that could depend on the position of the jet last used for a rho(jet) determination.

Implemented in fastjet::GridMedianBackgroundEstimator, and fastjet::JetMedianBackgroundEstimator.

◆ sigma() [2/2]

virtual double fastjet::BackgroundEstimatorBase::sigma ( const PseudoJet )
inlinevirtual

get sigma, the background fluctuations per unit area, locally at the position of a given jet.

As for rho(jet), it is non-const.

Reimplemented in fastjet::GridMedianBackgroundEstimator, and fastjet::JetMedianBackgroundEstimator.

Definition at line 252 of file BackgroundEstimatorBase.hh.

◆ has_sigma()

virtual bool fastjet::BackgroundEstimatorBase::has_sigma ( ) const
inlinevirtual

returns true if this background estimator has support for determination of sigma

Reimplemented in fastjet::GridMedianBackgroundEstimator, and fastjet::JetMedianBackgroundEstimator.

Definition at line 258 of file BackgroundEstimatorBase.hh.

◆ rho_m() [1/2]

virtual double fastjet::BackgroundEstimatorBase::rho_m ( ) const
inlinevirtual

returns rho_m, the purely longitudinal, particle-mass-induced component of the background density per unit area

Reimplemented in fastjet::GridMedianBackgroundEstimator, and fastjet::JetMedianBackgroundEstimator.

Definition at line 265 of file BackgroundEstimatorBase.hh.

◆ sigma_m() [1/2]

virtual double fastjet::BackgroundEstimatorBase::sigma_m ( ) const
inlinevirtual

returns sigma_m, a measure of the fluctuations in the purely longitudinal, particle-mass-induced component of the background density per unit square-root area; must be multipled by sqrt(area) to get fluctuations for a region of a given area.

Reimplemented in fastjet::GridMedianBackgroundEstimator, and fastjet::JetMedianBackgroundEstimator.

Definition at line 273 of file BackgroundEstimatorBase.hh.

◆ rho_m() [2/2]

virtual double fastjet::BackgroundEstimatorBase::rho_m ( const PseudoJet )
inlinevirtual

Returns rho_m locally at the jet position. As for rho(jet), it is non-const.

Reimplemented in fastjet::JetMedianBackgroundEstimator, and fastjet::GridMedianBackgroundEstimator.

Definition at line 278 of file BackgroundEstimatorBase.hh.

◆ sigma_m() [2/2]

virtual double fastjet::BackgroundEstimatorBase::sigma_m ( const PseudoJet )
inlinevirtual

Returns sigma_m locally at the jet position. As for rho(jet), it is non-const.

Reimplemented in fastjet::JetMedianBackgroundEstimator, and fastjet::GridMedianBackgroundEstimator.

Definition at line 283 of file BackgroundEstimatorBase.hh.

◆ has_rho_m()

virtual bool fastjet::BackgroundEstimatorBase::has_rho_m ( ) const
inlinevirtual

Returns true if this background estimator has support for determination of rho_m.

Note that support for sigma_m is automatic is one has sigma and rho_m support.

Reimplemented in fastjet::GridMedianBackgroundEstimator, and fastjet::JetMedianBackgroundEstimator.

Definition at line 292 of file BackgroundEstimatorBase.hh.

◆ set_rescaling_class()

virtual void fastjet::BackgroundEstimatorBase::set_rescaling_class ( const FunctionOfPseudoJet< double > *  rescaling_class_in)
inlinevirtual

Set a pointer to a class that calculates the rescaling factor as a function of the jet (position).

Note that the rescaling factor is used both in the determination of the "global" rho (the pt/A of each jet is divided by this factor) and when asking for a local rho (the result is multiplied by this factor).

The BackgroundRescalingYPolynomial class can be used to get a rescaling that depends just on rapidity.

There is currently no support for different rescaling classes for rho and rho_m determinations.

Reimplemented in fastjet::GridMedianBackgroundEstimator, and fastjet::JetMedianBackgroundEstimator.

Definition at line 311 of file BackgroundEstimatorBase.hh.

◆ rescaling_class()

const FunctionOfPseudoJet< double > * fastjet::BackgroundEstimatorBase::rescaling_class ( ) const
inline

return the pointer to the jet density class

Definition at line 314 of file BackgroundEstimatorBase.hh.

◆ description()

virtual std::string fastjet::BackgroundEstimatorBase::description ( ) const
pure virtual

returns a textual description of the background estimator

Implemented in fastjet::GridMedianBackgroundEstimator, and fastjet::JetMedianBackgroundEstimator.

◆ _median_and_stddev()

void fastjet::BackgroundEstimatorBase::_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
protected

given a quantity in a vector (e.g.

pt_over_area) and knowledge about the number of empty jets, calculate the median and stand_dev_if_gaussian (roughly from the 16th percentile)

If do_fj2_calculation is set to true then this performs FastJet 2.X estimation of the standard deviation, which has a spurious offset in the limit of a small number of jets.

Definition at line 57 of file BackgroundEstimatorBase.cc.

◆ _percentile()

double fastjet::BackgroundEstimatorBase::_percentile ( const std::vector< double > &  sorted_quantity_vector,
const double  percentile,
const double  nempty = 0.0,
const bool  do_fj2_calculation = false 
) const
protected

computes a percentile of a given sorted vector

Parameters
sorted_quantity_vectorthe vector contains the data sample
percentilethe percentile (defined between 0 and 1) to compute
nemptyan additional number of 0's (considered at the beginning of the quantity vector)
do_fj2_calculationcarry out the calculation as it was done in fj2 (suffers from "edge effects")

Definition at line 105 of file BackgroundEstimatorBase.cc.

◆ _set_cache_unavailable()

void fastjet::BackgroundEstimatorBase::_set_cache_unavailable ( )
inlineprotected

Definition at line 395 of file BackgroundEstimatorBase.hh.

◆ _lock_if_needed()

void fastjet::BackgroundEstimatorBase::_lock_if_needed ( ) const
protected

Definition at line 152 of file BackgroundEstimatorBase.cc.

◆ _unlock_if_needed()

void fastjet::BackgroundEstimatorBase::_unlock_if_needed ( ) const
protected

Definition at line 164 of file BackgroundEstimatorBase.cc.

Member Data Documentation

◆ _rescaling_class

const FunctionOfPseudoJet<double>* fastjet::BackgroundEstimatorBase::_rescaling_class
protected

Definition at line 368 of file BackgroundEstimatorBase.hh.

◆ _warnings_empty_area

LimitedWarning fastjet::BackgroundEstimatorBase::_warnings_empty_area
staticprotected

Definition at line 369 of file BackgroundEstimatorBase.hh.

◆ _cache_available

bool fastjet::BackgroundEstimatorBase::_cache_available
mutableprotected

Definition at line 373 of file BackgroundEstimatorBase.hh.

◆ _cached_estimate

BackgroundEstimate fastjet::BackgroundEstimatorBase::_cached_estimate
mutableprotected

all the info about what is computed

Definition at line 374 of file BackgroundEstimatorBase.hh.


The documentation for this class was generated from the following files: