fastjet 2.4.5
Public Member Functions | Private Attributes
fastjet::ATLASConePlugin Class Reference

ATLASConePlugin is a plugin for fastjet (v2.4 upwards) More...

#include <ATLASConePlugin.hh>

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

List of all members.

Public Member Functions

 ATLASConePlugin (double radius, double seedPt=2.0, double f=0.5)
 Main constructor for the ATLASCone Plugin class.
 ATLASConePlugin (const ATLASConePlugin &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 here we return the R of the last alg in the list
double seedPt () const
 seed threshold
double f () const
 split-merge overlap threshold

Private Attributes

double _radius
 the cone radius
double _seedPt
 the pt seed threshold used in stable-cone search
double _f
 the overlap thresholod used in the split-merge

Detailed Description

ATLASConePlugin is a plugin for fastjet (v2.4 upwards)

Definition at line 51 of file ATLASConePlugin.hh.


Constructor & Destructor Documentation

fastjet::ATLASConePlugin::ATLASConePlugin ( double  radius,
double  seedPt = 2.0,
double  f = 0.5 
) [inline]

Main constructor for the ATLASCone Plugin class.

Apparently the default parameters in the ATLAS software are the ones used here. SpartyJet uses a radius of 0.7, a seed threshold of 1 GeV and an overlap threshold of 0.75 For the ATLAS SW defaults, see http://atlas-sw.cern.ch/cgi-bin/viewcvs-atlas.cgi/groups/JetRoutines/SpartyJet/atlas/ in the JetdoneFinderTools.cxx (rev1.1) and JetSplitMergeTool.cxx (rev1.1) For SpartyJet, see atlas/ConeFinderTool.h

Finally, to agree with FastJet standards, we do not specify a default R, that in the ATLAS code is 0.7

Definition at line 65 of file ATLASConePlugin.hh.

    : _radius(radius), _seedPt(seedPt), _f(f){}
fastjet::ATLASConePlugin::ATLASConePlugin ( const ATLASConePlugin plugin) [inline]

copy constructor

Definition at line 69 of file ATLASConePlugin.hh.

                                                   {
    *this = plugin;
  }

Member Function Documentation

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

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

Implements fastjet::JetDefinition::Plugin.

Definition at line 48 of file ATLASConePlugin.cc.

                                           {
  ostringstream desc;
  desc << "ATLASCone plugin with R = "<< _radius 
       << ", seed threshold = " << _seedPt
       << ", overlap threshold f = " << _f;
  return desc.str();
}
double fastjet::ATLASConePlugin::f ( ) const [inline]

split-merge overlap threshold

Definition at line 86 of file ATLASConePlugin.hh.

{return _f;}
virtual double fastjet::ATLASConePlugin::R ( ) const [inline, virtual]

the plugin mechanism's standard way of accessing the jet radius here we return the R of the last alg in the list

Implements fastjet::JetDefinition::Plugin.

Definition at line 79 of file ATLASConePlugin.hh.

{return _radius;}
void fastjet::ATLASConePlugin::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 56 of file ATLASConePlugin.cc.

References fastjet::atlas::Jet::addConstituent(), fastjet::atlas::clear_list(), fastjet::PseudoJet::E(), fastjet::atlas::JetSplitMergeTool::execute(), fastjet::atlas::JetConeFinderTool::execute(), fastjet::atlas::Jet::index(), fastjet::ClusterSequence::jets(), fastjet::atlas::JetConeFinderTool::m_coneR, fastjet::atlas::JetSplitMergeTool::m_f, fastjet::atlas::JetConeFinderTool::m_seedPt, fastjet::ClusterSequence::plugin_record_iB_recombination(), fastjet::ClusterSequence::plugin_record_ij_recombination(), fastjet::PseudoJet::px(), fastjet::PseudoJet::py(), fastjet::PseudoJet::pz(), and fastjet::atlas::Jet::set_index().

                                                                      {

  // transfer the list of PseudoJet into a atlas::Jet::jet_list_t
  //  jet_list_t is a vector<Jet*>
  // We set the index of the 4-vect to trace the constituents at the end
  //------------------------------------------------------------------
  // cout << "ATLASConePlugin: transferring vectors from ClusterSequence" << endl;
  atlas::JetConeFinderTool::jetcollection_t jets_ptr;
  vector<atlas::Jet*> particles_ptr;

  for (unsigned int i=0 ; i<clust_seq.jets().size() ; i++) {
    const PseudoJet & mom = clust_seq.jets()[i];
    
    // first create the particle
    atlas::Jet *particle = new atlas::Jet(mom.px(), mom.py(), mom.pz(), mom.E(), i);
    particles_ptr.push_back(particle);

    // then add it to the list of particles we'll use for teh clustering
    atlas::Jet *jet = new atlas::Jet;
    jet->set_index(particle->index());
    jet->addConstituent(particle);

    // and finally add that one to the list of jets
    jets_ptr.push_back(jet);
  }
  // cout << "ATLASCone: " << jets_ptr.size() << " particles to cluster" << endl;

  // search the stable cones
  //------------------------------------------------------------------
  // cout << "ATLASConePlugin: searching for stable cones" << endl;
  atlas::JetConeFinderTool stable_cone_finder;

  // set the parameters
  stable_cone_finder.m_coneR  = _radius;
  stable_cone_finder.m_seedPt = _seedPt;

  // really do the search.
  // Note that the list of protocones is returned 
  // through the argument
  stable_cone_finder.execute(jets_ptr);
  // cout << "ATLASCone: " << jets_ptr.size() << " stable cones found" << endl;

  // perform the split-merge
  //------------------------------------------------------------------
  // cout << "ATLASConePlugin: running the split-merge" << endl;
  atlas::JetSplitMergeTool split_merge;
  
  // set the parameters
  split_merge.m_f = _f;

  // do the work
  // again, the list of jets is returned through the argument
  split_merge.execute(&jets_ptr);
  // cout << "ATLASCone: " << jets_ptr.size() << " jets after split--merge" << endl;

  // build the FastJet jets (a la SISConePlugin)
  //------------------------------------------------------------------
  // cout << "ATLASConePlugin: backporting jets to the ClusterSequence" << endl;
  for (atlas::Jet::jet_list_t::iterator jet_it = jets_ptr.begin() ;
         jet_it != jets_ptr.end(); jet_it++){
    // iterate over the constituents, starting from the first one
    // that we just take as a reference
    atlas::Jet::constit_vect_t::iterator constit_it = (*jet_it)->firstConstituent();
    // cout << " atlas: jet has " << (*jet_it)->getConstituentNum() << " constituents" << endl;
    int jet_k = (*constit_it)->index();
    constit_it++;
    
    // loop over the remaining particles
    while (constit_it != (*jet_it)->lastConstituent()){
      // take the last result of the merge
      int jet_i = jet_k;
      // and the next element of the jet
      int jet_j = (*constit_it)->index();
      // and merge them (with a fake dij)
      double dij = 0.0;

      // create the new jet by hand so that we can adjust its user index
      // Note again the use of the E-scheme recombination here!
      PseudoJet newjet = clust_seq.jets()[jet_i] + clust_seq.jets()[jet_j];
      clust_seq.plugin_record_ij_recombination(jet_i, jet_j, dij, newjet, jet_k);

      // jump to the next constituent
      constit_it++;
    }

    // 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);
  }
    
  // cout << "ATLASConePlugin: Bye" << endl;
  clear_list(particles_ptr);
  clear_list(jets_ptr);
}
double fastjet::ATLASConePlugin::seedPt ( ) const [inline]

seed threshold

Definition at line 83 of file ATLASConePlugin.hh.

{return _seedPt;}

Member Data Documentation

double fastjet::ATLASConePlugin::_f [private]

the overlap thresholod used in the split-merge

Definition at line 92 of file ATLASConePlugin.hh.

the cone radius

Definition at line 90 of file ATLASConePlugin.hh.

the pt seed threshold used in stable-cone search

Definition at line 91 of file ATLASConePlugin.hh.


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