Main Page | Namespace List | Class Hierarchy | Class List | Directories | File List | Namespace Members | Class Members | File Members

fastjet::CDFMidPointPlugin Class Reference

CDFMidPointPlugin is a plugin for fastjet (v2.1 upwards) that provides an interface to the CDF version of Run-II iterative cone algorithm with midpoint seeds (also known as the Iterative Legacy Cone Algorithm, ILCA). More...

#include <CDFMidPointPlugin.hh>

Inheritance diagram for fastjet::CDFMidPointPlugin:

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

Collaboration graph
[legend]
List of all members.

Public Types

enum  SplitMergeScale { SM_pt, SM_Et, SM_mt, SM_pttilde }
 the choice of scale to be used in the split-merge step More...

Public Member Functions

 CDFMidPointPlugin (double seed_threshold, double cone_radius, double cone_area_fraction, int max_pair_size, int max_iterations, double overlap_threshold, SplitMergeScale sm_scale=SM_pt)
 A CDFMidPointPlugin constructor that looks like the one provided by CDF.
 CDFMidPointPlugin (double cone_radius, double overlap_threshold=0.5, double seed_threshold=1.0, double cone_area_fraction=1.0)
 a compact constructor
double seed_threshold () const
double cone_radius () const
double cone_area_fraction () const
int max_pair_size () const
int max_iterations () const
double overlap_threshold () const
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(.


Private Attributes

double _seed_threshold
double _cone_radius
double _cone_area_fraction
int _max_pair_size
int _max_iterations
double _overlap_threshold
SplitMergeScale _sm_scale

Detailed Description

CDFMidPointPlugin is a plugin for fastjet (v2.1 upwards) that provides an interface to the CDF version of Run-II iterative cone algorithm with midpoint seeds (also known as the Iterative Legacy Cone Algorithm, ILCA).

The CDF code has been taken from Joey Huston's webpage http://www.pa.msu.edu/~huston/Les_Houches_2005/Les_Houches_SM.html

Note that the CDF midpoint code contains options that go beyond those described in the Tevatron run-II document (hep-ex/0005012), notably search-cones, as described in hep-ph/0111434, and midpoints bewteen multiplets of stable cones.

Additionally, the version of the CDF midpoint code distributed here has been modified by the FastJet authors, so as to allow one to choose the scale used in the split-merge step.

Definition at line 60 of file CDFMidPointPlugin.hh.


Member Enumeration Documentation

enum fastjet::CDFMidPointPlugin::SplitMergeScale
 

the choice of scale to be used in the split-merge step

Enumeration values:
SM_pt 
SM_Et 
SM_mt 
SM_pttilde 

Definition at line 64 of file CDFMidPointPlugin.hh.


Constructor & Destructor Documentation

fastjet::CDFMidPointPlugin::CDFMidPointPlugin double  seed_threshold,
double  cone_radius,
double  cone_area_fraction,
int  max_pair_size,
int  max_iterations,
double  overlap_threshold,
SplitMergeScale  sm_scale = SM_pt
[inline]
 

A CDFMidPointPlugin constructor that looks like the one provided by CDF.

Its arguments should have the following meaning:

  • seed_threshold: minimum pt for a particle to be considered a seed of the iteration.

  • cone_radius: standard meaning

  • cone_area_fraction: stable-cones are searched for with a radius Rsearch = R * sqrt(cone_area_fraction), and then expanded to size R afterwards; note (hep-ph/0610012) that this introduces IR unsafety at NLO for X+2-jet observables (where X any hard object).

  • max_pair_size: "midpoints" can be added between pairs of stable cones, triplets of stable cones, etc.; max_pair_size indicates the maximum number of stable cones that are assembled when adding midpoints.

  • max_iterations: the maximum number of iterations to carry out when looking for a stable cone.

  • overlap_threshold: if (overlapping_Et)/(Et_of_softer_protojet) < overlap_threshold, overlapping jets are split, otherwise they are merged.

  • sm_scale: a choice for the scale to be used in the split-merge step (both for ordering the momenta and quantifying the overlap); the three options are

. SM_pt: pt (default -- source of small IR safety issue in purely hadronic events)

. SM_Et: Et (not boost invariant, reduces to mt at zero rapidity and to pt and infinite rapidity)

. SM_mt: transverse mass = sqrt(m^2+pt^2)

Definition at line 105 of file CDFMidPointPlugin.hh.

fastjet::CDFMidPointPlugin::CDFMidPointPlugin double  cone_radius,
double  overlap_threshold = 0.5,
double  seed_threshold = 1.0,
double  cone_area_fraction = 1.0
[inline]
 

a compact constructor

Definition at line 122 of file CDFMidPointPlugin.hh.


Member Function Documentation

double fastjet::CDFMidPointPlugin::cone_area_fraction  )  const [inline]
 

Definition at line 138 of file CDFMidPointPlugin.hh.

Referenced by description().

00138 {return _cone_area_fraction ;}

double fastjet::CDFMidPointPlugin::cone_radius  )  const [inline]
 

Definition at line 137 of file CDFMidPointPlugin.hh.

Referenced by description().

00137 {return _cone_radius        ;}

string fastjet::CDFMidPointPlugin::description  )  const [virtual]
 

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

Implements fastjet::JetDefinition::Plugin.

Definition at line 45 of file CDFMidPointPlugin.cc.

References _sm_scale, cone_area_fraction(), cone_radius(), max_iterations(), max_pair_size(), overlap_threshold(), seed_threshold(), SM_Et, SM_mt, SM_pt, and SM_pttilde.

00045                                              {
00046   ostringstream desc;
00047   
00048   string sm_scale_string = "split-merge uses ";
00049   switch(_sm_scale) {
00050   case SM_pt:
00051     sm_scale_string += "pt";
00052     break;
00053   case SM_Et:
00054     sm_scale_string += "Et";
00055     break;
00056   case SM_mt:
00057     sm_scale_string += "mt";
00058     break;
00059   case SM_pttilde:
00060     sm_scale_string += "pttilde (scalar sum of pts)";
00061     break;
00062   default:
00063     ostringstream err;
00064     err << "Unrecognized split-merge scale choice = " << _sm_scale;
00065     throw Error(err.str());
00066   }
00067 
00068   
00069   if (cone_area_fraction() == 1) {
00070     desc << "CDF MidPoint jet algorithm, with " ;
00071   } else {
00072     desc << "CDF MidPoint+Searchcone jet algorithm, with ";
00073   }
00074   desc << "seed_threshold = "     << seed_threshold     () << ", "
00075        << "cone_radius = "        << cone_radius        () << ", "
00076        << "cone_area_fraction = " << cone_area_fraction () << ", " 
00077        << "max_pair_size = "      << max_pair_size      () << ", "
00078        << "max_iterations = "     << max_iterations     () << ", "
00079        << "overlap_threshold  = " << overlap_threshold  () << ", "
00080        << sm_scale_string ;
00081 
00082   return desc.str();
00083 }

int fastjet::CDFMidPointPlugin::max_iterations  )  const [inline]
 

Definition at line 140 of file CDFMidPointPlugin.hh.

Referenced by description().

00140 {return _max_iterations     ;}

int fastjet::CDFMidPointPlugin::max_pair_size  )  const [inline]
 

Definition at line 139 of file CDFMidPointPlugin.hh.

Referenced by description().

00139 {return _max_pair_size      ;}

double fastjet::CDFMidPointPlugin::overlap_threshold  )  const [inline]
 

Definition at line 141 of file CDFMidPointPlugin.hh.

Referenced by description().

00141 {return _overlap_threshold  ;}

void fastjet::CDFMidPointPlugin::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 86 of file CDFMidPointPlugin.cc.

References _cone_area_fraction, _cone_radius, _max_iterations, _max_pair_size, _overlap_threshold, _seed_threshold, _sm_scale, fastjet::ClusterSequence::jets(), fastjet::ClusterSequence::plugin_record_iB_recombination(), and fastjet::ClusterSequence::plugin_record_ij_recombination().

00086                                                                         {
00087  
00088   // create the physics towers needed by the CDF code
00089   vector<PhysicsTower> towers;
00090   towers.reserve(clust_seq.jets().size());
00091   for (unsigned i = 0; i < clust_seq.jets().size(); i++) {
00092     LorentzVector fourvect(clust_seq.jets()[i].px(),
00093                            clust_seq.jets()[i].py(),
00094                            clust_seq.jets()[i].pz(),
00095                            clust_seq.jets()[i].E());
00096     PhysicsTower tower(fourvect);
00097     // misuse one of the indices for tracking, since the MidPoint
00098     // implementation doesn't seem to make use of these indices
00099     tower.calTower.iEta = i;
00100     towers.push_back(tower);
00101   }
00102 
00103   // prepare the CDF algorithm
00104   MidPointAlgorithm m(_seed_threshold,_cone_radius,_cone_area_fraction,
00105                       _max_pair_size,_max_iterations,_overlap_threshold,
00106                       MidPointAlgorithm::SplitMergeScale(_sm_scale));
00107     
00108   // run the CDF algorithm
00109   std::vector<Cluster> jets;
00110   m.run(towers,jets);
00111 
00112 
00113   // now transfer the jets back into our own structure -- we will
00114   // mimic the cone code with a sequential recombination sequence in
00115   // which the jets are built up by adding one particle at a time
00116   for(vector<Cluster>::const_iterator jetIter = jets.begin(); 
00117                                       jetIter != jets.end(); jetIter++) {
00118     const vector<PhysicsTower> & tower_list = jetIter->towerList;
00119     int jet_k = tower_list[0].calTower.iEta;
00120   
00121     int ntow = int(jetIter->towerList.size());
00122     for (int itow = 1; itow < ntow; itow++) {
00123       int jet_i = jet_k;
00124       // retrieve our misappropriated index for the jet
00125       int jet_j = tower_list[itow].calTower.iEta;
00126       // do a fake recombination step with dij=0
00127       double dij = 0.0;
00128       clust_seq.plugin_record_ij_recombination(jet_i, jet_j, dij, jet_k);
00129     }
00130   
00131     // NB: put a sensible looking d_iB just to be nice...
00132     double d_iB = clust_seq.jets()[jet_k].perp2();
00133     clust_seq.plugin_record_iB_recombination(jet_k, d_iB);
00134   }
00135 
00136 
00137   // following code is for testing only
00138   //cout << endl;
00139   //for(vector<Cluster>::const_iterator jetIter = jets.begin(); 
00140   //                                    jetIter != jets.end(); jetIter++) {
00141   //  cout << jetIter->fourVector.pt() << " " << jetIter->fourVector.y() << endl;
00142   //}
00143   //cout << "-----------------------------------------------------\n";
00144   //vector<PseudoJet> ourjets(clust_seq.inclusive_jets());
00145   //for (vector<PseudoJet>::const_reverse_iterator ourjet = ourjets.rbegin();
00146   //     ourjet != ourjets.rend(); ourjet++) {
00147   //  cout << ourjet->perp() << " " << ourjet->rap() << endl;
00148   //}
00149   //cout << endl;
00150 }

double fastjet::CDFMidPointPlugin::seed_threshold  )  const [inline]
 

Definition at line 136 of file CDFMidPointPlugin.hh.

Referenced by description().

00136 {return _seed_threshold     ;}


Member Data Documentation

double fastjet::CDFMidPointPlugin::_cone_area_fraction [private]
 

Definition at line 152 of file CDFMidPointPlugin.hh.

Referenced by run_clustering().

double fastjet::CDFMidPointPlugin::_cone_radius [private]
 

Definition at line 151 of file CDFMidPointPlugin.hh.

Referenced by run_clustering().

int fastjet::CDFMidPointPlugin::_max_iterations [private]
 

Definition at line 154 of file CDFMidPointPlugin.hh.

Referenced by run_clustering().

int fastjet::CDFMidPointPlugin::_max_pair_size [private]
 

Definition at line 153 of file CDFMidPointPlugin.hh.

Referenced by run_clustering().

double fastjet::CDFMidPointPlugin::_overlap_threshold [private]
 

Definition at line 155 of file CDFMidPointPlugin.hh.

Referenced by run_clustering().

double fastjet::CDFMidPointPlugin::_seed_threshold [private]
 

Definition at line 150 of file CDFMidPointPlugin.hh.

Referenced by run_clustering().

SplitMergeScale fastjet::CDFMidPointPlugin::_sm_scale [private]
 

Definition at line 156 of file CDFMidPointPlugin.hh.

Referenced by description(), and run_clustering().


The documentation for this class was generated from the following files:
Generated on Mon Apr 2 20:58:24 2007 for fastjet by  doxygen 1.4.2