fastjet 2.4.5
|
JadePlugin is a plugin for fastjet (v2.4 upwards) More...
#include <JadePlugin.hh>
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: | |
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 |
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.
fastjet::JadePlugin::JadePlugin | ( | ) | [inline] |
fastjet::JadePlugin::JadePlugin | ( | const JadePlugin & | plugin | ) | [inline] |
string fastjet::JadePlugin::description | ( | ) | const [virtual] |
return a textual description of the jet-definition implemented in this plugin
Implements fastjet::JetDefinition::Plugin.
Definition at line 81 of file JadePlugin.cc.
{ ostringstream desc; desc << "e+e- JADE algorithm plugin"; return desc.str(); }
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 fastjet::JetDefinition::Plugin.
Definition at line 96 of file JadePlugin.hh.
{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 fastjet::JetDefinition::Plugin.
Definition at line 92 of file JadePlugin.hh.
{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:
Implements fastjet::JetDefinition::Plugin.
Definition at line 88 of file JadePlugin.cc.
References fastjet::NNH< BJ, I >::dij_min(), fastjet::ClusterSequence::jets(), fastjet::ClusterSequence::plugin_record_iB_recombination(), and fastjet::ClusterSequence::plugin_record_ij_recombination().
{ int njets = cs.jets().size(); NNH<JadeBriefJet> nnh(cs.jets()); // if testing against Hoeth's implementation, need to rescale the // dij by Q^2. //double Q2 = cs.Q2(); while (njets > 0) { int i, j, k; double dij = nnh.dij_min(i, j); if (j >= 0) { cs.plugin_record_ij_recombination(i, j, dij, k); nnh.merge_jets(i, j, cs.jets()[k], k); } else { double diB = cs.jets()[i].E()*cs.jets()[i].E(); // get new diB cs.plugin_record_iB_recombination(i, diB); nnh.remove_jet(i); } njets--; } }