CMSIterativeConePlugin is a plugin for fastjet (v2.4 upwards). More...
#include <CMSIterativeConePlugin.hh>
Public Member Functions | |
CMSIterativeConePlugin (double ConeRadius, double SeedThreshold=1.0) | |
Main constructor for the CMSIterativeCone Plugin class. | |
CMSIterativeConePlugin (const CMSIterativeConePlugin &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 | |
virtual double | seed_threshold () const |
get the seed threshold | |
Private Attributes | |
double | theConeRadius |
cone radius | |
double | theSeedThreshold |
seed threshold |
CMSIterativeConePlugin is a plugin for fastjet (v2.4 upwards).
Definition at line 47 of file CMSIterativeConePlugin.hh.
fastjet::CMSIterativeConePlugin::CMSIterativeConePlugin | ( | double | ConeRadius, | |
double | SeedThreshold = 1.0 | |||
) | [inline] |
Main constructor for the CMSIterativeCone Plugin class.
The arguments are ConeRadius the radius of the cone SeedThreshold a threshold for the seeds to iterate from
NOTE: to be more coherent with all other fastjet plugins, we've put the radius before the seed threshold. CMS does the opposite. In this way, we also put a default value of 0 for the seed threshold.
Definition at line 60 of file CMSIterativeConePlugin.hh.
00060 : 00061 theConeRadius(ConeRadius), theSeedThreshold(SeedThreshold){}
fastjet::CMSIterativeConePlugin::CMSIterativeConePlugin | ( | const CMSIterativeConePlugin & | plugin | ) | [inline] |
copy constructor
Definition at line 64 of file CMSIterativeConePlugin.hh.
string fastjet::CMSIterativeConePlugin::description | ( | ) | const [virtual] |
return a textual description of the jet-definition implemented in this plugin
Implements JetDefinition::Plugin.
Definition at line 68 of file CMSIterativeConePlugin.cc.
References theConeRadius, and theSeedThreshold.
00068 { 00069 ostringstream desc; 00070 desc << "CMSIterativeCone plugin with R = " << theConeRadius << " and seed threshold = " << theSeedThreshold; 00071 return desc.str(); 00072 }
virtual double fastjet::CMSIterativeConePlugin::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 74 of file CMSIterativeConePlugin.hh.
00074 {return theConeRadius;}
void fastjet::CMSIterativeConePlugin::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 74 of file CMSIterativeConePlugin.cc.
References fastjet::deltaR2(), PseudoJet::Et(), PseudoJet::eta(), ClusterSequence::jets(), M_PI, PseudoJet::phi(), fastjet::d0::inline_maths::phi(), ClusterSequence::plugin_record_iB_recombination(), ClusterSequence::plugin_record_ij_recombination(), theConeRadius, and theSeedThreshold.
00074 { 00075 00076 // This code is adapted from CMSIterativeConeAlgorithms.cc from the 00077 // CMSSW software. 00078 // The adaptation is just meant to use 00079 // - the FastJet 4-vectors instead of the CMS ones 00080 // - the FastJet clustering history structures instead of the 00081 // ProtoJet one used by CMS. 00082 00083 //make a list of input objects ordered by ET 00084 //cout << "copying the list of particles" << endl; 00085 list<PseudoJet> input; 00086 for (unsigned int i=0 ; i<clust_seq.jets().size() ; i++) { 00087 input.push_back(clust_seq.jets()[i]); 00088 } 00089 NumericSafeGreaterByEt<PseudoJet> compCandidate; 00090 //cout << "sorting" << endl; 00091 input.sort(compCandidate); 00092 00093 //find jets 00094 //cout << "launching the main loop" << endl; 00095 while( !input.empty() && (input.front().Et() > theSeedThreshold )) { 00096 //cone centre 00097 double eta0=input.front().eta(); 00098 double phi0=input.front().phi(); 00099 //protojet properties 00100 double eta=0; 00101 double phi=0; 00102 double et=0; 00103 //list of towers in cone 00104 list< list<PseudoJet>::iterator> cone; 00105 for(int iteration=0;iteration<100;iteration++){ 00106 //cout << "iterating" << endl; 00107 cone.clear(); 00108 eta=0; 00109 phi=0; 00110 et=0; 00111 for(list<PseudoJet>::iterator inp=input.begin(); 00112 inp!=input.end();inp++){ 00113 const PseudoJet tower = *inp; 00114 if( deltaR2(eta0,phi0,tower.eta(),tower.phi()) < 00115 theConeRadius*theConeRadius) { 00116 double tower_et = tower.Et(); 00117 cone.push_back(inp); 00118 eta+= tower_et*tower.eta(); 00119 double dphi=tower.phi()-phi0; 00120 if(dphi>M_PI) dphi-=2*M_PI; 00121 else if(dphi<=-M_PI) dphi+=2*M_PI; 00122 phi+=tower_et*dphi; 00123 et +=tower_et; 00124 } 00125 } 00126 eta=eta/et; 00127 phi=phi0+phi/et; 00128 if(phi>M_PI)phi-=2*M_PI; 00129 else if(phi<=-M_PI)phi+=2*M_PI; 00130 00131 if(fabs(eta-eta0)<.001 && fabs(phi-phi0)<.001) break;//stable cone found 00132 eta0=eta; 00133 phi0=phi; 00134 } 00135 00136 //cout << "make the jet final" << endl; 00137 00138 //make a final jet and remove the jet constituents from the input list 00139 // InputCollection jetConstituents; 00140 // list< list<InputItem>::iterator>::const_iterator inp; 00141 // for(inp=cone.begin();inp!=cone.end();inp++) { 00142 // jetConstituents.push_back(**inp); 00143 // input.erase(*inp); 00144 // } 00145 // fOutput->push_back (ProtoJet (jetConstituents)); 00146 // 00147 // IMPORTANT NOTE: 00148 // while the stability of the stable cone is tested using the Et 00149 // scheme recombination, the final jet uses E-scheme 00150 // recombination. 00151 // 00152 // The technique used here is the same as what we already used for 00153 // SISCone except that we act on the 'cone' list. 00154 // We successively merge the particles that make up the cone jet 00155 // until we have all particles in it. We start off with the zeroth 00156 // particle. 00157 list< list<PseudoJet>::iterator>::const_iterator inp; 00158 inp = cone.begin(); 00159 int jet_k = (*inp)->cluster_hist_index(); 00160 // gps tmp 00161 //float px=(*inp)->px(), py=(*inp)->py(), pz=(*inp)->pz(), E = (*inp)->E(); 00162 00163 // remove the particle from the list and jump to the next one 00164 input.erase(*inp); 00165 inp++; 00166 00167 // now loop over the remaining particles 00168 while (inp != cone.end()){ 00169 // take the last result of the merge 00170 int jet_i = jet_k; 00171 // and the next element of the jet 00172 int jet_j = (*inp)->cluster_hist_index(); 00173 // and merge them (with a fake dij) 00174 double dij = 0.0; 00175 00176 // create the new jet by hand so that we can adjust its user index 00177 // Note again the use of the E-scheme recombination here! 00178 PseudoJet newjet = clust_seq.jets()[jet_i] + clust_seq.jets()[jet_j]; 00179 00180 // gps tmp to try out floating issues 00181 //px+=(*inp)->px(), py+=(*inp)->py(), pz+=(*inp)->pz(), E += (*inp)->E(); 00182 //PseudoJet newjet(px,py,pz,E); 00183 00184 clust_seq.plugin_record_ij_recombination(jet_i, jet_j, dij, newjet, jet_k); 00185 00186 // remove the particle from the list and jump to the next one 00187 input.erase(*inp); 00188 inp++; 00189 } 00190 00191 // we have merged all the jet's particles into a single object, so now 00192 // "declare" it to be a beam (inclusive) jet. 00193 // [NB: put a sensible looking d_iB just to be nice...] 00194 double d_iB = clust_seq.jets()[jet_k].perp2(); 00195 clust_seq.plugin_record_iB_recombination(jet_k, d_iB); 00196 00197 00198 } //loop over seeds ended 00199 00200 00201 }
virtual double fastjet::CMSIterativeConePlugin::seed_threshold | ( | ) | const [inline, virtual] |
get the seed threshold
Definition at line 77 of file CMSIterativeConePlugin.hh.
00077 {return theSeedThreshold;}
double fastjet::CMSIterativeConePlugin::theConeRadius [private] |
cone radius
Definition at line 80 of file CMSIterativeConePlugin.hh.
Referenced by description(), and run_clustering().
double fastjet::CMSIterativeConePlugin::theSeedThreshold [private] |
seed threshold
Definition at line 81 of file CMSIterativeConePlugin.hh.
Referenced by description(), and run_clustering().