FastJet 3.4.1
|
Class to estimate the pt density of the background per unit area, using the median of the distribution of pt/area from jets that pass some selection criterion. More...
#include <fastjet/tools/JetMedianBackgroundEstimator.hh>
Classes | |
class | Extras |
an internal class to hold results of the calculation that are to be assigned to the "extras" part of a BackgroundEstimate More... | |
Public Member Functions | |
constructors and destructors | |
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 cluster the particles. More... | |
JetMedianBackgroundEstimator (const Selector &rho_range, const ClusterSequenceAreaBase &csa) | |
ctor from a ClusterSequenceAreaBase with area More... | |
JetMedianBackgroundEstimator (const Selector &rho_range=SelectorIdentity()) | |
Default constructor that optionally sets the rho range. More... | |
~JetMedianBackgroundEstimator () | |
default dtor More... | |
setting a new event | |
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. More... | |
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 ensure reproducibility of background estimators that rely on random numbers (specifically JetMedianBackgroundEstimator with ghosted areas) More... | |
void | set_cluster_sequence (const ClusterSequenceAreaBase &csa) |
(re)set the cluster sequence (with area support) to be used by future calls to rho() etc. More... | |
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. More... | |
void | set_selector (const Selector &rho_range_selector) |
(re)set the selector to be used for future calls to rho() etc. More... | |
void | set_compute_rho_m (bool enable) |
determine whether the automatic calculation of rho_m and sigma_m is enabled (by default true) More... | |
BackgroundEstimatorBase * | copy () const override |
return a pointer to a copy of this BGE; the user is responsible for eventually deleting the resulting object. More... | |
retrieving fundamental information | |
BackgroundEstimate | estimate () const override |
get the full set of background properties More... | |
BackgroundEstimate | estimate (const PseudoJet &jet) const override |
get the full set of background properties for a given reference jet More... | |
double | rho () const override |
get rho, the median background density per unit area More... | |
double | sigma () const override |
get sigma, the background fluctuations per unit area More... | |
double | rho (const PseudoJet &jet) override |
get rho, the median background density per unit area, locally at the position of a given jet. More... | |
double | sigma (const PseudoJet &jet) override |
get sigma, the background fluctuations per unit area, locally at the position of a given jet. More... | |
virtual bool | has_sigma () const override |
returns true if this background estimator has support for determination of sigma More... | |
virtual double | rho_m () const override |
returns rho_m, the purely longitudinal, particle-mass-induced component of the background density per unit area More... | |
virtual double | sigma_m () const override |
returns sigma_m, a measure of the fluctuations in the purely longitudinal, particle-mass-induced component of the background density per unit area; must be multipled by sqrt(area) to get fluctuations for a region of a given area. More... | |
virtual double | rho_m (const PseudoJet &) override |
Returns rho_m locally at the jet position. As for rho(jet), it is non-const. More... | |
virtual double | sigma_m (const PseudoJet &) override |
Returns sigma_m locally at the jet position. As for rho(jet), it is non-const. More... | |
virtual bool | has_rho_m () const override |
Returns true if this background estimator has support for determination of rho_m. More... | |
retrieving additional useful information | |
double | mean_area () const |
Returns the mean area of the jets used to actually compute the background properties in the last call of rho() or sigma() If the configuration has changed in the meantime, throw an error. More... | |
unsigned int | n_jets_used () const |
returns the number of jets used to actually compute the background properties in the last call of rho() or sigma() If the configuration has changed in the meantime, throw an error. More... | |
std::vector< PseudoJet > | jets_used () const |
returns the jets used to actually compute the background properties More... | |
double | empty_area () const |
Returns the estimate of the area (within the range defined by the selector) that is not occupied by jets. More... | |
double | n_empty_jets () const |
Returns the number of empty jets used when computing the background properties. More... | |
configuring behaviour | |
void | reset () |
Resets the class to its default state, including the choice to use 4-vector areas. More... | |
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 that is used in the ratiof . More... | |
bool | use_area_4vector () const |
check if the estimator uses the 4-vector area or the scalar area More... | |
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. More... | |
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 is null then pt/area is used (as occurs also if this function is not called). More... | |
const FunctionOfPseudoJet< double > * | jet_density_class () const |
return the pointer to the jet density class More... | |
virtual void | set_rescaling_class (const FunctionOfPseudoJet< double > *rescaling_class_in) override |
Set a pointer to a class that calculates the rescaling factor as a function of the jet (position). More... | |
Public Member Functions inherited from fastjet::BackgroundEstimatorBase | |
BackgroundEstimatorBase () | |
virtual | ~BackgroundEstimatorBase () |
a default virtual destructor that does nothing More... | |
const FunctionOfPseudoJet< double > * | rescaling_class () const |
return the pointer to the jet density class More... | |
description | |
std::string | description () const override |
returns a textual description of the background estimator More... | |
Additional Inherited Members | |
Protected Member Functions inherited from fastjet::BackgroundEstimatorBase | |
void | _set_cache_unavailable () |
void | _lock_if_needed () const |
void | _unlock_if_needed () const |
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... | |
Protected Attributes inherited from fastjet::BackgroundEstimatorBase | |
const FunctionOfPseudoJet< double > * | _rescaling_class |
bool | _cache_available |
BackgroundEstimate | _cached_estimate |
all the info about what is computed More... | |
Static Protected Attributes inherited from fastjet::BackgroundEstimatorBase | |
static LimitedWarning | _warnings_empty_area |
Class to estimate the pt density of the background per unit area, using the median of the distribution of pt/area from jets that pass some selection criterion.
Events are passed either in the form of the event particles (in which they're clustered by the class), a ClusterSequenceArea (in which case the jets used are those returned by "inclusive_jets()") or directly as a set of jets.
The selection criterion is typically a geometrical one (e.g. all jets with |y|<2) sometimes supplemented with some kinematical restriction (e.g. exclusion of the two hardest jets). It is passed to the class through a Selector.
Beware: by default, to correctly handle partially empty events, the class attempts to calculate an "empty area", based (schematically) on
range.total_area() - sum_{jets_in_range} jets.area()
For ranges with small areas, this can be inaccurate (particularly relevant in dense events where empty_area should be zero and ends up not being zero).
This calculation of empty area can be avoided if a ClusterSequenceArea class with explicit ghosts (ActiveAreaExplicitGhosts) is used. This is recommended unless speed requirements cause you to use Voronoi areas. For speedy background estimation you could also consider using GridMedianBackgroundEstimator.
Definition at line 86 of file JetMedianBackgroundEstimator.hh.
fastjet::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 cluster the particles.
Prior to the estimation of rho, one has to provide the particles to cluster using set_particles(...)
rho_range | the Selector specifying which jets will be considered |
jet_def | the jet definition to use for the clustering |
area_def | the area definition to use for the clustering |
Definition at line 84 of file JetMedianBackgroundEstimator.cc.
fastjet::JetMedianBackgroundEstimator::JetMedianBackgroundEstimator | ( | const Selector & | rho_range, |
const ClusterSequenceAreaBase & | csa | ||
) |
ctor from a ClusterSequenceAreaBase with area
rho_range | the Selector specifying which jets will be considered |
csa | the ClusterSequenceArea to use |
Pre-conditions:
Note that selectors with e.g. hardest-jets exclusion do not have a well-defined area. For this reasons, it is STRONGLY advised to use an area with explicit ghosts.
Definition at line 101 of file JetMedianBackgroundEstimator.cc.
|
inline |
Default constructor that optionally sets the rho range.
The configuration must be done later calling set_cluster_sequence(...) or set_jets(...).
rho_range | the Selector specifying which jets will be considered |
Definition at line 130 of file JetMedianBackgroundEstimator.hh.
|
inline |
default dtor
Definition at line 136 of file JetMedianBackgroundEstimator.hh.
|
overridevirtual |
tell the background estimator that it has a new event, composed of the specified particles.
Implements fastjet::BackgroundEstimatorBase.
Definition at line 117 of file JetMedianBackgroundEstimator.cc.
|
overridevirtual |
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 from fastjet::BackgroundEstimatorBase.
Definition at line 128 of file JetMedianBackgroundEstimator.cc.
void fastjet::JetMedianBackgroundEstimator::set_cluster_sequence | ( | const ClusterSequenceAreaBase & | csa | ) |
(re)set the cluster sequence (with area support) to be used by future calls to rho() etc.
csa | the cluster sequence area |
Pre-conditions:
Note that selectors with e.g. hardest-jets exclusion do not have a well-defined area. For this reasons, it is STRONGLY advised to use an area with explicit ghosts.
Definition at line 183 of file JetMedianBackgroundEstimator.cc.
void fastjet::JetMedianBackgroundEstimator::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.
; for the conditions that must be satisfied by the jets, see the Constructor that takes jets.
Definition at line 214 of file JetMedianBackgroundEstimator.cc.
|
inline |
(re)set the selector to be used for future calls to rho() etc.
Definition at line 179 of file JetMedianBackgroundEstimator.hh.
|
inline |
determine whether the automatic calculation of rho_m and sigma_m is enabled (by default true)
Definition at line 186 of file JetMedianBackgroundEstimator.hh.
|
inlineoverridevirtual |
return a pointer to a copy of this BGE; the user is responsible for eventually deleting the resulting object.
Implements fastjet::BackgroundEstimatorBase.
Definition at line 196 of file JetMedianBackgroundEstimator.hh.
|
overridevirtual |
get the full set of background properties
For background estimators using a local ranges, this throws an error (use estimate(jet) instead) In the presence of a rescaling, the rescaling is not included
Implements fastjet::BackgroundEstimatorBase.
Definition at line 270 of file JetMedianBackgroundEstimator.cc.
|
overridevirtual |
get the full set of background properties for a given reference jet
Implements fastjet::BackgroundEstimatorBase.
Definition at line 280 of file JetMedianBackgroundEstimator.cc.
|
overridevirtual |
get rho, the median background density per unit area
Implements fastjet::BackgroundEstimatorBase.
Definition at line 303 of file JetMedianBackgroundEstimator.cc.
|
overridevirtual |
get sigma, the background fluctuations per unit area
Reimplemented from fastjet::BackgroundEstimatorBase.
Definition at line 314 of file JetMedianBackgroundEstimator.cc.
|
overridevirtual |
get rho, the median background density per unit area, locally at the position of a given jet.
If the Selector associated with the range takes a reference jet (i.e. is relocatable), then for subsequent operations the Selector has that jet set as its reference.
Implements fastjet::BackgroundEstimatorBase.
Definition at line 327 of file JetMedianBackgroundEstimator.cc.
|
overridevirtual |
get sigma, the background fluctuations per unit area, locally at the position of a given jet.
If the Selector associated with the range takes a reference jet (i.e. is relocatable), then for subsequent operations the Selector has that jet set as its reference.
Reimplemented from fastjet::BackgroundEstimatorBase.
Definition at line 351 of file JetMedianBackgroundEstimator.cc.
|
inlineoverridevirtual |
returns true if this background estimator has support for determination of sigma
Reimplemented from fastjet::BackgroundEstimatorBase.
Definition at line 238 of file JetMedianBackgroundEstimator.hh.
|
overridevirtual |
returns rho_m, the purely longitudinal, particle-mass-induced component of the background density per unit area
Reimplemented from fastjet::BackgroundEstimatorBase.
Definition at line 369 of file JetMedianBackgroundEstimator.cc.
|
overridevirtual |
returns sigma_m, a measure of the fluctuations in the purely longitudinal, particle-mass-induced component of the background density per unit area; must be multipled by sqrt(area) to get fluctuations for a region of a given area.
Reimplemented from fastjet::BackgroundEstimatorBase.
Definition at line 385 of file JetMedianBackgroundEstimator.cc.
|
overridevirtual |
Returns rho_m locally at the jet position. As for rho(jet), it is non-const.
Reimplemented from fastjet::BackgroundEstimatorBase.
Definition at line 398 of file JetMedianBackgroundEstimator.cc.
|
overridevirtual |
Returns sigma_m locally at the jet position. As for rho(jet), it is non-const.
Reimplemented from fastjet::BackgroundEstimatorBase.
Definition at line 417 of file JetMedianBackgroundEstimator.cc.
|
inlineoverridevirtual |
Returns true if this background estimator has support for determination of rho_m.
In te presence of a density class, support for rho_m is automatically disabled
Note that support for sigma_m is automatic is one has sigma and rho_m support.
Reimplemented from fastjet::BackgroundEstimatorBase.
Definition at line 267 of file JetMedianBackgroundEstimator.hh.
double fastjet::JetMedianBackgroundEstimator::mean_area | ( | ) | const |
Returns the mean area of the jets used to actually compute the background properties in the last call of rho() or sigma() If the configuration has changed in the meantime, throw an error.
Definition at line 437 of file JetMedianBackgroundEstimator.cc.
unsigned int fastjet::JetMedianBackgroundEstimator::n_jets_used | ( | ) | const |
returns the number of jets used to actually compute the background properties in the last call of rho() or sigma() If the configuration has changed in the meantime, throw an error.
Definition at line 457 of file JetMedianBackgroundEstimator.cc.
std::vector< PseudoJet > fastjet::JetMedianBackgroundEstimator::jets_used | ( | ) | const |
returns the jets used to actually compute the background properties
Definition at line 476 of file JetMedianBackgroundEstimator.cc.
double fastjet::JetMedianBackgroundEstimator::empty_area | ( | ) | const |
Returns the estimate of the area (within the range defined by the selector) that is not occupied by jets.
The value is that for the last call of rho() or sigma() If the configuration has changed in the meantime, throw an error.
The answer is defined to be zero if the area calculation involved explicit ghosts; if the area calculation was an active area, then use is made of the active area's internal list of pure ghost jets (taking those that pass the selector); otherwise it is based on the difference between the selector's total area and the area of the jets that pass the selector.
The result here is just the cached result of the corresponding call to the ClusterSequenceAreaBase function.
Definition at line 518 of file JetMedianBackgroundEstimator.cc.
double fastjet::JetMedianBackgroundEstimator::n_empty_jets | ( | ) | const |
Returns the number of empty jets used when computing the background properties.
The value is that for the last call of rho() or sigma(). If the configuration has changed in the meantime, throw an error.
If the area has explicit ghosts the result is zero; for active areas it is the number of internal pure ghost jets that pass the selector; otherwise it is deduced from the empty area, divided by (the average pure-ghost-jet area).
The result here is just the cached result of the corresponding call to the ClusterSequenceAreaBase function.
Definition at line 547 of file JetMedianBackgroundEstimator.cc.
void fastjet::JetMedianBackgroundEstimator::reset | ( | ) |
Resets the class to its default state, including the choice to use 4-vector areas.
Definition at line 571 of file JetMedianBackgroundEstimator.cc.
|
inline |
By default when calculating pt/Area for a jet, it is the transverse component of the 4-vector area that is used in the ratiof .
Calling this function with a "false" argument causes the scalar area to be used instead.
While the difference between the two choices is usually small, for high-precision work it is usually the 4-vector area that is to be preferred.
use_it | whether one uses the 4-vector area or not (true by default) |
Definition at line 339 of file JetMedianBackgroundEstimator.hh.
|
inline |
check if the estimator uses the 4-vector area or the scalar area
Definition at line 345 of file JetMedianBackgroundEstimator.hh.
|
inline |
The FastJet v2.X sigma calculation had a small spurious offset in the limit of a small number of jets.
This is fixed by default in versions 3 upwards. The old behaviour can be obtained with a call to this function.
Definition at line 351 of file JetMedianBackgroundEstimator.hh.
void fastjet::JetMedianBackgroundEstimator::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 is null then pt/area is used (as occurs also if this function is not called).
Note that this is still preliminary in FastJet 3.0 and that backward compatibility is not guaranteed in future releases of FastJet
Definition at line 591 of file JetMedianBackgroundEstimator.cc.
|
inline |
return the pointer to the jet density class
Definition at line 366 of file JetMedianBackgroundEstimator.hh.
|
inlineoverridevirtual |
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.
Reimplemented from fastjet::BackgroundEstimatorBase.
Definition at line 378 of file JetMedianBackgroundEstimator.hh.
|
overridevirtual |
returns a textual description of the background estimator
Implements fastjet::BackgroundEstimatorBase.
Definition at line 602 of file JetMedianBackgroundEstimator.cc.