fastjet::SISConeSphericalPlugin Class Reference

SISConeSphericalPlugin 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 <SISConeSphericalPlugin.hh>

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

List of all members.

Public Types

enum  SplitMergeScale { SM_E, SM_Etilde }
 

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

More...

Public Member Functions

 SISConeSphericalPlugin (double cone_radius, double overlap_threshold, int n_pass_max=0, double protojet_Emin=0.0, bool caching=false, SplitMergeScale split_merge_scale=SM_Etilde, double split_merge_stopping_scale=0.0)
 Main constructor for the SISConeSpherical Plugin class.
double protojet_Emin () const
 minimum energy for a protojet to be considered in the split-merge step of the algorithm
double protojet_or_ghost_Emin () const
 return the scale to be passed to SISCone as the protojet_Emin -- if we have a ghost separation scale that is above the protojet_ptmin, then the ghost_separation_scale becomes the relevant one to use here
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_use_E_weighted_splitting () const
 indicate if the splittings are done using the anti-kt distance
void set_split_merge_use_E_weighted_splitting (bool val)
virtual bool supports_ghosted_passive_areas () const
 overload the default as we don't provide support for passive areas.
virtual std::string description () const
 plugin description
virtual void run_clustering (ClusterSequence &) const
 really do the clustering work

Protected Member Functions

virtual void reset_stored_plugin () const
 call the re-clustering itself

Private Attributes

double _protojet_Emin
SplitMergeScale _split_merge_scale
bool _use_E_weighted_splitting

Static Private Attributes

static std::auto_ptr
< SISConeSphericalPlugin
stored_plugin
static std::auto_ptr
< std::vector< PseudoJet > > 
stored_particles
static std::auto_ptr
< siscone_spherical::CSphsiscone > 
stored_siscone

Detailed Description

SISConeSphericalPlugin 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.

This is the version of SISCone using spherical coordinates. Compared to the original cylindrical version:

For moderate values of R the problem should not be too severe (or may even be absent for some values of the overlap parameter), however the user should be aware of the issue.

The default split-merge scale may change at a later date to resolve this issue.

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, arXiv:0704.0292 [JHEP 0705:086,2007].

For documentation about the implementation, see the siscone/doc/html/index.html file.

Definition at line 89 of file SISConeSphericalPlugin.hh.


Member Enumeration Documentation

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

Enumerator:
SM_E 

Energy (IR unsafe with momentum conservation).

SM_Etilde 

sum_{i jet} E_i [1+sin^2(theta_iJ)]

Definition at line 94 of file SISConeSphericalPlugin.hh.

00094                        {SM_E,        
00095                         SM_Etilde    
00096   };


Constructor & Destructor Documentation

fastjet::SISConeSphericalPlugin::SISConeSphericalPlugin ( double  cone_radius,
double  overlap_threshold,
int  n_pass_max = 0,
double  protojet_Emin = 0.0,
bool  caching = false,
SplitMergeScale  split_merge_scale = SM_Etilde,
double  split_merge_stopping_scale = 0.0 
) [inline]

Main constructor for the SISConeSpherical Plugin class.

Definition at line 102 of file SISConeSphericalPlugin.hh.

Referenced by reset_stored_plugin().


Member Function Documentation

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

plugin description

Implements fastjet::SISConeBasePlugin.

Definition at line 35 of file SISConeSphericalPlugin.cc.

References fastjet::SISConeBasePlugin::_split_merge_stopping_scale, _use_E_weighted_splitting, fastjet::SISConeBasePlugin::_use_jet_def_recombiner, fastjet::SISConeBasePlugin::caching(), fastjet::SISConeBasePlugin::cone_radius(), fastjet::SISConeBasePlugin::n_pass_max(), fastjet::SISConeBasePlugin::overlap_threshold(), protojet_Emin(), and split_merge_scale().

00035                                                   {
00036   ostringstream desc;
00037   
00038   const string on = "on";
00039   const string off = "off";
00040 
00041   string sm_scale_string = "split-merge uses " + 
00042     split_merge_scale_name(Esplit_merge_scale(split_merge_scale()));
00043 
00044   desc << "Spherical SISCone jet algorithm with " ;
00045   desc << "cone_radius = "       << cone_radius        () << ", ";
00046   desc << "overlap_threshold = " << overlap_threshold  () << ", ";
00047   desc << "n_pass_max = "        << n_pass_max         () << ", ";
00048   desc << "protojet_Emin = "     << protojet_Emin()      << ", ";
00049   desc <<  sm_scale_string                                << ", ";
00050   desc << "caching turned "      << (caching() ? on : off);
00051   desc << ", SM stop scale = "     << _split_merge_stopping_scale;
00052 
00053   // add a note to the description if we use the pt-weighted splitting
00054   if (_use_E_weighted_splitting){
00055     desc << ", using E-weighted splitting";
00056   }
00057 
00058   if (_use_jet_def_recombiner){
00059     desc << ", using jet-definition's own recombiner";
00060   }
00061 
00062   // create a fake siscone object so that we can find out more about it
00063   CSphsiscone siscone;
00064   if (siscone.merge_identical_protocones) {
00065     desc << ", and (IR unsafe) merge_indentical_protocones=true" ;
00066   }
00067 
00068   desc << ", SISCone code v" << siscone_version();
00069 
00070   return desc.str();
00071 }

double fastjet::SISConeSphericalPlugin::protojet_Emin (  )  const [inline]

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

Definition at line 122 of file SISConeSphericalPlugin.hh.

Referenced by description().

00122 {return _protojet_Emin  ;}

double fastjet::SISConeSphericalPlugin::protojet_or_ghost_Emin (  )  const [inline]

return the scale to be passed to SISCone as the protojet_Emin -- if we have a ghost separation scale that is above the protojet_ptmin, then the ghost_separation_scale becomes the relevant one to use here

Definition at line 128 of file SISConeSphericalPlugin.hh.

Referenced by run_clustering().

00128                                           {return std::max(_protojet_Emin,
00129                                                            _ghost_sep_scale);}

void fastjet::SISConeSphericalPlugin::reset_stored_plugin (  )  const [protected, virtual]

call the re-clustering itself

Implements fastjet::SISConeBasePlugin.

Definition at line 209 of file SISConeSphericalPlugin.cc.

References SISConeSphericalPlugin(), and stored_plugin.

Referenced by run_clustering().

00209                                                       {
00210   stored_plugin.reset( new SISConeSphericalPlugin(*this));
00211 }

void fastjet::SISConeSphericalPlugin::run_clustering ( ClusterSequence  )  const [virtual]

really do the clustering work

Implements fastjet::SISConeBasePlugin.

Definition at line 75 of file SISConeSphericalPlugin.cc.

References fastjet::SISConeBaseExtras::_jet_def_plugin, fastjet::SISConeBaseExtras::_most_ambiguous_split, fastjet::SISConeBaseExtras::_pass, fastjet::SISConeBaseExtras::_protocones, fastjet::SISConeBasePlugin::_split_merge_stopping_scale, _use_E_weighted_splitting, fastjet::SISConeBasePlugin::_use_jet_def_recombiner, fastjet::SISConeBasePlugin::caching(), fastjet::SISConeBasePlugin::cone_radius(), PseudoJet::E(), fastjet::SISConeBasePlugin::ghost_separation_scale(), have_same_momentum(), ClusterSequence::jets(), fastjet::SISConeBasePlugin::n_pass_max(), fastjet::SISConeBasePlugin::overlap_threshold(), ClusterSequence::plugin_associate_extras(), ClusterSequence::plugin_record_iB_recombination(), ClusterSequence::plugin_record_ij_recombination(), protojet_or_ghost_Emin(), PseudoJet::px(), PseudoJet::py(), PseudoJet::pz(), reset_stored_plugin(), PseudoJet::set_user_index(), split_merge_scale(), stored_particles, stored_plugin, and stored_siscone.

00075                                                                              {
00076 
00077   CSphsiscone   local_siscone;
00078   CSphsiscone * siscone;
00079 
00080   unsigned n = clust_seq.jets().size();
00081 
00082   bool new_siscone = true; // by default we'll be running it
00083 
00084   if (caching()) {
00085 
00086     // Establish if we have a cached run with the same R, npass and
00087     // particles. If not then do any tidying up / reallocation that's
00088     // necessary for the next round of caching, otherwise just set
00089     // relevant pointers so that we can reuse and old run.
00090     if (stored_siscone.get() != 0) {
00091       new_siscone = !(stored_plugin->cone_radius()   == cone_radius()
00092                       && stored_plugin->n_pass_max() == n_pass_max()
00093                       && stored_particles->size()    == n);
00094       if (!new_siscone) {
00095         for(unsigned i = 0; i < n; i++) {
00096           // only check momentum because indices will be correctly dealt
00097           // with anyway when extracting the clustering order.
00098           new_siscone |= !have_same_momentum(clust_seq.jets()[i],
00099                                              (*stored_particles)[i]);
00100         }
00101       }
00102     }
00103 
00104     // allocate the new siscone, etc., if need be
00105     if (new_siscone) {
00106       stored_siscone  .reset( new CSphsiscone );
00107       stored_particles.reset( new std::vector<PseudoJet>(clust_seq.jets()));
00108       reset_stored_plugin();
00109     }
00110 
00111     siscone = stored_siscone.get();
00112   } else {
00113     siscone = &local_siscone;
00114   }
00115 
00116   // make sure stopping scale is set in siscone
00117   siscone->SM_var2_hardest_cut_off = _split_merge_stopping_scale*_split_merge_stopping_scale;
00118 
00119   // set the specific parameters
00120   // when running with ghosts for passive areas, do not put the
00121   // ghosts into the stable-cone search (not relevant)
00122   siscone->stable_cone_soft_E2_cutoff = ghost_separation_scale()
00123                                       * ghost_separation_scale();
00124   // set the type of splitting we want (default=std one, true->pt-weighted split)
00125   siscone->set_E_weighted_splitting(_use_E_weighted_splitting);
00126 
00127 
00128   if (new_siscone) {
00129     // transfer fastjet initial particles into the siscone type
00130     std::vector<CSphmomentum> siscone_momenta(n);
00131     for(unsigned i = 0; i < n; i++) {
00132       const PseudoJet & p = clust_seq.jets()[i]; // shorthand
00133       siscone_momenta[i] = CSphmomentum(p.px(), p.py(), p.pz(), p.E());
00134     }
00135 
00136     // run the jet finding
00137     //cout << "plg sms: " << split_merge_scale() << endl;
00138     siscone->compute_jets(siscone_momenta, cone_radius(), overlap_threshold(),
00139                           n_pass_max(), protojet_or_ghost_Emin(), 
00140                           Esplit_merge_scale(split_merge_scale()));
00141   } else {
00142     // rerun the jet finding
00143     // just run the overlap part of the jets.
00144     //cout << "plg rcmp sms: " << split_merge_scale() << endl;
00145     siscone->recompute_jets(overlap_threshold(), protojet_or_ghost_Emin(), 
00146                             Esplit_merge_scale(split_merge_scale()));
00147   }
00148 
00149   // extract the jets [in reverse order -- to get nice ordering in pt at end]
00150   int njet = siscone->jets.size();
00151 
00152   // allocate space for the extras object
00153   SISConeSphericalExtras * extras = new SISConeSphericalExtras(n);
00154 
00155   for (int ijet = njet-1; ijet >= 0; ijet--) {
00156     const CSphjet & jet = siscone->jets[ijet]; // shorthand
00157 
00158     // Successively merge the particles that make up the cone jet
00159     // until we have all particles in it.  Start off with the zeroth
00160     // particle.
00161     int jet_k = jet.contents[0];
00162     for (unsigned ipart = 1; ipart < jet.contents.size(); ipart++) {
00163       // take the last result of the merge
00164       int jet_i = jet_k;
00165       // and the next element of the jet
00166       int jet_j = jet.contents[ipart];
00167       // and merge them (with a fake dij)
00168       double dij = 0.0;
00169 
00170       if (_use_jet_def_recombiner) {
00171        clust_seq.plugin_record_ij_recombination(jet_i, jet_j, dij, jet_k);
00172       } else {
00173        // create the new jet by hand so that we can adjust its user index
00174        PseudoJet newjet = clust_seq.jets()[jet_i] + clust_seq.jets()[jet_j];
00175        // set the user index to be the pass in which the jet was discovered
00176        newjet.set_user_index(jet.pass);
00177        clust_seq.plugin_record_ij_recombination(jet_i, jet_j, dij, newjet, jet_k);
00178       }
00179 
00180     }
00181 
00182     // we have merged all the jet's particles into a single object, so now
00183     // "declare" it to be a beam (inclusive) jet.
00184     // [NB: put a sensible looking d_iB just to be nice...]
00185     double d_iB = clust_seq.jets()[jet_k].perp2();
00186     clust_seq.plugin_record_iB_recombination(jet_k, d_iB);
00187 
00188     // now record the pass of the jet in the extras object
00189     extras->_pass[clust_seq.jets()[jet_k].cluster_hist_index()] = jet.pass;
00190   }
00191 
00192   // now copy the list of protocones into an "extras" objects
00193   for (unsigned ipass = 0; ipass < siscone->protocones_list.size(); ipass++) {
00194     for (unsigned ipc = 0; ipc < siscone->protocones_list[ipass].size(); ipc++) {
00195       PseudoJet protocone(siscone->protocones_list[ipass][ipc]);
00196       protocone.set_user_index(ipass);
00197       extras->_protocones.push_back(protocone);
00198     }
00199   }
00200   extras->_most_ambiguous_split = siscone->most_ambiguous_split;
00201 
00202   // tell it what the jet definition was
00203   extras->_jet_def_plugin = this;
00204 
00205   // give the extras object to the cluster sequence.
00206   clust_seq.plugin_associate_extras(std::auto_ptr<ClusterSequence::Extras>(extras));
00207 }

void fastjet::SISConeSphericalPlugin::set_split_merge_scale ( SplitMergeScale  sms  )  [inline]

sets scale used in split-merge

Definition at line 134 of file SISConeSphericalPlugin.hh.

00134 {_split_merge_scale = sms;}

void fastjet::SISConeSphericalPlugin::set_split_merge_use_E_weighted_splitting ( bool  val  )  [inline]

Definition at line 138 of file SISConeSphericalPlugin.hh.

00138                                                           {
00139     _use_E_weighted_splitting = val;}

SplitMergeScale fastjet::SISConeSphericalPlugin::split_merge_scale (  )  const [inline]

indicates scale used in split-merge

Definition at line 132 of file SISConeSphericalPlugin.hh.

Referenced by description(), and run_clustering().

00132 {return _split_merge_scale;}

bool fastjet::SISConeSphericalPlugin::split_merge_use_E_weighted_splitting (  )  const [inline]

indicate if the splittings are done using the anti-kt distance

Definition at line 137 of file SISConeSphericalPlugin.hh.

00137 {return _use_E_weighted_splitting;}

virtual bool fastjet::SISConeSphericalPlugin::supports_ghosted_passive_areas (  )  const [inline, virtual]

overload the default as we don't provide support for passive areas.

Reimplemented from fastjet::SISConeBasePlugin.

Definition at line 143 of file SISConeSphericalPlugin.hh.

00143 {return true;}


Member Data Documentation

Definition at line 153 of file SISConeSphericalPlugin.hh.

Definition at line 154 of file SISConeSphericalPlugin.hh.

Definition at line 155 of file SISConeSphericalPlugin.hh.

Referenced by description(), and run_clustering().

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

Definition at line 160 of file SISConeSphericalPlugin.hh.

Referenced by run_clustering().

Definition at line 159 of file SISConeSphericalPlugin.hh.

Referenced by reset_stored_plugin(), and run_clustering().

std::auto_ptr< CSphsiscone > fastjet::SISConeSphericalPlugin::stored_siscone [static, private]

Definition at line 161 of file SISConeSphericalPlugin.hh.

Referenced by run_clustering().


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

Generated on 26 Feb 2010 for fastjet by  doxygen 1.6.1