FastJet 3.0.6
Public Member Functions | Protected Member Functions | Protected Attributes | Static Protected Attributes
fastjet::D0RunIBaseConePlugin Class Reference

D0RunIBaseConePlugin is base class for a plugin for FastJet (v3.0 or later) that provides an interface to the D0 version of Run-I cone algorithm. More...

#include <fastjet/D0RunIBaseConePlugin.hh>

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

List of all members.

Public Member Functions

 D0RunIBaseConePlugin (double CONErad_in, double JETmne_in, double SPLifr_in=_DEFAULT_SPLifr)
 A D0RunIConePlugin constructor which sets the "free" parameters of the algorithm:
double CONErad () const
double JETmne () const
double SPLifr () const
double TWOrad () const
bool D0_Angle () const
bool Increase_Delta_R () const
bool Kill_Far_Clusters () const
bool Jet_Et_Min_On_Iter () const
double Far_Ratio () const
double Eitem_Negdrop () const
double Et_Min_Ratio () const
double Thresh_Diff_Et () const
double overlap_threshold () const
 access the split_ratio() also by the name overlap_threshold()
virtual std::string description () const =0
 return a textual description of the jet-definition implemented in this plugin
virtual void run_clustering (ClusterSequence &) const =0
 given a ClusterSequence that has been filled up with initial particles, the following function should fill up the rest of the ClusterSequence, using the following member functions of ClusterSequence:
virtual double R () const
 the plugin mechanism's standard way of accessing the jet radius

Protected Member Functions

template<typename HepEntityType >
void run_clustering_worker (ClusterSequence &) const

Protected Attributes

double _CONErad
double _JETmne
double _SPLifr
double _TWOrad
bool _D0_Angle
bool _Increase_Delta_R
bool _Kill_Far_Clusters
bool _Jet_Et_Min_On_Iter
double _Far_Ratio
double _Eitem_Negdrop
double _Et_Min_Ratio
double _Thresh_Diff_Et

Static Protected Attributes

static const double _DEFAULT_SPLifr = 0.5
static const double _DEFAULT_TWOrad = 0.
static const bool _DEFAULT_D0_Angle = false
static const bool _DEFAULT_Increase_Delta_R = true
static const bool _DEFAULT_Kill_Far_Clusters = true
static const bool _DEFAULT_Jet_Et_Min_On_Iter = true
static const double _DEFAULT_Far_Ratio = 0.5
static const double _DEFAULT_Eitem_Negdrop = -1.0
static const double _DEFAULT_Et_Min_Ratio = 0.5
static const double _DEFAULT_Thresh_Diff_Et = 0.01

Detailed Description

D0RunIBaseConePlugin is base class for a plugin for FastJet (v3.0 or later) that provides an interface to the D0 version of Run-I cone algorithm.

Note that this base class is purely virtual and thus needs to be overloaded. In practice this means that you should use one of D0RunIConePlugin or D0RunIpre96ConePlugin.

The D0 code has been obtained from Lars Sonnenschein's web-space http://www-d0.fnal.gov/~sonne/D0RunIcone.tgz

The version of the D0 Run I code distributed here has been modified by the FastJet authors, so as to provide access to the contents of the jets (as is necessary for the plugin). This does not modify the results of the clustering.

Definition at line 59 of file D0RunIBaseConePlugin.hh.


Constructor & Destructor Documentation

fastjet::D0RunIBaseConePlugin::D0RunIBaseConePlugin ( double  CONErad_in,
double  JETmne_in,
double  SPLifr_in = _DEFAULT_SPLifr 
) [inline]

A D0RunIConePlugin constructor which sets the "free" parameters of the algorithm:

Parameters:
CONEradis the cone radius
JETmneis a minimum ET requirement on every iteration (jet dropped if Et < JETmne * Et_min_ratio ). The value that has been used by D0 for JETmne: 8 GeV (and Et_min_ratio is 0.5)
SPlifris the shared Et fraction splitting threshold, and a value of 0.5 was usually used by D0

The remaining parameters of the algorithm are not to be modified if the algorithm is to correspond to the one actually used by D0.

Definition at line 76 of file D0RunIBaseConePlugin.hh.


Member Function Documentation

virtual void fastjet::D0RunIBaseConePlugin::run_clustering ( ClusterSequence ) const [pure virtual]

given a ClusterSequence that has been filled up with initial particles, the following function should fill up the rest of the ClusterSequence, using the following member functions of ClusterSequence:

  • plugin_do_ij_recombination(...)
  • plugin_do_iB_recombination(...)

Implements fastjet::JetDefinition::Plugin.

Implemented in fastjet::D0RunIConePlugin, and fastjet::D0RunIpre96ConePlugin.


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