fastjet 2.4.5
Public Member Functions | Private Member Functions | Private Attributes | Static Private Attributes
fastjet::ClusterSequenceArea Class Reference

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.

Private Member Functions

void _warn_if_range_unsuitable (const RangeDefinition &range) const
 print a warning if the range is unsuitable for the current calculation of the area (e.g.
template<class L >
void initialize_and_run_cswa (const std::vector< L > &pseudojets, const JetDefinition &jet_def)

Private Attributes

std::auto_ptr
< ClusterSequenceAreaBase
_area_base
AreaDefinition _area_def

Static Private Attributes

static LimitedWarning _range_warnings
static LimitedWarning _explicit_ghosts_repeats_warnings

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 49 of file ClusterSequenceArea.hh.


Constructor & Destructor Documentation

template<class L >
fastjet::ClusterSequenceArea::ClusterSequenceArea ( const std::vector< L > &  pseudojets,
const JetDefinition jet_def,
const AreaDefinition area_def_in 
) [inline]

main constructor

Definition at line 53 of file ClusterSequenceArea.hh.

                                               : _area_def(area_def_in) {
    initialize_and_run_cswa(pseudojets, jet_def);
  }
template<class L >
fastjet::ClusterSequenceArea::ClusterSequenceArea ( const std::vector< L > &  pseudojets,
const JetDefinition jet_def,
const GhostedAreaSpec ghost_spec 
) [inline]

constructor with a GhostedAreaSpec

Definition at line 61 of file ClusterSequenceArea.hh.

                                                : _area_def(ghost_spec){
    initialize_and_run_cswa(pseudojets, jet_def);
  }
template<class L >
fastjet::ClusterSequenceArea::ClusterSequenceArea ( const std::vector< L > &  pseudojets,
const JetDefinition jet_def,
const VoronoiAreaSpec voronoi_spec 
) [inline]

constructor with a VoronoiAreaSpec

Definition at line 69 of file ClusterSequenceArea.hh.

                                                  : _area_def(voronoi_spec){
    initialize_and_run_cswa(pseudojets, jet_def);
  }

Member Function Documentation

void fastjet::ClusterSequenceArea::_warn_if_range_unsuitable ( const RangeDefinition range) const [private]

print a warning if the range is unsuitable for the current calculation of the area (e.g.

because ghosts do not extend far enough).

Definition at line 11 of file ClusterSequenceArea.cc.

References fastjet::RangeDefinition::get_rap_limits(), fastjet::kt_algorithm, fastjet::passive_area, and fastjet::voronoi_area.

                                                                                       {
  bool no_ghosts = (_area_def.area_type() == voronoi_area)
    || (_area_def.area_type() == passive_area
        && jet_def().jet_algorithm() == kt_algorithm);
  if (! no_ghosts) {
    double rapmin, rapmax;
    range.get_rap_limits(rapmin, rapmax);
    if (rapmin < -_area_def.ghost_spec().ghost_maxrap()+0.95*jet_def().R() ||
        rapmax >  _area_def.ghost_spec().ghost_maxrap()-0.95*jet_def().R()) {
      _range_warnings.warn("rapidity range for median (rho) extends beyond +-(ghost_maxrap - 0.95*R); this is likely to cause the results to be unreliable; safest option is to increase ghost_maxrap in the area definition");
    }
  }
}
virtual double fastjet::ClusterSequenceArea::area ( const PseudoJet jet) const [inline, virtual]

return the area associated with the given jet

Reimplemented from fastjet::ClusterSequenceAreaBase.

Definition at line 80 of file ClusterSequenceArea.hh.

                                                          {
    return _area_base->area(jet);}
virtual PseudoJet fastjet::ClusterSequenceArea::area_4vector ( const PseudoJet jet) const [inline, virtual]

return the 4-vector area

Reimplemented from fastjet::ClusterSequenceAreaBase.

Definition at line 89 of file ClusterSequenceArea.hh.

                                                              {
    return _area_base->area_4vector(jet);}
const AreaDefinition& fastjet::ClusterSequenceArea::area_def ( ) const [inline]

return a reference to the area definition

Definition at line 76 of file ClusterSequenceArea.hh.

{return _area_def;}
virtual double fastjet::ClusterSequenceArea::area_error ( const PseudoJet jet) const [inline, virtual]

return the error (uncertainty) associated with the determination of the area of this jet

Reimplemented from fastjet::ClusterSequenceAreaBase.

Definition at line 85 of file ClusterSequenceArea.hh.

                                                          {
    return _area_base->area_error(jet);}
virtual double fastjet::ClusterSequenceArea::empty_area ( const RangeDefinition range) const [inline, virtual]

return the total area, in the given rap-phi range, that is free of jets

Reimplemented from fastjet::ClusterSequenceAreaBase.

Definition at line 105 of file ClusterSequenceArea.hh.

                                                                 {
    return _area_base->empty_area(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 131 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 147 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 158 of file ClusterSequenceArea.hh.

                                                                  {
    ClusterSequenceAreaBase::get_median_rho_and_sigma(range,use_area_4vector,
                                                      median,sigma, mean_area);
  }
virtual bool fastjet::ClusterSequenceArea::has_explicit_ghosts ( ) const [inline, virtual]

true if this ClusterSequence has explicit ghosts

Reimplemented from fastjet::ClusterSequenceAreaBase.

Definition at line 123 of file ClusterSequenceArea.hh.

                                           {
    return _area_base->has_explicit_ghosts();
  }
template<class L >
void fastjet::ClusterSequenceArea::initialize_and_run_cswa ( const std::vector< L > &  pseudojets,
const JetDefinition jet_def 
) [private]

Definition at line 199 of file ClusterSequenceArea.hh.

References fastjet::active_area, fastjet::active_area_explicit_ghosts, fastjet::one_ghost_passive_area, fastjet::passive_area, and fastjet::voronoi_area.

 {
  
  ClusterSequenceAreaBase * _area_base_ptr;
  switch(_area_def.area_type()) {
  case active_area:
    _area_base_ptr = new ClusterSequenceActiveArea(pseudojets, 
                                                   jet_def, 
                                                   _area_def.ghost_spec());
    break;
  case active_area_explicit_ghosts:
    if (_area_def.ghost_spec().repeat() != 1) 
      _explicit_ghosts_repeats_warnings.warn("Requested active area with explicit ghosts with repeat != 1; only 1 set of ghosts will be used");
    _area_base_ptr = new ClusterSequenceActiveAreaExplicitGhosts(pseudojets, 
                                                   jet_def, 
                                                   _area_def.ghost_spec());
    break;
  case voronoi_area:
    _area_base_ptr = new ClusterSequenceVoronoiArea(pseudojets, 
                                                   jet_def, 
                                                   _area_def.voronoi_spec());
    break;
  case one_ghost_passive_area:
    _area_base_ptr = new ClusterSequence1GhostPassiveArea(pseudojets, 
                                                    jet_def, 
                                                    _area_def.ghost_spec());
    break;
  case passive_area:
    _area_base_ptr = new ClusterSequencePassiveArea(pseudojets, 
                                                    jet_def, 
                                                    _area_def.ghost_spec());
    break;
  default:
    std::cerr << "Error: unrecognized area_type in ClusterSequenceArea:" 
              << _area_def.area_type() << std::endl;
    exit(-1);
  }
  // now copy across the information from the area base class
  _area_base = std::auto_ptr<ClusterSequenceAreaBase>(_area_base_ptr);
  transfer_from_sequence(*_area_base);
}
virtual bool fastjet::ClusterSequenceArea::is_pure_ghost ( const PseudoJet jet) const [inline, virtual]

true if a jet is made exclusively of ghosts

Reimplemented from fastjet::ClusterSequenceAreaBase.

Definition at line 118 of file ClusterSequenceArea.hh.

                                                          {
    return _area_base->is_pure_ghost(jet);
  }
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 113 of file ClusterSequenceArea.hh.

                                                                   {
    return _area_base->n_empty_jets(range);
  }
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 170 of file ClusterSequenceArea.hh.

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

Member Data Documentation

Definition at line 191 of file ClusterSequenceArea.hh.

Definition at line 192 of file ClusterSequenceArea.hh.

Definition at line 194 of file ClusterSequenceArea.hh.

Definition at line 193 of file ClusterSequenceArea.hh.


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