Public Member Functions

fastjet::ClusterSequenceArea Class Reference
[Area-related classes]

General class for user to obtain ClusterSequence with additional area information. More...

#include <ClusterSequenceArea.hh>

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

List of all members.

Public Member Functions

template<class L >
 ClusterSequenceArea (const std::vector< L > &pseudojets, const JetDefinition &jet_def, const AreaDefinition &area_def_in)
 main constructor
template<class L >
 ClusterSequenceArea (const std::vector< L > &pseudojets, const JetDefinition &jet_def, const GhostedAreaSpec &ghost_spec)
 constructor with a GhostedAreaSpec
template<class L >
 ClusterSequenceArea (const std::vector< L > &pseudojets, const JetDefinition &jet_def, const VoronoiAreaSpec &voronoi_spec)
 constructor with a VoronoiAreaSpec
const AreaDefinitionarea_def () const
 return a reference to the area definition
virtual double area (const PseudoJet &jet) const
 return the area associated with the given jet
virtual double area_error (const PseudoJet &jet) const
 return the error (uncertainty) associated with the determination of the area of this jet
virtual PseudoJet area_4vector (const PseudoJet &jet) const
 return the 4-vector area
virtual double empty_area (const RangeDefinition &range) const
 return the total area, in the given rap-phi range, that is free of jets
virtual double n_empty_jets (const RangeDefinition &range) const
 return something similar to the number of pure ghost jets in the given rap-phi range in an active area case.
virtual bool is_pure_ghost (const PseudoJet &jet) const
 true if a jet is made exclusively of ghosts
virtual bool has_explicit_ghosts () const
 true if this ClusterSequence has explicit ghosts
virtual void get_median_rho_and_sigma (const std::vector< PseudoJet > &all_jets, const RangeDefinition &range, bool use_area_4vector, double &median, double &sigma, double &mean_area, bool all_are_incl=false) const
 overload version of what's in the ClusterSequenceAreaBase class, which additionally checks compatibility between "range" and region in which ghosts are thrown.
virtual void get_median_rho_and_sigma (const RangeDefinition &range, bool use_area_4vector, double &median, double &sigma) const
 overload version of what's in the ClusterSequenceAreaBase class, which actually just does the same thing as the base version (but since we've overridden the 5-argument version above, we have to override the 4-argument version too.
virtual void get_median_rho_and_sigma (const RangeDefinition &range, bool use_area_4vector, double &median, double &sigma, double &mean_area) const
 overload version of what's in the ClusterSequenceAreaBase class, which actually just does the same thing as the base version (but since we've overridden the multi-argument version above, we have to override the 5-argument version too.
virtual void parabolic_pt_per_unit_area (double &a, double &b, const RangeDefinition &range, double exclude_above=-1.0, bool use_area_4vector=false) const
 overload version of what's in the ClusterSequenceAreaBase class, which additionally checks compatibility between "range" and region in which ghosts are thrown.

Detailed Description

General class for user to obtain ClusterSequence with additional area information.

Based on the area_def, it automatically dispatches the work to the appropriate actual ClusterSequenceAreaBase-derived-class to do the real work.

Definition at line 51 of file ClusterSequenceArea.hh.


Member Function Documentation

virtual double fastjet::ClusterSequenceArea::n_empty_jets ( const RangeDefinition range  )  const [inline, virtual]

return something similar to the number of pure ghost jets in the given rap-phi range in an active area case.

For the local implementation we return empty_area/(0.55 pi R^2), based on measured properties of ghost jets with kt and cam. Note that the number returned is a double.

Reimplemented from fastjet::ClusterSequenceAreaBase.

Definition at line 115 of file ClusterSequenceArea.hh.

                                                                   {
    return _area_base->n_empty_jets(range);
  }

virtual void fastjet::ClusterSequenceArea::get_median_rho_and_sigma ( const std::vector< PseudoJet > &  all_jets,
const RangeDefinition range,
bool  use_area_4vector,
double &  median,
double &  sigma,
double &  mean_area,
bool  all_are_incl = false 
) const [inline, virtual]

overload version of what's in the ClusterSequenceAreaBase class, which additionally checks compatibility between "range" and region in which ghosts are thrown.

Reimplemented from fastjet::ClusterSequenceAreaBase.

Definition at line 133 of file ClusterSequenceArea.hh.

                                                                         {
    _warn_if_range_unsuitable(range);
    ClusterSequenceAreaBase::get_median_rho_and_sigma(
                                 all_jets, range, use_area_4vector,
                                 median, sigma, mean_area, all_are_incl);
  }

virtual void fastjet::ClusterSequenceArea::get_median_rho_and_sigma ( const RangeDefinition range,
bool  use_area_4vector,
double &  median,
double &  sigma 
) const [inline, virtual]

overload version of what's in the ClusterSequenceAreaBase class, which actually just does the same thing as the base version (but since we've overridden the 5-argument version above, we have to override the 4-argument version too.

Reimplemented from fastjet::ClusterSequenceAreaBase.

Definition at line 149 of file ClusterSequenceArea.hh.

                                                                               {
    ClusterSequenceAreaBase::get_median_rho_and_sigma(range,use_area_4vector,
                                                      median,sigma);
  }

virtual void fastjet::ClusterSequenceArea::get_median_rho_and_sigma ( const RangeDefinition range,
bool  use_area_4vector,
double &  median,
double &  sigma,
double &  mean_area 
) const [inline, virtual]

overload version of what's in the ClusterSequenceAreaBase class, which actually just does the same thing as the base version (but since we've overridden the multi-argument version above, we have to override the 5-argument version too.

Reimplemented from fastjet::ClusterSequenceAreaBase.

Definition at line 160 of file ClusterSequenceArea.hh.

                                                                  {
    ClusterSequenceAreaBase::get_median_rho_and_sigma(range,use_area_4vector,
                                                      median,sigma, mean_area);
  }

virtual void fastjet::ClusterSequenceArea::parabolic_pt_per_unit_area ( double &  a,
double &  b,
const RangeDefinition range,
double  exclude_above = -1.0,
bool  use_area_4vector = false 
) const [inline, virtual]

overload version of what's in the ClusterSequenceAreaBase class, which additionally checks compatibility between "range" and region in which ghosts are thrown.

Reimplemented from fastjet::ClusterSequenceAreaBase.

Definition at line 172 of file ClusterSequenceArea.hh.

                                                                             {
    _warn_if_range_unsuitable(range);
    ClusterSequenceAreaBase::parabolic_pt_per_unit_area(
                                a,b,range, exclude_above, use_area_4vector);
  }


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