fastjet 2.4.5
Classes | Public Member Functions | Private Member Functions | Private Attributes
fastjet::ClusterSequenceVoronoiArea Class Reference

Handle the computation of Voronoi jet area. More...

#include <ClusterSequenceVoronoiArea.hh>

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

List of all members.

Classes

class  VoronoiAreaCalc
 class for carrying out a voronoi area calculation on a set of initial vectors More...

Public Member Functions

template<class L >
 ClusterSequenceVoronoiArea (const std::vector< L > &pseudojets, const JetDefinition &jet_def, const VoronoiAreaSpec &spec=VoronoiAreaSpec(), const bool &writeout_combinations=false)
 template ctor
 ~ClusterSequenceVoronoiArea ()
 default dtor
virtual double area (const PseudoJet &jet) const
 return the area associated with the given jet
virtual PseudoJet area_4vector (const PseudoJet &jet) const
 return a 4-vector area associated with the given jet -- stricly this is not the exact 4-vector area, but rather an approximation made of sums of centres of all Voronoi cells in jet, each contributing with a normalisation equal to the area of the cell
virtual double area_error (const PseudoJet &jet) const
 return the error of the area associated with the given jet (0 by definition for a voronoi area)

Private Member Functions

void _initializeVA ()
 initialisation of the Voronoi Area

Private Attributes

std::vector< double > _voronoi_area
 vector containing the result
std::vector< PseudoJet_voronoi_area_4vector
 vector containing approx 4-vect areas
VoronoiAreaCalc_pa_calc
 area calculator
double _effective_Rfact
 effective radius

Detailed Description

Handle the computation of Voronoi jet area.

Definition at line 46 of file ClusterSequenceVoronoiArea.hh.


Constructor & Destructor Documentation

template<class L >
fastjet::ClusterSequenceVoronoiArea::ClusterSequenceVoronoiArea ( const std::vector< L > &  pseudojets,
const JetDefinition jet_def,
const VoronoiAreaSpec spec = VoronoiAreaSpec(),
const bool &  writeout_combinations = false 
)

template ctor

template constructor need to be specified in the header!

Parameters:
pseudojetlist of jets (template type)
jet_defjet definition
effective_Rfacteffective radius
writeout_combinations??????

Definition at line 99 of file ClusterSequenceVoronoiArea.hh.

                                     :
  _effective_Rfact(spec.effective_Rfact()) {

  // transfer the initial jets (type L) into our own array
  _transfer_input_jets(pseudojets);

  // run the clustering
  _initialise_and_run(jet_def,writeout_combinations);

  // the jet clustering's already been done, now worry about areas...
  _initializeVA();
}
fastjet::ClusterSequenceVoronoiArea::~ClusterSequenceVoronoiArea ( )

default dtor

Definition at line 313 of file ClusterSequenceVoronoiArea.cc.

                                                        {
  delete _pa_calc;
}

Member Function Documentation

void fastjet::ClusterSequenceVoronoiArea::_initializeVA ( ) [private]

initialisation of the Voronoi Area

Definition at line 271 of file ClusterSequenceVoronoiArea.cc.

                                                {
  // run the VAC on our original particles
  _pa_calc = new VAC(_jets.begin(), 
                     _jets.begin()+n_particles(),
                     _effective_Rfact*_jet_def.R());

  // transfer the areas to our local structure
  //  -- first the initial ones
  _voronoi_area.reserve(2*n_particles());
  for (unsigned int i=0; i<n_particles(); i++) {
    _voronoi_area.push_back(_pa_calc->area(i));
    // make a stab at a 4-vector area
    if (_jets[i].perp2() > 0) {
      _voronoi_area_4vector.push_back((_pa_calc->area(i)/_jets[i].perp())
                                      * _jets[i]);
    } else {
      // not sure what to do here -- just put zero (it won't be meaningful
      // anyway)
      _voronoi_area_4vector.push_back(PseudoJet(0.0,0.0,0.0,0.0));
    }
  }
           
  //  -- then the combined areas that arise from the clustering
  for (unsigned int i = n_particles(); i < _history.size(); i++) {
    double area;
    PseudoJet area_4vect;
    if (_history[i].parent2 >= 0) {
      area = _voronoi_area[_history[i].parent1] + 
             _voronoi_area[_history[i].parent2];
      area_4vect = _voronoi_area_4vector[_history[i].parent1] + 
                   _voronoi_area_4vector[_history[i].parent2];
    } else {
      area = _voronoi_area[_history[i].parent1];
      area_4vect = _voronoi_area_4vector[_history[i].parent1];
    }
    _voronoi_area.push_back(area);
    _voronoi_area_4vector.push_back(area_4vect);
  }

}
virtual double fastjet::ClusterSequenceVoronoiArea::area ( const PseudoJet jet) const [inline, virtual]

return the area associated with the given jet

Reimplemented from fastjet::ClusterSequenceAreaBase.

Definition at line 63 of file ClusterSequenceVoronoiArea.hh.

References fastjet::PseudoJet::cluster_hist_index().

Referenced by fastjet::ClusterSequencePassiveArea::_initialise_and_run_PA().

                                                          {
    return _voronoi_area[jet.cluster_hist_index()];};
virtual PseudoJet fastjet::ClusterSequenceVoronoiArea::area_4vector ( const PseudoJet jet) const [inline, virtual]

return a 4-vector area associated with the given jet -- stricly this is not the exact 4-vector area, but rather an approximation made of sums of centres of all Voronoi cells in jet, each contributing with a normalisation equal to the area of the cell

Reimplemented from fastjet::ClusterSequenceAreaBase.

Definition at line 70 of file ClusterSequenceVoronoiArea.hh.

References fastjet::PseudoJet::cluster_hist_index().

Referenced by fastjet::ClusterSequencePassiveArea::_initialise_and_run_PA().

                                                                     {
    return _voronoi_area_4vector[jet.cluster_hist_index()];};
virtual double fastjet::ClusterSequenceVoronoiArea::area_error ( const PseudoJet jet) const [inline, virtual]

return the error of the area associated with the given jet (0 by definition for a voronoi area)

Reimplemented from fastjet::ClusterSequenceAreaBase.

Definition at line 75 of file ClusterSequenceVoronoiArea.hh.

                                                                {
    return 0.0;};

Member Data Documentation

effective radius

Definition at line 90 of file ClusterSequenceVoronoiArea.hh.

area calculator

Definition at line 89 of file ClusterSequenceVoronoiArea.hh.

vector containing the result

Definition at line 87 of file ClusterSequenceVoronoiArea.hh.

vector containing approx 4-vect areas

Definition at line 88 of file ClusterSequenceVoronoiArea.hh.


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