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.
00064 {SM_pt, SM_Et, SM_mt, SM_pttilde};
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:
. 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.
00112 : 00113 _seed_threshold (seed_threshold ), 00114 _cone_radius (cone_radius ), 00115 _cone_area_fraction (cone_area_fraction ), 00116 _max_pair_size (max_pair_size ), 00117 _max_iterations (max_iterations ), 00118 _overlap_threshold (overlap_threshold ), 00119 _sm_scale (sm_scale) {}
fastjet::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.
00131 : 00132 _seed_threshold (seed_threshold ), 00133 _cone_radius (cone_radius ), 00134 _cone_area_fraction (cone_area_fraction ), 00135 _max_pair_size (2 ), 00136 _max_iterations (100 ), 00137 _overlap_threshold (overlap_threshold ), 00138 _sm_scale (SM_pt) {}
double fastjet::CDFMidPointPlugin::cone_area_fraction | ( | ) | const [inline] |
Definition at line 144 of file CDFMidPointPlugin.hh.
Referenced by description().
00144 {return _cone_area_fraction ;}
double fastjet::CDFMidPointPlugin::cone_radius | ( | ) | const [inline] |
Definition at line 143 of file CDFMidPointPlugin.hh.
Referenced by description().
00143 {return _cone_radius ;}
string fastjet::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.
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.
00046 { 00047 ostringstream desc; 00048 00049 string sm_scale_string = "split-merge uses "; 00050 switch(_sm_scale) { 00051 case SM_pt: 00052 sm_scale_string += "pt"; 00053 break; 00054 case SM_Et: 00055 sm_scale_string += "Et"; 00056 break; 00057 case SM_mt: 00058 sm_scale_string += "mt"; 00059 break; 00060 case SM_pttilde: 00061 sm_scale_string += "pttilde (scalar sum of pts)"; 00062 break; 00063 default: 00064 ostringstream err; 00065 err << "Unrecognized split-merge scale choice = " << _sm_scale; 00066 throw Error(err.str()); 00067 } 00068 00069 00070 if (cone_area_fraction() == 1) { 00071 desc << "CDF MidPoint jet algorithm, with " ; 00072 } else { 00073 desc << "CDF MidPoint+Searchcone jet algorithm, with "; 00074 } 00075 desc << "seed_threshold = " << seed_threshold () << ", " 00076 << "cone_radius = " << cone_radius () << ", " 00077 << "cone_area_fraction = " << cone_area_fraction () << ", " 00078 << "max_pair_size = " << max_pair_size () << ", " 00079 << "max_iterations = " << max_iterations () << ", " 00080 << "overlap_threshold = " << overlap_threshold () << ", " 00081 << sm_scale_string ; 00082 00083 return desc.str(); 00084 }
int fastjet::CDFMidPointPlugin::max_iterations | ( | ) | const [inline] |
Definition at line 146 of file CDFMidPointPlugin.hh.
Referenced by description().
00146 {return _max_iterations ;}
int fastjet::CDFMidPointPlugin::max_pair_size | ( | ) | const [inline] |
Definition at line 145 of file CDFMidPointPlugin.hh.
Referenced by description().
00145 {return _max_pair_size ;}
double fastjet::CDFMidPointPlugin::overlap_threshold | ( | ) | const [inline] |
Definition at line 147 of file CDFMidPointPlugin.hh.
Referenced by description().
00147 {return _overlap_threshold ;}
virtual double fastjet::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.
00154 {return cone_radius();}
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:
..)
Implements JetDefinition::Plugin.
Definition at line 87 of file CDFMidPointPlugin.cc.
References _cone_area_fraction, _cone_radius, _max_iterations, _max_pair_size, _overlap_threshold, _seed_threshold, _sm_scale, ClusterSequence::jets(), ClusterSequence::plugin_record_iB_recombination(), and ClusterSequence::plugin_record_ij_recombination().
00087 { 00088 00089 // create the physics towers needed by the CDF code 00090 vector<PhysicsTower> towers; 00091 towers.reserve(clust_seq.jets().size()); 00092 for (unsigned i = 0; i < clust_seq.jets().size(); i++) { 00093 LorentzVector fourvect(clust_seq.jets()[i].px(), 00094 clust_seq.jets()[i].py(), 00095 clust_seq.jets()[i].pz(), 00096 clust_seq.jets()[i].E()); 00097 PhysicsTower tower(fourvect); 00098 // misuse one of the indices for tracking, since the MidPoint 00099 // implementation doesn't seem to make use of these indices 00100 tower.calTower.iEta = i; 00101 towers.push_back(tower); 00102 } 00103 00104 // prepare the CDF algorithm 00105 MidPointAlgorithm m(_seed_threshold,_cone_radius,_cone_area_fraction, 00106 _max_pair_size,_max_iterations,_overlap_threshold, 00107 MidPointAlgorithm::SplitMergeScale(_sm_scale)); 00108 00109 // run the CDF algorithm 00110 std::vector<Cluster> jets; 00111 m.run(towers,jets); 00112 00113 00114 // now transfer the jets back into our own structure -- we will 00115 // mimic the cone code with a sequential recombination sequence in 00116 // which the jets are built up by adding one particle at a time 00117 for(vector<Cluster>::const_iterator jetIter = jets.begin(); 00118 jetIter != jets.end(); jetIter++) { 00119 const vector<PhysicsTower> & tower_list = jetIter->towerList; 00120 int jet_k = tower_list[0].calTower.iEta; 00121 00122 int ntow = int(jetIter->towerList.size()); 00123 for (int itow = 1; itow < ntow; itow++) { 00124 int jet_i = jet_k; 00125 // retrieve our misappropriated index for the jet 00126 int jet_j = tower_list[itow].calTower.iEta; 00127 // do a fake recombination step with dij=0 00128 double dij = 0.0; 00129 clust_seq.plugin_record_ij_recombination(jet_i, jet_j, dij, jet_k); 00130 } 00131 00132 // NB: put a sensible looking d_iB just to be nice... 00133 double d_iB = clust_seq.jets()[jet_k].perp2(); 00134 clust_seq.plugin_record_iB_recombination(jet_k, d_iB); 00135 } 00136 00137 00138 // following code is for testing only 00139 //cout << endl; 00140 //for(vector<Cluster>::const_iterator jetIter = jets.begin(); 00141 // jetIter != jets.end(); jetIter++) { 00142 // cout << jetIter->fourVector.pt() << " " << jetIter->fourVector.y() << endl; 00143 //} 00144 //cout << "-----------------------------------------------------\n"; 00145 //vector<PseudoJet> ourjets(clust_seq.inclusive_jets()); 00146 //for (vector<PseudoJet>::const_reverse_iterator ourjet = ourjets.rbegin(); 00147 // ourjet != ourjets.rend(); ourjet++) { 00148 // cout << ourjet->perp() << " " << ourjet->rap() << endl; 00149 //} 00150 //cout << endl; 00151 }
double fastjet::CDFMidPointPlugin::seed_threshold | ( | ) | const [inline] |
Definition at line 142 of file CDFMidPointPlugin.hh.
Referenced by description().
00142 {return _seed_threshold ;}
double fastjet::CDFMidPointPlugin::_cone_area_fraction [private] |
Definition at line 160 of file CDFMidPointPlugin.hh.
Referenced by run_clustering().
double fastjet::CDFMidPointPlugin::_cone_radius [private] |
Definition at line 159 of file CDFMidPointPlugin.hh.
Referenced by run_clustering().
int fastjet::CDFMidPointPlugin::_max_iterations [private] |
Definition at line 162 of file CDFMidPointPlugin.hh.
Referenced by run_clustering().
int fastjet::CDFMidPointPlugin::_max_pair_size [private] |
Definition at line 161 of file CDFMidPointPlugin.hh.
Referenced by run_clustering().
double fastjet::CDFMidPointPlugin::_overlap_threshold [private] |
Definition at line 163 of file CDFMidPointPlugin.hh.
Referenced by run_clustering().
double fastjet::CDFMidPointPlugin::_seed_threshold [private] |
Definition at line 158 of file CDFMidPointPlugin.hh.
Referenced by run_clustering().
Definition at line 164 of file CDFMidPointPlugin.hh.
Referenced by description(), and run_clustering().