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:

  • plugin_do_ij_recombination(.

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.

00066     : _radius(radius), _seedPt(seedPt), _f(f){}

fastjet::ATLASConePlugin::ATLASConePlugin ( const ATLASConePlugin plugin  )  [inline]

copy constructor

Definition at line 69 of file ATLASConePlugin.hh.

00069                                                    {
00070     *this = plugin;
00071   }


Member Function Documentation

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

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

Implements JetDefinition::Plugin.

Definition at line 48 of file ATLASConePlugin.cc.

References _f, _radius, and _seedPt.

00048                                            {
00049   ostringstream desc;
00050   desc << "ATLASCone plugin with R = "<< _radius 
00051        << ", seed threshold = " << _seedPt
00052        << ", overlap threshold f = " << _f;
00053   return desc.str();
00054 }

double fastjet::ATLASConePlugin::f (  )  const [inline]

split-merge overlap threshold

Definition at line 86 of file ATLASConePlugin.hh.

00086 {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 JetDefinition::Plugin.

Definition at line 79 of file ATLASConePlugin.hh.

00079 {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 JetDefinition::Plugin.

Definition at line 56 of file ATLASConePlugin.cc.

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

00056                                                                       {
00057 
00058   // transfer the list of PseudoJet into a atlas::Jet::jet_list_t
00059   //  jet_list_t is a vector<Jet*>
00060   // We set the index of the 4-vect to trace the constituents at the end
00061   //------------------------------------------------------------------
00062   // cout << "ATLASConePlugin: transferring vectors from ClusterSequence" << endl;
00063   atlas::JetConeFinderTool::jetcollection_t jets_ptr;
00064   vector<atlas::Jet*> particles_ptr;
00065 
00066   for (unsigned int i=0 ; i<clust_seq.jets().size() ; i++) {
00067     const PseudoJet & mom = clust_seq.jets()[i];
00068     
00069     // first create the particle
00070     atlas::Jet *particle = new atlas::Jet(mom.px(), mom.py(), mom.pz(), mom.E(), i);
00071     particles_ptr.push_back(particle);
00072 
00073     // then add it to the list of particles we'll use for teh clustering
00074     atlas::Jet *jet = new atlas::Jet;
00075     jet->set_index(particle->index());
00076     jet->addConstituent(particle);
00077 
00078     // and finally add that one to the list of jets
00079     jets_ptr.push_back(jet);
00080   }
00081   // cout << "ATLASCone: " << jets_ptr.size() << " particles to cluster" << endl;
00082 
00083   // search the stable cones
00084   //------------------------------------------------------------------
00085   // cout << "ATLASConePlugin: searching for stable cones" << endl;
00086   atlas::JetConeFinderTool stable_cone_finder;
00087 
00088   // set the parameters
00089   stable_cone_finder.m_coneR  = _radius;
00090   stable_cone_finder.m_seedPt = _seedPt;
00091 
00092   // really do the search.
00093   // Note that the list of protocones is returned 
00094   // through the argument
00095   stable_cone_finder.execute(jets_ptr);
00096   // cout << "ATLASCone: " << jets_ptr.size() << " stable cones found" << endl;
00097 
00098   // perform the split-merge
00099   //------------------------------------------------------------------
00100   // cout << "ATLASConePlugin: running the split-merge" << endl;
00101   atlas::JetSplitMergeTool split_merge;
00102   
00103   // set the parameters
00104   split_merge.m_f = _f;
00105 
00106   // do the work
00107   // again, the list of jets is returned through the argument
00108   split_merge.execute(&jets_ptr);
00109   // cout << "ATLASCone: " << jets_ptr.size() << " jets after split--merge" << endl;
00110 
00111   // build the FastJet jets (a la SISConePlugin)
00112   //------------------------------------------------------------------
00113   // cout << "ATLASConePlugin: backporting jets to the ClusterSequence" << endl;
00114   for (atlas::Jet::jet_list_t::iterator jet_it = jets_ptr.begin() ;
00115          jet_it != jets_ptr.end(); jet_it++){
00116     // iterate over the constituents, starting from the first one
00117     // that we just take as a reference
00118     atlas::Jet::constit_vect_t::iterator constit_it = (*jet_it)->firstConstituent();
00119     // cout << " atlas: jet has " << (*jet_it)->getConstituentNum() << " constituents" << endl;
00120     int jet_k = (*constit_it)->index();
00121     constit_it++;
00122     
00123     // loop over the remaining particles
00124     while (constit_it != (*jet_it)->lastConstituent()){
00125       // take the last result of the merge
00126       int jet_i = jet_k;
00127       // and the next element of the jet
00128       int jet_j = (*constit_it)->index();
00129       // and merge them (with a fake dij)
00130       double dij = 0.0;
00131 
00132       // create the new jet by hand so that we can adjust its user index
00133       // Note again the use of the E-scheme recombination here!
00134       PseudoJet newjet = clust_seq.jets()[jet_i] + clust_seq.jets()[jet_j];
00135       clust_seq.plugin_record_ij_recombination(jet_i, jet_j, dij, newjet, jet_k);
00136 
00137       // jump to the next constituent
00138       constit_it++;
00139     }
00140 
00141     // we have merged all the jet's particles into a single object, so now
00142     // "declare" it to be a beam (inclusive) jet.
00143     // [NB: put a sensible looking d_iB just to be nice...]
00144     double d_iB = clust_seq.jets()[jet_k].perp2();
00145     clust_seq.plugin_record_iB_recombination(jet_k, d_iB);
00146   }
00147     
00148   // cout << "ATLASConePlugin: Bye" << endl;
00149   clear_list(particles_ptr);
00150 }

double fastjet::ATLASConePlugin::seedPt (  )  const [inline]

seed threshold

Definition at line 83 of file ATLASConePlugin.hh.

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

Referenced by description(), and run_clustering().

the cone radius

Definition at line 90 of file ATLASConePlugin.hh.

Referenced by description(), and run_clustering().

the pt seed threshold used in stable-cone search

Definition at line 91 of file ATLASConePlugin.hh.

Referenced by description(), and 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