fastjet::JadePlugin Class Reference

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

#include <JadePlugin.hh>

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

List of all members.

Public Member Functions

 JadePlugin ()
 Main constructor for the Jade Plugin class.
 JadePlugin (const JadePlugin &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(.

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

Detailed Description

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

It implements the JADE algorithm, which is an e+e- sequential recombination algorithm with interparticle distance

dij = 2 E_i E_j (1 - cos theta_ij)

or equivalently

yij = dij/E_{vis}^2

This corresponds to the distance measured used in

"Experimental Investigation of the Energy Dependence of the Strong Coupling Strength." JADE Collaboration (S. Bethke et al.) Phys.Lett.B213:235,1988

The JADE article carries out particle recombinations in the E-scheme (4-vector recombination), which is the default procedure for this plugin.

NOTE: other widely used schemes include E0, P, P0; however they also involve modifications to the distance measure. Be sure of what you're doing before running a JADE type algorithm.

To access the jets with a given ycut value (clustering stops once all yij > ycut), use

vector<PseudoJet> jets = cluster_sequence.exclusive_jets_ycut(ycut);

and related routines.

Definition at line 75 of file JadePlugin.hh.


Constructor & Destructor Documentation

fastjet::JadePlugin::JadePlugin (  )  [inline]

Main constructor for the Jade Plugin class.

Definition at line 78 of file JadePlugin.hh.

00078 {}

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

copy constructor

Definition at line 81 of file JadePlugin.hh.

00081                                          {
00082     *this = plugin;
00083   }


Member Function Documentation

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

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

Implements JetDefinition::Plugin.

Definition at line 81 of file JadePlugin.cc.

00081                                       {
00082   ostringstream desc;
00083   desc << "e+e- JADE algorithm plugin";
00084   return desc.str();
00085 }

virtual bool fastjet::JadePlugin::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 96 of file JadePlugin.hh.

00096 {return true;}

virtual double fastjet::JadePlugin::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 92 of file JadePlugin.hh.

00092 {return 1.0;}

void fastjet::JadePlugin::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 88 of file JadePlugin.cc.

References ClusterSequence::jets(), ClusterSequence::plugin_record_iB_recombination(), and ClusterSequence::plugin_record_ij_recombination().

00088                                                           {
00089   int njets = cs.jets().size();
00090   NNH<JadeBriefJet> nnh(cs.jets());
00091 
00092   // if testing against Hoeth's implementation, need to rescale the
00093   // dij by Q^2.
00094   //double Q2 = cs.Q2(); 
00095 
00096   while (njets > 0) {
00097     int i, j, k;
00098     double dij = nnh.dij_min(i, j);
00099 
00100     if (j >= 0) {
00101       cs.plugin_record_ij_recombination(i, j, dij, k);
00102       nnh.merge_jets(i, j, cs.jets()[k], k);
00103     } else {
00104       double diB = cs.jets()[i].E()*cs.jets()[i].E(); // get new diB
00105       cs.plugin_record_iB_recombination(i, diB);
00106       nnh.remove_jet(i);
00107     }
00108     njets--;
00109   }
00110 }


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

Generated on 26 Feb 2010 for fastjet by  doxygen 1.6.1