FastJet
3.1.3
|
Background Estimator based on the median pt/area of a set of grid cells. More...
#include <fastjet/tools/GridMedianBackgroundEstimator.hh>
Public Member Functions | |
constructors and destructors | |
GridMedianBackgroundEstimator (double ymax, double requested_grid_spacing) | |
GridMedianBackgroundEstimator (const RectangularGrid &grid) | |
Constructor based on a user's fully specified RectangularGrid. | |
GridMedianBackgroundEstimator (double rapmin_in, double rapmax_in, double drap_in, double dphi_in, Selector tile_selector=Selector()) | |
Constructor with the explicit parameters for the underlying RectangularGrid. More... | |
setting a new event | |
void | set_particles (const std::vector< PseudoJet > &particles) |
tell the background estimator that it has a new event, composed of the specified particles. More... | |
void | set_compute_rho_m (bool enable) |
determine whether the automatic calculation of rho_m and sigma_m is enabled (by default true) | |
retrieving fundamental information | |
double | rho () const |
returns rho, the median background density per unit area | |
double | sigma () const |
returns sigma, the background fluctuations per unit area; must be multipled by sqrt(area) to get fluctuations for a region of a given area. More... | |
double | rho (const PseudoJet &jet) |
returns rho, the background density per unit area, locally at the position of a given jet. More... | |
double | sigma (const PseudoJet &jet) |
returns sigma, the background fluctuations per unit area, locally at the position of a given jet. More... | |
bool | has_sigma () |
returns true if this background estimator has support for determination of sigma | |
double | rho_m () const |
Returns rho_m, the purely longitudinal, particle-mass-induced component of the background density per unit area. | |
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 area; must be multipled by sqrt(area) to get fluctuations for a region of a given area. More... | |
double | rho_m (const PseudoJet &jet) |
Returns rho_m locally at the jet position. As for rho(jet), it is non-const. | |
double | sigma_m (const PseudoJet &jet) |
Returns sigma_m locally at the jet position. As for rho(jet), it is non-const. | |
bool | has_rho_m () const |
Returns true if this background estimator has support for determination of rho_m. More... | |
double | mean_area () const |
returns the area of the grid cells (all identical, but referred to as "mean" area for uniformity with JetMedianBGE). More... | |
configuring the behaviour | |
virtual void | set_rescaling_class (const FunctionOfPseudoJet< double > *rescaling_class) |
Set a pointer to a class that calculates the rescaling factor as a function of the jet (position). More... | |
description | |
std::string | description () const |
returns a textual description of the background estimator | |
Public Member Functions inherited from fastjet::BackgroundEstimatorBase | |
virtual | ~BackgroundEstimatorBase () |
a default virtual destructor that does nothing | |
BackgroundEstimatorBase () | |
const FunctionOfPseudoJet < double > * | rescaling_class () const |
return the pointer to the jet density class | |
Public Member Functions inherited from fastjet::RectangularGrid | |
RectangularGrid (double rapmax_in, double cell_size) | |
ctor with simple initialisation More... | |
RectangularGrid (double rapmin_in, double rapmax_in, double drap_in, double dphi_in, Selector tile_selector=Selector()) | |
ctor with more control over initialisation More... | |
RectangularGrid () | |
dummy ctor (will give an unusable grid) | |
virtual int | n_tiles () const |
returns the total number of tiles in the tiling; valid tile indices run from 0 ... More... | |
virtual int | n_good_tiles () const |
returns the number of tiles that are "good"; i.e. More... | |
virtual int | tile_index (const PseudoJet &p) const |
returns the index of the tile in which p is located, or -1 if p is outside the tiling region | |
virtual bool | tile_is_good (int itile) const |
returns whether a given tile is good | |
virtual double | tile_area (int) const |
returns the area of tile itile. | |
virtual double | mean_tile_area () const |
returns the mean area of tiles. | |
double | rapmin () const |
returns the minimum rapidity extent of the grid | |
double | rapmax () const |
returns the maxmium rapidity extent of the grid | |
double | drap () const |
returns the spacing of the grid in rapidity | |
double | dphi () const |
returns the spacing of the grid in azimuth | |
virtual bool | is_initialised () const |
returns true if the grid is in a suitably initialised state | |
Public Member Functions inherited from fastjet::TilingBase | |
virtual bool | all_tiles_good () const |
returns whether all tiles are good | |
virtual bool | all_tiles_equal_area () const |
returns true if all tiles have the same area | |
bool | is_initialized () const |
virtual | ~TilingBase () |
virtual destructor | |
Additional Inherited Members | |
Protected Member Functions inherited from fastjet::BackgroundEstimatorBase | |
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 |
Static Protected Attributes inherited from fastjet::BackgroundEstimatorBase | |
static LimitedWarning | _warnings_empty_area |
Background Estimator based on the median pt/area of a set of grid cells.
Description of the method: This background estimator works by projecting the event onto a grid in rapidity and azimuth. In each grid cell, the scalar pt sum of the particles in the cell is computed. The background density is then estimated by the median of (scalar pt sum/cell area) for all cells.
Parameters: The class takes 2 arguments: the absolute rapidity extent of the cells and the size of the grid cells. Note that the size of the cell will be adjusted in azimuth to satisfy the 2pi periodicity and in rapidity to match the requested rapidity extent.
Rescaling: It is possible to use a rescaling profile. In this case, the profile needs to be set before setting the particles and it will be applied to each particle (i.e. not to each cell). Note also that in this case one needs to call rho(jet) instead of rho() [Without rescaling, they are identical]
Definition at line 77 of file GridMedianBackgroundEstimator.hh.
|
inline |
ymax | maximal absolute rapidity extent of the grid |
requested_grid_spacing | size of the grid cell. The "real" cell size could differ due e.g. to the 2pi periodicity in azimuthal angle (size, not area) |
Definition at line 92 of file GridMedianBackgroundEstimator.hh.
|
inline |
Constructor with the explicit parameters for the underlying RectangularGrid.
rapmin | the minimum rapidity extent of the grid |
rapmax | the maximum rapidity extent of the grid |
drap | the grid spacing in rapidity |
dphi | the grid spacing in azimuth |
tile_selector | optional (geometric) selector to specify which tiles are good; a tile is good if a massless 4-vector at the center of the tile passes the selection |
Definition at line 117 of file GridMedianBackgroundEstimator.hh.
|
virtual |
tell the background estimator that it has a new event, composed of the specified particles.
Implements fastjet::BackgroundEstimatorBase.
Definition at line 43 of file GridMedianBackgroundEstimator.cc.
|
virtual |
returns sigma, the background fluctuations 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 146 of file GridMedianBackgroundEstimator.cc.
|
virtual |
returns 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.
Implements fastjet::BackgroundEstimatorBase.
Definition at line 157 of file GridMedianBackgroundEstimator.cc.
|
virtual |
returns sigma, the background fluctuations per unit area, locally at the position of a given jet.
As for rho(jet), it is non-const.
Reimplemented from fastjet::BackgroundEstimatorBase.
Definition at line 167 of file GridMedianBackgroundEstimator.cc.
|
virtual |
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 188 of file GridMedianBackgroundEstimator.cc.
|
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 from fastjet::BackgroundEstimatorBase.
Definition at line 198 of file GridMedianBackgroundEstimator.hh.
|
inline |
returns the area of the grid cells (all identical, but referred to as "mean" area for uniformity with JetMedianBGE).
Definition at line 203 of file GridMedianBackgroundEstimator.hh.
|
virtual |
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.
Note that this has to be called BEFORE any attempt to do an actual computation
The same profile will be used for both pt and mt (this is probabaly a good approximation since the particle density changes is what dominates the rapidity profile)
Reimplemented from fastjet::BackgroundEstimatorBase.
Definition at line 252 of file GridMedianBackgroundEstimator.cc.