fastjet 2.4.3
|
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>
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, 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:
| |
virtual double | R () const |
the plugin mechanism's standard way of accessing the jet radius | |
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 |
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.
the choice of scale to be used in the split-merge step
Definition at line 64 of file CDFMidPointPlugin.hh.
{SM_pt, SM_Et, SM_mt, SM_pttilde};
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:
. 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.
: _seed_threshold (seed_threshold ), _cone_radius (cone_radius ), _cone_area_fraction (cone_area_fraction ), _max_pair_size (max_pair_size ), _max_iterations (max_iterations ), _overlap_threshold (overlap_threshold ), _sm_scale (sm_scale) {}
CDFMidPointPlugin::CDFMidPointPlugin | ( | double | cone_radius, |
double | overlap_threshold, | ||
double | seed_threshold = 1.0 , |
||
double | cone_area_fraction = 1.0 |
||
) | [inline] |
a compact constructor
NB: as of version 2.4, the default value for the overlap_threshold threshold has been removed, to avoid misleading people into using the value of 0.5 without thinking, which is known to have adverse effects in high-noise environments. A recommended value is 0.75.
Definition at line 128 of file CDFMidPointPlugin.hh.
: _seed_threshold (seed_threshold ), _cone_radius (cone_radius ), _cone_area_fraction (cone_area_fraction ), _max_pair_size (2 ), _max_iterations (100 ), _overlap_threshold (overlap_threshold ), _sm_scale (SM_pt) {}
double CDFMidPointPlugin::cone_area_fraction | ( | ) | const [inline] |
Definition at line 144 of file CDFMidPointPlugin.hh.
References _cone_area_fraction.
{return _cone_area_fraction ;}
double CDFMidPointPlugin::cone_radius | ( | ) | const [inline] |
Definition at line 143 of file CDFMidPointPlugin.hh.
References _cone_radius.
Referenced by R().
{return _cone_radius ;}
string CDFMidPointPlugin::description | ( | ) | const [virtual] |
return a textual description of the jet-definition implemented in this plugin
Implements JetDefinition::Plugin.
Definition at line 46 of file CDFMidPointPlugin.cc.
{ ostringstream desc; string sm_scale_string = "split-merge uses "; switch(_sm_scale) { case SM_pt: sm_scale_string += "pt"; break; case SM_Et: sm_scale_string += "Et"; break; case SM_mt: sm_scale_string += "mt"; break; case SM_pttilde: sm_scale_string += "pttilde (scalar sum of pts)"; break; default: ostringstream err; err << "Unrecognized split-merge scale choice = " << _sm_scale; throw Error(err.str()); } if (cone_area_fraction() == 1) { desc << "CDF MidPoint jet algorithm, with " ; } else { desc << "CDF MidPoint+Searchcone jet algorithm, with "; } desc << "seed_threshold = " << seed_threshold () << ", " << "cone_radius = " << cone_radius () << ", " << "cone_area_fraction = " << cone_area_fraction () << ", " << "max_pair_size = " << max_pair_size () << ", " << "max_iterations = " << max_iterations () << ", " << "overlap_threshold = " << overlap_threshold () << ", " << sm_scale_string ; return desc.str(); }
int CDFMidPointPlugin::max_iterations | ( | ) | const [inline] |
Definition at line 146 of file CDFMidPointPlugin.hh.
References _max_iterations.
{return _max_iterations ;}
int CDFMidPointPlugin::max_pair_size | ( | ) | const [inline] |
Definition at line 145 of file CDFMidPointPlugin.hh.
References _max_pair_size.
{return _max_pair_size ;}
double CDFMidPointPlugin::overlap_threshold | ( | ) | const [inline] |
Definition at line 147 of file CDFMidPointPlugin.hh.
References _overlap_threshold.
{return _overlap_threshold ;}
virtual double CDFMidPointPlugin::R | ( | ) | const [inline, virtual] |
the plugin mechanism's standard way of accessing the jet radius
Implements JetDefinition::Plugin.
Definition at line 154 of file CDFMidPointPlugin.hh.
References cone_radius().
{return cone_radius();}
void 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:
Implements JetDefinition::Plugin.
Definition at line 87 of file CDFMidPointPlugin.cc.
References ClusterSequence::jets(), ClusterSequence::plugin_record_iB_recombination(), and ClusterSequence::plugin_record_ij_recombination().
{ // create the physics towers needed by the CDF code vector<PhysicsTower> towers; towers.reserve(clust_seq.jets().size()); for (unsigned i = 0; i < clust_seq.jets().size(); i++) { LorentzVector fourvect(clust_seq.jets()[i].px(), clust_seq.jets()[i].py(), clust_seq.jets()[i].pz(), clust_seq.jets()[i].E()); PhysicsTower tower(fourvect); // misuse one of the indices for tracking, since the MidPoint // implementation doesn't seem to make use of these indices tower.calTower.iEta = i; towers.push_back(tower); } // prepare the CDF algorithm MidPointAlgorithm m(_seed_threshold,_cone_radius,_cone_area_fraction, _max_pair_size,_max_iterations,_overlap_threshold, MidPointAlgorithm::SplitMergeScale(_sm_scale)); // run the CDF algorithm std::vector<Cluster> jets; m.run(towers,jets); // now transfer the jets back into our own structure -- we will // mimic the cone code with a sequential recombination sequence in // which the jets are built up by adding one particle at a time for(vector<Cluster>::const_iterator jetIter = jets.begin(); jetIter != jets.end(); jetIter++) { const vector<PhysicsTower> & tower_list = jetIter->towerList; int jet_k = tower_list[0].calTower.iEta; int ntow = int(jetIter->towerList.size()); for (int itow = 1; itow < ntow; itow++) { int jet_i = jet_k; // retrieve our misappropriated index for the jet int jet_j = tower_list[itow].calTower.iEta; // do a fake recombination step with dij=0 double dij = 0.0; clust_seq.plugin_record_ij_recombination(jet_i, jet_j, dij, jet_k); } // 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); } // following code is for testing only //cout << endl; //for(vector<Cluster>::const_iterator jetIter = jets.begin(); // jetIter != jets.end(); jetIter++) { // cout << jetIter->fourVector.pt() << " " << jetIter->fourVector.y() << endl; //} //cout << "-----------------------------------------------------\n"; //vector<PseudoJet> ourjets(clust_seq.inclusive_jets()); //for (vector<PseudoJet>::const_reverse_iterator ourjet = ourjets.rbegin(); // ourjet != ourjets.rend(); ourjet++) { // cout << ourjet->perp() << " " << ourjet->rap() << endl; //} //cout << endl; }
double CDFMidPointPlugin::seed_threshold | ( | ) | const [inline] |
Definition at line 142 of file CDFMidPointPlugin.hh.
References _seed_threshold.
{return _seed_threshold ;}
double CDFMidPointPlugin::_cone_area_fraction [private] |
Definition at line 160 of file CDFMidPointPlugin.hh.
Referenced by cone_area_fraction().
double CDFMidPointPlugin::_cone_radius [private] |
Definition at line 159 of file CDFMidPointPlugin.hh.
Referenced by cone_radius().
int CDFMidPointPlugin::_max_iterations [private] |
Definition at line 162 of file CDFMidPointPlugin.hh.
Referenced by max_iterations().
int CDFMidPointPlugin::_max_pair_size [private] |
Definition at line 161 of file CDFMidPointPlugin.hh.
Referenced by max_pair_size().
double CDFMidPointPlugin::_overlap_threshold [private] |
Definition at line 163 of file CDFMidPointPlugin.hh.
Referenced by overlap_threshold().
double CDFMidPointPlugin::_seed_threshold [private] |
Definition at line 158 of file CDFMidPointPlugin.hh.
Referenced by seed_threshold().
SplitMergeScale CDFMidPointPlugin::_sm_scale [private] |
Definition at line 164 of file CDFMidPointPlugin.hh.