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:

  • plugin_do_ij_recombination(.

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.

00062 : _ycut(ycut) {}

fastjet::EECambridgePlugin::EECambridgePlugin ( const EECambridgePlugin plugin  )  [inline]

copy constructor

Definition at line 65 of file EECambridgePlugin.hh.

00065                                                        {
00066     *this = plugin;
00067   }


Member Function Documentation

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

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

Implements JetDefinition::Plugin.

Definition at line 71 of file EECambridgePlugin.cc.

References ycut().

00071                                              {
00072   ostringstream desc;
00073   desc << "EECambridge plugin with ycut = " << ycut() ;
00074   return desc.str();
00075 }

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 JetDefinition::Plugin.

Definition at line 82 of file EECambridgePlugin.hh.

00082 {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 JetDefinition::Plugin.

Definition at line 78 of file EECambridgePlugin.hh.

00078 {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 JetDefinition::Plugin.

Definition at line 77 of file EECambridgePlugin.cc.

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

00077                                                                  {
00078   int njets = cs.jets().size();
00079   NNH<EECamBriefJet> nnh(cs.jets());
00080 
00081   double Q2 = cs.Q2(); 
00082 
00083   while (njets > 0) {
00084     int i, j, k;
00085     // here we get a minimum based on the purely angular variable from the NNH class
00086     // (called dij there, but vij in the Cambridge article (which uses dij for 
00087     // a kt distance...)
00088     double vij = nnh.dij_min(i, j); // i,j are return values...
00089 
00090     // next we work out the dij (ee kt distance), and based on its
00091     // value decide whether we have soft-freezing (represented here by
00092     // a "Beam" clustering) or not
00093     double dij;
00094     if (j >= 0) {
00095       double scale = min(cs.jets()[i].E(), cs.jets()[j].E());
00096       dij = 2 * vij * scale * scale;
00097       if (dij > Q2 * ycut()) {
00098         // we'll call the softer partner a "beam" jet
00099         if (cs.jets()[i].E() > cs.jets()[j].E()) swap(i,j);
00100         j = -1;
00101       }
00102     } else {
00103       // for the last particle left, just use yij = 1
00104       dij = Q2;
00105     }
00106     
00107     if (j >= 0) {
00108       cs.plugin_record_ij_recombination(i, j, dij, k);
00109       nnh.merge_jets(i, j, cs.jets()[k], k);
00110     } else {
00111       cs.plugin_record_iB_recombination(i, dij);
00112       nnh.remove_jet(i);
00113     }
00114     njets--;
00115   }
00116     
00117 }

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

Definition at line 73 of file EECambridgePlugin.hh.

Referenced by description(), and run_clustering().

00073 {return _ycut;}


Member Data Documentation

Definition at line 85 of file EECambridgePlugin.hh.


The documentation for this class was generated from the following files:

Generated on 26 Feb 2010 for fastjet by  doxygen 1.6.1