fastjet 2.4.5
Public Member Functions | Private Attributes
fastjet::EECambridgePlugin Class Reference

EECambridgePlugin is a plugin for fastjet (v2.4 upwards) More...

#include <EECambridgePlugin.hh>

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

List of all members.

Public Member Functions

 EECambridgePlugin (double ycut)
 Main constructor for the EECambridge Plugin class.
 EECambridgePlugin (const EECambridgePlugin &plugin)
 copy constructor
virtual std::string description () const
 return a textual description of the jet-definition implemented in this plugin
virtual void run_clustering (ClusterSequence &) const
 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:
double ycut () const
virtual double R () const
 the plugin mechanism's standard way of accessing the jet radius.
virtual bool exclusive_sequence_meaningful () const
 avoid the warning whenever the user requests "exclusive" jets from the cluster sequence

Private Attributes

double _ycut

Detailed Description

EECambridgePlugin is a plugin for fastjet (v2.4 upwards)

It implements the Cambridge algorithm, as defined in

Better jet clustering algorithms Yuri Dokshitzer, Garth Leder, Stefano Moretti, Bryan Webber JHEP 9708 (1997) 001 http://www-spires.slac.stanford.edu/spires/find/hep/www?rawcmd=FIND+j+JHEPA%2C9708%2C001

On construction one must supply a ycut value.

To get the jets at the end call ClusterSequence::inclusive_jets();

Definition at line 57 of file EECambridgePlugin.hh.


Constructor & Destructor Documentation

fastjet::EECambridgePlugin::EECambridgePlugin ( double  ycut) [inline]

Main constructor for the EECambridge Plugin class.

It takes the dimensionless parameter ycut (the Q value for normalisation of the kt-distances is taken from the sum of all particle energies).

Definition at line 62 of file EECambridgePlugin.hh.

: _ycut(ycut) {}
fastjet::EECambridgePlugin::EECambridgePlugin ( const EECambridgePlugin plugin) [inline]

copy constructor

Definition at line 65 of file EECambridgePlugin.hh.

                                                       {
    *this = plugin;
  }

Member Function Documentation

string fastjet::EECambridgePlugin::description ( ) const [virtual]

return a textual description of the jet-definition implemented in this plugin

Implements fastjet::JetDefinition::Plugin.

Definition at line 71 of file EECambridgePlugin.cc.

                                             {
  ostringstream desc;
  desc << "EECambridge plugin with ycut = " << ycut() ;
  return desc.str();
}
virtual bool fastjet::EECambridgePlugin::exclusive_sequence_meaningful ( ) const [inline, virtual]

avoid the warning whenever the user requests "exclusive" jets from the cluster sequence

Reimplemented from fastjet::JetDefinition::Plugin.

Definition at line 82 of file EECambridgePlugin.hh.

{return true;}
virtual double fastjet::EECambridgePlugin::R ( ) const [inline, virtual]

the plugin mechanism's standard way of accessing the jet radius.

This must be set to return something sensible, even if R does not make sense for this algorithm!

Implements fastjet::JetDefinition::Plugin.

Definition at line 78 of file EECambridgePlugin.hh.

{return 1.0;}
void fastjet::EECambridgePlugin::run_clustering ( ClusterSequence ) const [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.

Definition at line 77 of file EECambridgePlugin.cc.

References fastjet::ClusterSequence::jets(), fastjet::d0::inline_maths::min(), fastjet::ClusterSequence::plugin_record_iB_recombination(), fastjet::ClusterSequence::plugin_record_ij_recombination(), and fastjet::ClusterSequence::Q2().

                                                                 {
  int njets = cs.jets().size();
  NNH<EECamBriefJet> nnh(cs.jets());

  double Q2 = cs.Q2(); 

  while (njets > 0) {
    int i, j, k;
    // here we get a minimum based on the purely angular variable from the NNH class
    // (called dij there, but vij in the Cambridge article (which uses dij for 
    // a kt distance...)
    double vij = nnh.dij_min(i, j); // i,j are return values...

    // next we work out the dij (ee kt distance), and based on its
    // value decide whether we have soft-freezing (represented here by
    // a "Beam" clustering) or not
    double dij;
    if (j >= 0) {
      double scale = min(cs.jets()[i].E(), cs.jets()[j].E());
      dij = 2 * vij * scale * scale;
      if (dij > Q2 * ycut()) {
        // we'll call the softer partner a "beam" jet
        if (cs.jets()[i].E() > cs.jets()[j].E()) swap(i,j);
        j = -1;
      }
    } else {
      // for the last particle left, just use yij = 1
      dij = Q2;
    }
    
    if (j >= 0) {
      cs.plugin_record_ij_recombination(i, j, dij, k);
      nnh.merge_jets(i, j, cs.jets()[k], k);
    } else {
      cs.plugin_record_iB_recombination(i, dij);
      nnh.remove_jet(i);
    }
    njets--;
  }
    
}
double fastjet::EECambridgePlugin::ycut ( ) const [inline]

Definition at line 73 of file EECambridgePlugin.hh.

{return _ycut;}

Member Data Documentation

Definition at line 85 of file EECambridgePlugin.hh.


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