Main Page | Namespace List | Class Hierarchy | Class List | Directories | File List | Namespace Members | Class Members | File Members

fastjet::SISConePlugin Class Reference

SISConePlugin is a plugin for fastjet (v2.1 upwards) that provides an interface to the seedless infrared safe cone jet finder by Gregory Soyez and Gavin Salam. More...

#include <SISConePlugin.hh>

Inheritance diagram for fastjet::SISConePlugin:

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

Collaboration graph
[legend]
List of all members.

Public Types

enum  SplitMergeScale { SM_pt, SM_Et, SM_mt, SM_pttilde }
 enum for the different split-merge scale choices; Note that order _must_ be the same as in siscone More...

Public Member Functions

 SISConePlugin (double cone_radius, double overlap_threshold=0.5, int n_pass_max=0, double protojet_ptmin=0.0, bool caching=false, SplitMergeScale split_merge_scale=SM_pttilde)
 Constructor for the SISCone Plugin class.
 SISConePlugin (double cone_radius, double overlap_threshold, int n_pass_max, double protojet_ptmin, bool caching, bool split_merge_on_transverse_mass)
 Backwards compatible constructor for the SISCone Plugin class.
 SISConePlugin (double cone_radius, double overlap_threshold, int n_pass_max, bool caching)
 backwards compatible constructor for the SISCone Plugin class (avoid using this in future).
 SISConePlugin (const SISConePlugin &plugin)
 copy constructor
double cone_radius () const
 the cone radius
double overlap_threshold () const
 Fraction of overlap energy in a jet above which jets are merged and below which jets are split.
int n_pass_max () const
 the maximum number of passes of stable-cone searching (<=0 is same as infinity).
double protojet_ptmin () const
 minimum pt for a protojet to be considered in the split-merge step of the algorithm
SplitMergeScale split_merge_scale () const
 indicates scale used in split-merge
void set_split_merge_scale (SplitMergeScale sms)
 sets scale used in split-merge
bool split_merge_on_transverse_mass () const
 indicates whether the split-merge orders on transverse mass or not.
void set_split_merge_on_transverse_mass (bool val)
bool caching () const
 indicates whether caching is turned on or not.
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(.


Private Attributes

double _cone_radius
double _overlap_threshold
int _n_pass_max
double _protojet_ptmin
bool _caching
SplitMergeScale _split_merge_scale

Static Private Attributes

static std::auto_ptr< SISConePluginstored_plugin
static std::auto_ptr< std::vector<
PseudoJet > > 
stored_particles
static std::auto_ptr< siscone::Csiscone > stored_siscone

Detailed Description

SISConePlugin is a plugin for fastjet (v2.1 upwards) that provides an interface to the seedless infrared safe cone jet finder by Gregory Soyez and Gavin Salam.

As of 2006-12-26, this plugin is beta, as is the SISCone code itself.

SISCone uses geometrical techniques to exhaustively consider all possible distinct cones. It then finds out which ones are stable and sends the result to the Tevatron Run-II type split-merge procedure for overlapping cones.

Four parameters govern the "physics" of the algorithm:

One parameter governs some internal algorithmic shortcuts:

The final jets can be accessed by requestion the inclusive_jets(...) from the ClusterSequence object. Note that these PseudoJets have their user_index() set to the index of the pass in which they were found (first pass = 0). NB: This does not currently work for jets that consist of a single particle.

For further information on the details of the algorithm see the SISCone paper; for documentation about the implementation, see the siscone/doc/html/index.html file.

Definition at line 75 of file SISConePlugin.hh.


Member Enumeration Documentation

enum fastjet::SISConePlugin::SplitMergeScale
 

enum for the different split-merge scale choices; Note that order _must_ be the same as in siscone

Enumeration values:
SM_pt  transverse momentum (E-scheme), IR unsafe
SM_Et  transverse energy (E-scheme), not long.

boost invariant original run-II choice [may not be implemented]

SM_mt  transverse mass (E-scheme), IR safe except in decays of two identical narrow heavy particles
SM_pttilde  pt-scheme pt = {i in jet} |p_{ti}|, should be IR safe in all cases

Definition at line 80 of file SISConePlugin.hh.

00080                        {SM_pt,     
00081                         SM_Et,     
00082 
00083                         SM_mt,     
00084 
00085                         SM_pttilde 
00086 
00087   };


Constructor & Destructor Documentation

fastjet::SISConePlugin::SISConePlugin double  cone_radius,
double  overlap_threshold = 0.5,
int  n_pass_max = 0,
double  protojet_ptmin = 0.0,
bool  caching = false,
SplitMergeScale  split_merge_scale = SM_pttilde
[inline]
 

Constructor for the SISCone Plugin class.

Definition at line 91 of file SISConePlugin.hh.

Referenced by run_clustering().

fastjet::SISConePlugin::SISConePlugin double  cone_radius,
double  overlap_threshold,
int  n_pass_max,
double  protojet_ptmin,
bool  caching,
bool  split_merge_on_transverse_mass
[inline]
 

Backwards compatible constructor for the SISCone Plugin class.

Definition at line 105 of file SISConePlugin.hh.

fastjet::SISConePlugin::SISConePlugin double  cone_radius,
double  overlap_threshold,
int  n_pass_max,
bool  caching
[inline]
 

backwards compatible constructor for the SISCone Plugin class (avoid using this in future).

Definition at line 120 of file SISConePlugin.hh.

fastjet::SISConePlugin::SISConePlugin const SISConePlugin plugin  )  [inline]
 

copy constructor

Definition at line 132 of file SISConePlugin.hh.

00132                                                {
00133     *this = plugin;
00134   }


Member Function Documentation

bool fastjet::SISConePlugin::caching  )  const [inline]
 

indicates whether caching is turned on or not.

Definition at line 165 of file SISConePlugin.hh.

Referenced by description(), and run_clustering().

00165 {return _caching ;}

double fastjet::SISConePlugin::cone_radius  )  const [inline]
 

the cone radius

Definition at line 137 of file SISConePlugin.hh.

Referenced by description(), and run_clustering().

00137 {return _cone_radius        ;}

string fastjet::SISConePlugin::description  )  const [virtual]
 

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

Implements fastjet::JetDefinition::Plugin.

Definition at line 23 of file SISConePlugin.cc.

References caching(), cone_radius(), n_pass_max(), overlap_threshold(), protojet_ptmin(), and split_merge_scale().

00023                                          {
00024   ostringstream desc;
00025   
00026   const string on = "on";
00027   const string off = "off";
00028 
00029   string sm_scale_string = "split-merge uses " + 
00030     split_merge_scale_name(Esplit_merge_scale(split_merge_scale()));
00031 
00032   desc << "SISCone jet algorithm with " ;
00033   desc << "cone_radius = "       << cone_radius        () << ", ";
00034   desc << "overlap_threshold = " << overlap_threshold  () << ", ";
00035   desc << "n_pass_max = "        << n_pass_max         () << ", ";
00036   desc << "protojet_ptmin = "    << protojet_ptmin()      << ", ";
00037   desc <<  sm_scale_string                                << ", ";
00038   desc << "caching turned "      << (caching() ? on : off);
00039 
00040   // create a fake scones object so that we can find out more about it
00041   Csiscone siscone;
00042   if (siscone.merge_identical_protocones) {
00043     desc << ", and (IR unsafe) merge_indentical_protocones=true" ;
00044   }
00045 
00046   return desc.str();
00047 }

int fastjet::SISConePlugin::n_pass_max  )  const [inline]
 

the maximum number of passes of stable-cone searching (<=0 is same as infinity).

Definition at line 145 of file SISConePlugin.hh.

Referenced by description(), and run_clustering().

00145 {return _n_pass_max  ;}

double fastjet::SISConePlugin::overlap_threshold  )  const [inline]
 

Fraction of overlap energy in a jet above which jets are merged and below which jets are split.

Definition at line 141 of file SISConePlugin.hh.

Referenced by description(), and run_clustering().

00141 {return _overlap_threshold  ;}

double fastjet::SISConePlugin::protojet_ptmin  )  const [inline]
 

minimum pt for a protojet to be considered in the split-merge step of the algorithm

Definition at line 149 of file SISConePlugin.hh.

Referenced by description(), and run_clustering().

00149 {return _protojet_ptmin  ;}

void fastjet::SISConePlugin::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 57 of file SISConePlugin.cc.

References fastjet::SISConeExtras::_jet_def_plugin, fastjet::SISConeExtras::_most_ambiguous_split, fastjet::SISConeExtras::_protocones, caching(), cone_radius(), fastjet::PseudoJet::E(), fastjet::have_same_momentum(), fastjet::ClusterSequence::jets(), n_pass_max(), overlap_threshold(), fastjet::ClusterSequence::plugin_associate_extras(), fastjet::ClusterSequence::plugin_record_iB_recombination(), fastjet::ClusterSequence::plugin_record_ij_recombination(), protojet_ptmin(), fastjet::PseudoJet::px(), fastjet::PseudoJet::py(), fastjet::PseudoJet::pz(), fastjet::PseudoJet::set_user_index(), SISConePlugin(), split_merge_scale(), stored_particles, stored_plugin, and stored_siscone.

00057                                                                     {
00058 
00059 
00060   Csiscone * siscone;
00061   Csiscone   local_siscone;
00062 
00063   unsigned n = clust_seq.jets().size();
00064 
00065   bool new_siscone = true; // by default we'll be running it
00066 
00067   if (caching()) {
00068 
00069     // Establish if we have a cached run with the same R, npass and
00070     // particles. If not then do any tidying up / reallocation that's
00071     // necessary for the next round of caching, otherwise just set
00072     // relevant pointers so that we can reuse and old run.
00073     if (stored_siscone.get() != 0) {
00074       new_siscone = !(stored_plugin->cone_radius()   == cone_radius()
00075                       && stored_plugin->n_pass_max() == n_pass_max()  
00076                       && stored_particles->size()    == n);
00077       if (!new_siscone) {
00078         for(unsigned i = 0; i < n; i++) {
00079           // only check momentum because indices will be correctly dealt
00080           // with anyway when extracting the clustering order.
00081           new_siscone |= !have_same_momentum(clust_seq.jets()[i], 
00082                                              (*stored_particles)[i]);
00083         }
00084       }
00085     } 
00086       
00087     // allocate the new siscone, etc., if need be
00088     if (new_siscone) {
00089       stored_siscone  .reset( new Csiscone                           );
00090       stored_particles.reset( new vector<PseudoJet>(clust_seq.jets()));
00091       stored_plugin   .reset( new SISConePlugin(*this)               );
00092     }
00093 
00094     siscone = stored_siscone.get();
00095   } else {
00096     siscone = &local_siscone;
00097   }
00098 
00099   if (new_siscone) {
00100     // transfer fastjet initial particles into the siscone type
00101     vector<Cmomentum> siscone_momenta(n);
00102     for(unsigned i = 0; i < n; i++) {
00103       const PseudoJet & p = clust_seq.jets()[i]; // shorthand
00104       siscone_momenta[i] = Cmomentum(p.px(), p.py(), p.pz(), p.E());
00105     }
00106     
00107     // run the jet finding
00108     //cout << "plg sms: " << split_merge_scale() << endl;
00109     siscone->compute_jets(siscone_momenta, cone_radius(), overlap_threshold(),
00110                           n_pass_max(), protojet_ptmin(), 
00111                           Esplit_merge_scale(split_merge_scale()));
00112   } else {
00113     // just run the overlap part of the jets.
00114     //cout << "plg rcmp sms: " << split_merge_scale() << endl;
00115     siscone->recompute_jets(overlap_threshold(), protojet_ptmin(), 
00116                             Esplit_merge_scale(split_merge_scale()));
00117   }
00118 
00119   // extract the jets [in reverse order -- to get nice ordering in pt at end]
00120   int njet = siscone->jets.size();
00121 
00122   for (int ijet = njet-1; ijet >= 0; ijet--) {
00123     const Cjet & jet = siscone->jets[ijet]; // shorthand
00124     
00125     // Successively merge the particles that make up the cone jet
00126     // until we have all particles in it.  Start off with the zeroth
00127     // particle.
00128     int jet_k = jet.contents[0];
00129     for (unsigned ipart = 1; ipart < jet.contents.size(); ipart++) {
00130       // take the last result of the merge
00131       int jet_i = jet_k;
00132       // and the next element of the jet
00133       int jet_j = jet.contents[ipart];
00134       // and merge them (with a fake dij)
00135       double dij = 0.0;
00136 
00137       // create the new jet by hand so that we can adjust its user index
00138       PseudoJet newjet = clust_seq.jets()[jet_i] + clust_seq.jets()[jet_j];
00139 
00140       // set the user index to be the pass in which the jet was discovered
00141       newjet.set_user_index(jet.pass);
00142         
00143       clust_seq.plugin_record_ij_recombination(jet_i, jet_j, dij, newjet, jet_k);
00144     }
00145     // we have merged all the jet's particles into a single object, so now
00146     // "declare" it to be a beam (inclusive) jet.
00147     // [NB: put a sensible looking d_iB just to be nice...]
00148     double d_iB = clust_seq.jets()[jet_k].perp2();
00149     clust_seq.plugin_record_iB_recombination(jet_k, d_iB);
00150   }
00151 
00152   // now copy the list of protocones into an "extras" objects
00153   SISConeExtras * extras = new SISConeExtras;
00154   for (unsigned ipass = 0; ipass < siscone->protocones_list.size(); ipass++) {
00155     for (unsigned ipc = 0; ipc < siscone->protocones_list[ipass].size(); ipc++) {
00156       double rap = siscone->protocones_list[ipass][ipc].eta;
00157       double phi = siscone->protocones_list[ipass][ipc].phi;
00158       PseudoJet protocone(cos(phi),sin(phi),sinh(rap),cosh(rap));
00159       //PseudoJet protocone(siscone->protocones_list[ipass][ipc]);
00160       protocone.set_user_index(ipass);
00161       extras->_protocones.push_back(protocone);
00162     }
00163   }
00164   extras->_most_ambiguous_split = siscone->most_ambiguous_split;
00165 
00166   // tell it what the jet definition was
00167   extras->_jet_def_plugin = this;
00168 
00169   // give the extras object to the cluster sequence.
00170   clust_seq.plugin_associate_extras(auto_ptr<ClusterSequence::Extras>(extras));
00171 }

void fastjet::SISConePlugin::set_split_merge_on_transverse_mass bool  val  )  [inline]
 

Definition at line 160 of file SISConePlugin.hh.

00160                                                     {
00161     _split_merge_scale = val  ? SM_mt : SM_pt;}

void fastjet::SISConePlugin::set_split_merge_scale SplitMergeScale  sms  )  [inline]
 

sets scale used in split-merge

Definition at line 155 of file SISConePlugin.hh.

00155 {_split_merge_scale = sms;}

bool fastjet::SISConePlugin::split_merge_on_transverse_mass  )  const [inline]
 

indicates whether the split-merge orders on transverse mass or not.

retained for backwards compatibility with 2.1.0b3

Definition at line 159 of file SISConePlugin.hh.

00159 {return _split_merge_scale == SM_mt ;}

SplitMergeScale fastjet::SISConePlugin::split_merge_scale  )  const [inline]
 

indicates scale used in split-merge

Definition at line 153 of file SISConePlugin.hh.

Referenced by description(), and run_clustering().

00153 {return _split_merge_scale;}


Member Data Documentation

bool fastjet::SISConePlugin::_caching [private]
 

Definition at line 175 of file SISConePlugin.hh.

double fastjet::SISConePlugin::_cone_radius [private]
 

Definition at line 172 of file SISConePlugin.hh.

int fastjet::SISConePlugin::_n_pass_max [private]
 

Definition at line 173 of file SISConePlugin.hh.

double fastjet::SISConePlugin::_overlap_threshold [private]
 

Definition at line 172 of file SISConePlugin.hh.

double fastjet::SISConePlugin::_protojet_ptmin [private]
 

Definition at line 174 of file SISConePlugin.hh.

SplitMergeScale fastjet::SISConePlugin::_split_merge_scale [private]
 

Definition at line 176 of file SISConePlugin.hh.

auto_ptr< vector< PseudoJet > > fastjet::SISConePlugin::stored_particles [static, private]
 

Definition at line 20 of file SISConePlugin.cc.

Referenced by run_clustering().

auto_ptr< SISConePlugin > fastjet::SISConePlugin::stored_plugin [static, private]
 

Definition at line 19 of file SISConePlugin.cc.

Referenced by run_clustering().

auto_ptr< Csiscone > fastjet::SISConePlugin::stored_siscone [static, private]
 

Definition at line 21 of file SISConePlugin.cc.

Referenced by run_clustering().


The documentation for this class was generated from the following files:
Generated on Mon Apr 2 20:58:23 2007 for fastjet by  doxygen 1.4.2