FastJet 3.0.2
|
Like ClusterSequence with computation of the Voronoi jet area. More...
#include <fastjet/ClusterSequenceVoronoiArea.hh>
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 &) const |
return the error of the area associated with the given jet (0 by definition for a voronoi area) |
Like ClusterSequence with computation of the Voronoi jet area.
Handle the computation of Voronoi jet area.
This class should not be used directly. Rather use ClusterSequenceArea with the appropriate AreaDefinition
Definition at line 48 of file ClusterSequenceVoronoiArea.hh.
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!
pseudojet | list of jets (template type) |
jet_def | jet definition |
effective_Rfact | effective radius |
writeout_combinations | ?????? |
Definition at line 101 of file ClusterSequenceVoronoiArea.hh.