fastjet 2.4.5
Public Types | Public Member Functions | Protected Member Functions | Private Attributes | Static Private Attributes
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
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.


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.

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

                                         {
  ostringstream desc;
  
  const string on = "on";
  const string off = "off";

  string sm_scale_string = "split-merge uses " + 
    split_merge_scale_name(Esplit_merge_scale(split_merge_scale()));

  desc << "SISCone jet algorithm with " ;
  desc << "cone_radius = "       << cone_radius        () << ", ";
  desc << "overlap_threshold = " << overlap_threshold  () << ", ";
  desc << "n_pass_max = "        << n_pass_max         () << ", ";
  desc << "protojet_ptmin = "    << protojet_ptmin()      << ", ";
  desc <<  sm_scale_string                                << ", ";
  desc << "caching turned "      << (caching() ? on : off);
  desc << ", SM stop scale = "     << _split_merge_stopping_scale;

  // add a note to the description if we use the pt-weighted splitting
  if (_use_pt_weighted_splitting){
    desc << ", using pt-weighted splitting";
  }

  if (_use_jet_def_recombiner){
    desc << ", using jet-definition's own recombiner";
  }

  // create a fake siscone object so that we can find out more about it
  Csiscone siscone;
  if (siscone.merge_identical_protocones) {
    desc << ", and (IR unsafe) merge_indentical_protocones=true" ;
  }

  desc << ", SISCone code v" << siscone_version();

  return desc.str();
}
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.

                                           {return std::max(_protojet_ptmin,
                                                            _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.

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

                                             {
  stored_plugin.reset( new SISConePlugin(*this));
}
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::PseudoJet::E(), fastjet::have_same_momentum(), fastjet::ClusterSequence::jets(), fastjet::ClusterSequence::plugin_associate_extras(), fastjet::ClusterSequence::plugin_record_iB_recombination(), fastjet::ClusterSequence::plugin_record_ij_recombination(), fastjet::PseudoJet::px(), fastjet::PseudoJet::py(), fastjet::PseudoJet::pz(), and fastjet::PseudoJet::set_user_index().

                                                                    {

  Csiscone   local_siscone;
  Csiscone * siscone;
  
  unsigned n = clust_seq.jets().size();

  bool new_siscone = true; // by default we'll be running it

  if (caching()) {

    // Establish if we have a cached run with the same R, npass and
    // particles. If not then do any tidying up / reallocation that's
    // necessary for the next round of caching, otherwise just set
    // relevant pointers so that we can reuse and old run.
    if (stored_siscone.get() != 0) {
      new_siscone = !(stored_plugin->cone_radius()   == cone_radius()
                      && stored_plugin->n_pass_max() == n_pass_max()  
                      && stored_particles->size()    == n);
      if (!new_siscone) {
        for(unsigned i = 0; i < n; i++) {
          // only check momentum because indices will be correctly dealt
          // with anyway when extracting the clustering order.
          new_siscone |= !have_same_momentum(clust_seq.jets()[i], 
                                             (*stored_particles)[i]);
        }
      }
    } 
      
    // allocate the new siscone, etc., if need be
    if (new_siscone) {
      stored_siscone  .reset( new Csiscone );
      stored_particles.reset( new std::vector<PseudoJet>(clust_seq.jets()));
      reset_stored_plugin();
    }

    siscone = stored_siscone.get();
  } else {
    siscone = &local_siscone;
  }

  // make sure stopping scale is set in siscone
  siscone->SM_var2_hardest_cut_off = _split_merge_stopping_scale*_split_merge_stopping_scale;

  // set the specific parameters
  // when running with ghosts for passive areas, do not put the
  // ghosts into the stable-cone search (not relevant)
  siscone->stable_cone_soft_pt2_cutoff = ghost_separation_scale()
                                         * ghost_separation_scale();
  // set the type of splitting we want (default=std one, true->pt-weighted split)
  siscone->set_pt_weighted_splitting(_use_pt_weighted_splitting);

  if (new_siscone) {
    // transfer fastjet initial particles into the siscone type
    std::vector<Cmomentum> siscone_momenta(n);
    for(unsigned i = 0; i < n; i++) {
      const PseudoJet & p = clust_seq.jets()[i]; // shorthand
      siscone_momenta[i] = Cmomentum(p.px(), p.py(), p.pz(), p.E());
    }
    
    // run the jet finding
    //cout << "plg sms: " << split_merge_scale() << endl;
    siscone->compute_jets(siscone_momenta, cone_radius(), overlap_threshold(),
                          n_pass_max(), protojet_or_ghost_ptmin(), 
                          Esplit_merge_scale(split_merge_scale()));
  } else {
    // rerun the jet finding
    // just run the overlap part of the jets.
    //cout << "plg rcmp sms: " << split_merge_scale() << endl;
    siscone->recompute_jets(overlap_threshold(), protojet_or_ghost_ptmin(), 
                            Esplit_merge_scale(split_merge_scale()));
  }

  // extract the jets [in reverse order -- to get nice ordering in pt at end]
  int njet = siscone->jets.size();

  // allocate space for the extras object
  SISConeExtras * extras = new SISConeExtras(n);

  for (int ijet = njet-1; ijet >= 0; ijet--) {
    const Cjet & jet = siscone->jets[ijet]; // shorthand
    
    // Successively merge the particles that make up the cone jet
    // until we have all particles in it.  Start off with the zeroth
    // particle.
    int jet_k = jet.contents[0];
    for (unsigned ipart = 1; ipart < jet.contents.size(); ipart++) {
      // take the last result of the merge
      int jet_i = jet_k;
      // and the next element of the jet
      int jet_j = jet.contents[ipart];
      // and merge them (with a fake dij)
      double dij = 0.0;

      if (_use_jet_def_recombiner) {
        clust_seq.plugin_record_ij_recombination(jet_i, jet_j, dij, jet_k);
      } else {
        // create the new jet by hand so that we can adjust its user index
        PseudoJet newjet = clust_seq.jets()[jet_i] + clust_seq.jets()[jet_j];
        // set the user index to be the pass in which the jet was discovered
        newjet.set_user_index(jet.pass);
        clust_seq.plugin_record_ij_recombination(jet_i, jet_j, dij, newjet, jet_k);
      }

    }

    // we have merged all the jet's particles into a single object, so now
    // "declare" it to be a beam (inclusive) jet.
    // [NB: put a sensible looking d_iB just to be nice...]
    double d_iB = clust_seq.jets()[jet_k].perp2();
    clust_seq.plugin_record_iB_recombination(jet_k, d_iB);

    // now record the pass of the jet in the extras object
    extras->_pass[clust_seq.jets()[jet_k].cluster_hist_index()] = jet.pass;
  }

  // now copy the list of protocones into an "extras" objects
  for (unsigned ipass = 0; ipass < siscone->protocones_list.size(); ipass++) {
    for (unsigned ipc = 0; ipc < siscone->protocones_list[ipass].size(); ipc++) {
      PseudoJet protocone(siscone->protocones_list[ipass][ipc]);
      protocone.set_user_index(ipass);
      extras->_protocones.push_back(protocone);
    }
  }
  extras->_most_ambiguous_split = siscone->most_ambiguous_split;

  // tell it what the jet definition was
  extras->_jet_def_plugin = this;

  // give the extras object to the cluster sequence.
  clust_seq.plugin_associate_extras(std::auto_ptr<ClusterSequence::Extras>(extras));
}
void fastjet::SISConePlugin::set_split_merge_on_transverse_mass ( bool  val) [inline]

Definition at line 160 of file SISConePlugin.hh.

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

sets scale used in split-merge

Definition at line 155 of file SISConePlugin.hh.

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

Definition at line 166 of file SISConePlugin.hh.

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.

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

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.

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

Definition at line 185 of file SISConePlugin.hh.

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

Definition at line 184 of file SISConePlugin.hh.

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

Definition at line 186 of file SISConePlugin.hh.


The documentation for this class was generated from the following files:
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines