FastJet 3.0.2
fastjet::GridMedianBackgroundEstimator Class Reference

Background Estimator based on the median pt/area of a set of grid cells. More...

#include <fastjet/tools/GridMedianBackgroundEstimator.hh>

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

List of all members.

Public Member Functions

constructors and destructors
 GridMedianBackgroundEstimator (double ymax, double requested_grid_spacing)
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.
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.
double rho (const PseudoJet &jet)
 returns rho, the background density per unit area, locally at the position of a given jet.
double sigma (const PseudoJet &jet)
 returns sigma, the background fluctuations per unit area, locally at the position of a given jet.
bool has_sigma ()
 returns true if this background estimator has support for determination of sigma
double mean_area () const
 returns the area of the grid cells (all identical, but referred to as "mean" area for uniformity with JetMedianBGE).
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).
description
std::string description () const
 returns a textual description of the background estimator

Detailed Description

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 63 of file GridMedianBackgroundEstimator.hh.


Constructor & Destructor Documentation

fastjet::GridMedianBackgroundEstimator::GridMedianBackgroundEstimator ( double  ymax,
double  requested_grid_spacing 
) [inline]
Parameters:
ymaxmaximal absolute rapidity extent of the grid
requested_grid_spacingsize 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 72 of file GridMedianBackgroundEstimator.hh.


Member Function Documentation

void fastjet::GridMedianBackgroundEstimator::set_particles ( const std::vector< PseudoJet > &  particles) [virtual]

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

Implements fastjet::BackgroundEstimatorBase.

double fastjet::GridMedianBackgroundEstimator::sigma ( ) const [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 71 of file GridMedianBackgroundEstimator.cc.

double fastjet::GridMedianBackgroundEstimator::rho ( const PseudoJet jet) [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 86 of file GridMedianBackgroundEstimator.cc.

double fastjet::GridMedianBackgroundEstimator::sigma ( const PseudoJet jet) [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 96 of file GridMedianBackgroundEstimator.cc.

double fastjet::GridMedianBackgroundEstimator::mean_area ( ) const [inline]

returns the area of the grid cells (all identical, but referred to as "mean" area for uniformity with JetMedianBGE).

Definition at line 118 of file GridMedianBackgroundEstimator.hh.

void fastjet::GridMedianBackgroundEstimator::set_rescaling_class ( const FunctionOfPseudoJet< double > *  rescaling_class) [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

Reimplemented from fastjet::BackgroundEstimatorBase.

Definition at line 134 of file GridMedianBackgroundEstimator.cc.


The documentation for this class was generated from the following files:
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends