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, int n_pass_max=0, double protojet_ptmin=0.0, bool caching=false, SplitMergeScale split_merge_scale=SM_pttilde, double split_merge_stopping_scale=0.0)
 Main 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).
double protojet_ptmin () const
 minimum pt for a protojet to be considered in the split-merge step of the algorithm
double protojet_or_ghost_ptmin () const
 return the scale to be passed to SISCone as the protojet_ptmin -- 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_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 split_merge_use_pt_weighted_splitting () const
 indicates whether the split-merge orders on transverse mass or not.
void set_split_merge_use_pt_weighted_splitting (bool val)
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_ptmin
SplitMergeScale _split_merge_scale
bool _use_pt_weighted_splitting

Static Private Attributes

static std::auto_ptr
< SISConePlugin
stored_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.

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 67 of file SISConePlugin.hh.


Member Enumeration Documentation

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

Enumerator:
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 72 of file SISConePlugin.hh.

00072                        {SM_pt,     
00073                         SM_Et,     
00074 
00075                         SM_mt,     
00076 
00077                         SM_pttilde 
00078 
00079   };


Constructor & Destructor Documentation

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

Main constructor for the SISCone Plugin class.

Note: wrt version prior to 2.4 this constructor differs in that a the default value has been removed for overlap_threshold. The former has been removed because the old default of 0.5 was found to be unsuitable in high-noise environments; so the user should now explicitly think about the value for this -- we recommend 0.75.

Definition at line 91 of file SISConePlugin.hh.

Referenced by reset_stored_plugin().

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 110 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 127 of file SISConePlugin.hh.


Member Function Documentation

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

plugin description

Implements fastjet::SISConeBasePlugin.

Definition at line 36 of file SISConePlugin.cc.

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

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

double fastjet::SISConePlugin::protojet_or_ghost_ptmin (  )  const [inline]

return the scale to be passed to SISCone as the protojet_ptmin -- 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 149 of file SISConePlugin.hh.

Referenced by run_clustering().

00149                                            {return std::max(_protojet_ptmin,
00150                                                             _ghost_sep_scale);}

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 143 of file SISConePlugin.hh.

Referenced by description().

00143 {return _protojet_ptmin  ;}

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

call the re-clustering itself

Implements fastjet::SISConeBasePlugin.

Definition at line 210 of file SISConePlugin.cc.

References SISConePlugin(), and stored_plugin.

Referenced by run_clustering().

00210                                              {
00211   stored_plugin.reset( new SISConePlugin(*this));
00212 }

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

really do the clustering work

Implements fastjet::SISConeBasePlugin.

Definition at line 76 of file SISConePlugin.cc.

References fastjet::SISConeBaseExtras::_jet_def_plugin, fastjet::SISConeBaseExtras::_most_ambiguous_split, fastjet::SISConeBaseExtras::_pass, fastjet::SISConeBaseExtras::_protocones, fastjet::SISConeBasePlugin::_split_merge_stopping_scale, fastjet::SISConeBasePlugin::_use_jet_def_recombiner, _use_pt_weighted_splitting, 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_ptmin(), PseudoJet::px(), PseudoJet::py(), PseudoJet::pz(), reset_stored_plugin(), PseudoJet::set_user_index(), split_merge_scale(), stored_particles, stored_plugin, and stored_siscone.

00076                                                                     {
00077 
00078   Csiscone   local_siscone;
00079   Csiscone * siscone;
00080   
00081   unsigned n = clust_seq.jets().size();
00082 
00083   bool new_siscone = true; // by default we'll be running it
00084 
00085   if (caching()) {
00086 
00087     // Establish if we have a cached run with the same R, npass and
00088     // particles. If not then do any tidying up / reallocation that's
00089     // necessary for the next round of caching, otherwise just set
00090     // relevant pointers so that we can reuse and old run.
00091     if (stored_siscone.get() != 0) {
00092       new_siscone = !(stored_plugin->cone_radius()   == cone_radius()
00093                       && stored_plugin->n_pass_max() == n_pass_max()  
00094                       && stored_particles->size()    == n);
00095       if (!new_siscone) {
00096         for(unsigned i = 0; i < n; i++) {
00097           // only check momentum because indices will be correctly dealt
00098           // with anyway when extracting the clustering order.
00099           new_siscone |= !have_same_momentum(clust_seq.jets()[i], 
00100                                              (*stored_particles)[i]);
00101         }
00102       }
00103     } 
00104       
00105     // allocate the new siscone, etc., if need be
00106     if (new_siscone) {
00107       stored_siscone  .reset( new Csiscone );
00108       stored_particles.reset( new std::vector<PseudoJet>(clust_seq.jets()));
00109       reset_stored_plugin();
00110     }
00111 
00112     siscone = stored_siscone.get();
00113   } else {
00114     siscone = &local_siscone;
00115   }
00116 
00117   // make sure stopping scale is set in siscone
00118   siscone->SM_var2_hardest_cut_off = _split_merge_stopping_scale*_split_merge_stopping_scale;
00119 
00120   // set the specific parameters
00121   // when running with ghosts for passive areas, do not put the
00122   // ghosts into the stable-cone search (not relevant)
00123   siscone->stable_cone_soft_pt2_cutoff = ghost_separation_scale()
00124                                          * ghost_separation_scale();
00125   // set the type of splitting we want (default=std one, true->pt-weighted split)
00126   siscone->set_pt_weighted_splitting(_use_pt_weighted_splitting);
00127 
00128   if (new_siscone) {
00129     // transfer fastjet initial particles into the siscone type
00130     std::vector<Cmomentum> siscone_momenta(n);
00131     for(unsigned i = 0; i < n; i++) {
00132       const PseudoJet & p = clust_seq.jets()[i]; // shorthand
00133       siscone_momenta[i] = Cmomentum(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_ptmin(), 
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_ptmin(), 
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   SISConeExtras * extras = new SISConeExtras(n);
00154 
00155   for (int ijet = njet-1; ijet >= 0; ijet--) {
00156     const Cjet & 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::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;}

void fastjet::SISConePlugin::set_split_merge_use_pt_weighted_splitting ( bool  val  )  [inline]

Definition at line 166 of file SISConePlugin.hh.

00166                                                            {
00167     _use_pt_weighted_splitting = val;}

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;}

bool fastjet::SISConePlugin::split_merge_use_pt_weighted_splitting (  )  const [inline]

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

retained for backwards compatibility with 2.1.0b3

Definition at line 165 of file SISConePlugin.hh.


Member Data Documentation

Definition at line 177 of file SISConePlugin.hh.

Definition at line 178 of file SISConePlugin.hh.

Definition at line 180 of file SISConePlugin.hh.

Referenced by description(), and run_clustering().

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

Definition at line 185 of file SISConePlugin.hh.

Referenced by run_clustering().

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

Definition at line 184 of file SISConePlugin.hh.

Referenced by reset_stored_plugin(), and run_clustering().

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

Definition at line 186 of file SISConePlugin.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