FastJets.cc
Go to the documentation of this file.
00001 // -*- C++ -*- 00002 #include "Rivet/Config/RivetCommon.hh" 00003 #include "Rivet/Tools/Logging.hh" 00004 #include "Rivet/Projections/FastJets.hh" 00005 #include "Rivet/Projections/HeavyHadrons.hh" 00006 #include "Rivet/Projections/TauFinder.hh" 00007 00008 namespace Rivet { 00009 00010 00011 void FastJets::_init1(JetAlgName alg, double rparameter, double seed_threshold) { 00012 setName("FastJets"); 00013 MSG_DEBUG("JetAlg = " << alg); 00014 MSG_DEBUG("R parameter = " << rparameter); 00015 MSG_DEBUG("Seed threshold = " << seed_threshold); 00016 if (alg == KT) { 00017 _jdef = fastjet::JetDefinition(fastjet::kt_algorithm, rparameter, fastjet::E_scheme); 00018 } else if (alg == CAM) { 00019 _jdef = fastjet::JetDefinition(fastjet::cambridge_algorithm, rparameter, fastjet::E_scheme); 00020 } else if (alg == ANTIKT) { 00021 _jdef = fastjet::JetDefinition(fastjet::antikt_algorithm, rparameter, fastjet::E_scheme); 00022 } else if (alg == DURHAM) { 00023 _jdef = fastjet::JetDefinition(fastjet::ee_kt_algorithm, fastjet::E_scheme); 00024 } else { 00025 // Plugins: 00026 if (alg == SISCONE) { 00027 const double OVERLAP_THRESHOLD = 0.75; 00028 _plugin.reset(new fastjet::SISConePlugin(rparameter, OVERLAP_THRESHOLD)); 00029 } else if (alg == PXCONE) { 00030 string msg = "PxCone currently not supported, since FastJet doesn't install it by default. "; 00031 msg += "Please notify the Rivet authors if this behaviour should be changed."; 00032 throw Error(msg); 00033 //_plugin.reset(new fastjet::PxConePlugin(rparameter)); 00034 } else if (alg == ATLASCONE) { 00035 const double OVERLAP_THRESHOLD = 0.5; 00036 _plugin.reset(new fastjet::ATLASConePlugin(rparameter, seed_threshold, OVERLAP_THRESHOLD)); 00037 } else if (alg == CMSCONE) { 00038 _plugin.reset(new fastjet::CMSIterativeConePlugin(rparameter, seed_threshold)); 00039 } else if (alg == CDFJETCLU) { 00040 const double OVERLAP_THRESHOLD = 0.75; 00041 _plugin.reset(new fastjet::CDFJetCluPlugin(rparameter, OVERLAP_THRESHOLD, seed_threshold)); 00042 } else if (alg == CDFMIDPOINT) { 00043 const double OVERLAP_THRESHOLD = 0.5; 00044 _plugin.reset(new fastjet::CDFMidPointPlugin(rparameter, OVERLAP_THRESHOLD, seed_threshold)); 00045 } else if (alg == D0ILCONE) { 00046 const double min_jet_Et = 6.0; 00047 _plugin.reset(new fastjet::D0RunIIConePlugin(rparameter, min_jet_Et)); 00048 } else if (alg == JADE) { 00049 _plugin.reset(new fastjet::JadePlugin()); 00050 } else if (alg == TRACKJET) { 00051 _plugin.reset(new fastjet::TrackJetPlugin(rparameter)); 00052 } 00053 _jdef = fastjet::JetDefinition(_plugin.get()); 00054 } 00055 addProjection(HeavyHadrons(), "HFHadrons"); 00056 addProjection(TauFinder(TauFinder::HADRONIC), "Taus"); 00057 } 00058 00059 void FastJets::_init2(fastjet::JetAlgorithm type, 00060 fastjet::RecombinationScheme recom, double rparameter) { 00061 setName("FastJets"); 00062 _jdef = fastjet::JetDefinition(type, rparameter, recom); 00063 addProjection(HeavyHadrons(), "HFHadrons"); 00064 addProjection(TauFinder(TauFinder::HADRONIC), "Taus"); 00065 } 00066 00067 void FastJets::_init3(fastjet::JetDefinition::Plugin* plugin) { 00068 setName("FastJets"); 00069 /// @todo Should we be copying the plugin? 00070 _plugin.reset(plugin); 00071 _jdef = fastjet::JetDefinition(_plugin.get()); 00072 addProjection(HeavyHadrons(), "HFHadrons"); 00073 addProjection(TauFinder(TauFinder::HADRONIC), "Taus"); 00074 } 00075 00076 00077 00078 int FastJets::compare(const Projection& p) const { 00079 const FastJets& other = dynamic_cast<const FastJets&>(p); 00080 return \ 00081 cmp(_useInvisibles, other._useInvisibles) || 00082 (_useInvisibles ? mkNamedPCmp(other, "FS") : mkNamedPCmp(other, "VFS")) || 00083 cmp(_jdef.jet_algorithm(), other._jdef.jet_algorithm()) || 00084 cmp(_jdef.recombination_scheme(), other._jdef.recombination_scheme()) || 00085 cmp(_jdef.plugin(), other._jdef.plugin()) || 00086 cmp(_jdef.R(), other._jdef.R()) || 00087 cmp(_adef, other._adef); 00088 } 00089 00090 00091 00092 void FastJets::project(const Event& e) { 00093 // Final state particles 00094 const string fskey = _useInvisibles ? "FS" : "VFS"; 00095 const Particles fsparticles = applyProjection<FinalState>(e, fskey).particles(); 00096 // Tagging particles 00097 const Particles chadrons = applyProjection<HeavyHadrons>(e, "HFHadrons").cHadrons(); 00098 const Particles bhadrons = applyProjection<HeavyHadrons>(e, "HFHadrons").bHadrons(); 00099 const Particles taus = applyProjection<FinalState>(e, "Taus").particles(); 00100 calc(fsparticles, chadrons+bhadrons+taus); 00101 } 00102 00103 00104 void FastJets::calc(const Particles& fsparticles, const Particles& tagparticles) { 00105 _particles.clear(); 00106 vector<fastjet::PseudoJet> pjs; 00107 00108 /// @todo Use FastJet3's UserInfo system 00109 00110 // Store 4 vector data about each particle into FastJet's PseudoJets 00111 int counter = 1; 00112 foreach (const Particle& p, fsparticles) { 00113 const FourMomentum fv = p.momentum(); 00114 fastjet::PseudoJet pj(fv.px(), fv.py(), fv.pz(), fv.E()); ///< @todo Eliminate? 00115 pj.set_user_index(counter); 00116 pjs.push_back(pj); 00117 _particles[counter] = p; 00118 counter += 1; 00119 } 00120 // And the same for ghost tagging particles (with negative user indices) 00121 counter = 1; 00122 foreach (const Particle& p, tagparticles) { 00123 const FourMomentum fv = 1e-20 * p.momentum(); ///< Ghostify the momentum 00124 fastjet::PseudoJet pj(fv.px(), fv.py(), fv.pz(), fv.E()); ///< @todo Eliminate? 00125 pj.set_user_index(-counter); 00126 pjs.push_back(pj); 00127 _particles[-counter] = p; 00128 counter += 1; 00129 } 00130 00131 MSG_DEBUG("Running FastJet ClusterSequence construction"); 00132 // Choose cseq as basic or area-calculating 00133 if (_adef == NULL) { 00134 _cseq.reset(new fastjet::ClusterSequence(pjs, _jdef)); 00135 } else { 00136 _cseq.reset(new fastjet::ClusterSequenceArea(pjs, _jdef, *_adef)); 00137 } 00138 } 00139 00140 00141 void FastJets::reset() { 00142 _yscales.clear(); 00143 _particles.clear(); 00144 /// @todo _cseq = fastjet::ClusterSequence(); 00145 } 00146 00147 00148 size_t FastJets::numJets(double ptmin) const { 00149 if (_cseq.get() == NULL) return 0; 00150 return _cseq->inclusive_jets(ptmin).size(); 00151 } 00152 00153 00154 Jets FastJets::_jets(double ptmin) const { 00155 Jets rtn; 00156 foreach (const fastjet::PseudoJet& pj, pseudojets()) { 00157 assert(clusterSeq() != NULL); 00158 const PseudoJets parts = clusterSeq()->constituents(pj); 00159 vector<Particle> constituents, tags; 00160 constituents.reserve(parts.size()); 00161 foreach (const fastjet::PseudoJet& p, parts) { 00162 map<int, Particle>::const_iterator found = _particles.find(p.user_index()); 00163 assert(found != _particles.end()); 00164 assert(found->first != 0); 00165 if (found->first > 0) constituents.push_back(found->second); 00166 else if (found->first < 0) tags.push_back(found->second); 00167 } 00168 Jet j(pj, constituents, tags); 00169 rtn.push_back(j); 00170 } 00171 /// @todo Cache? 00172 return rtn; 00173 } 00174 00175 00176 PseudoJets FastJets::pseudoJets(double ptmin) const { 00177 return (clusterSeq() != NULL) ? clusterSeq()->inclusive_jets(ptmin) : PseudoJets(); 00178 } 00179 00180 00181 00182 00183 00184 vector<double> FastJets::ySubJet(const fastjet::PseudoJet& jet) const { 00185 assert(clusterSeq()); 00186 fastjet::ClusterSequence subjet_cseq(clusterSeq()->constituents(jet), _jdef); 00187 vector<double> yMergeVals; 00188 for (int i = 1; i < 4; ++i) { 00189 // Multiply the dmerge value by R^2 so that it corresponds to a 00190 // relative k_T (fastjet has 1/R^2 in the d_ij distance by default) 00191 const double ktmerge = subjet_cseq.exclusive_dmerge(i) * _jdef.R()*_jdef.R(); 00192 yMergeVals.push_back(ktmerge/jet.perp2()); 00193 } 00194 _yscales.insert(make_pair( jet.cluster_hist_index(), yMergeVals )); 00195 return yMergeVals; 00196 } 00197 00198 00199 00200 fastjet::PseudoJet FastJets::splitJet(fastjet::PseudoJet jet, double& last_R) const { 00201 // Sanity cuts 00202 if (jet.E() <= 0 || _cseq->constituents(jet).size() <= 1) { 00203 return jet; 00204 } 00205 00206 // Build a new cluster sequence just using the consituents of this jet. 00207 assert(clusterSeq()); 00208 fastjet::ClusterSequence cs(clusterSeq()->constituents(jet), _jdef); 00209 00210 // Get the jet back again 00211 fastjet::PseudoJet remadeJet = cs.inclusive_jets()[0]; 00212 MSG_DEBUG("Jet2:" << remadeJet.m() << "," << remadeJet.e()); 00213 00214 fastjet::PseudoJet parent1, parent2; 00215 fastjet::PseudoJet split(0.0, 0.0, 0.0, 0.0); 00216 while (cs.has_parents(remadeJet, parent1, parent2)) { 00217 MSG_DEBUG("Parents:" << parent1.m() << "," << parent2.m()); 00218 if (parent1.m2() < parent2.m2()) { 00219 fastjet::PseudoJet tmp; 00220 tmp = parent1; parent1 = parent2; parent2 = tmp; 00221 } 00222 00223 double ktdist = parent1.kt_distance(parent2); 00224 double rtycut2 = 0.3*0.3; 00225 if (parent1.m() < ((2.0*remadeJet.m())/3.0) && ktdist > rtycut2*remadeJet.m2()) { 00226 break; 00227 } else { 00228 remadeJet = parent1; 00229 } 00230 } 00231 00232 last_R = 0.5 * sqrt(parent1.squared_distance(parent2)); 00233 split.reset(remadeJet.px(), remadeJet.py(), remadeJet.pz(), remadeJet.E()); 00234 return split; 00235 } 00236 00237 00238 00239 fastjet::PseudoJet FastJets::filterJet(fastjet::PseudoJet jet, 00240 double& stingy_R, const double def_R) const { 00241 assert(clusterSeq()); 00242 00243 if (jet.E() <= 0.0 || clusterSeq()->constituents(jet).size() == 0) { 00244 return jet; 00245 } 00246 if (stingy_R == 0.0) { 00247 stingy_R = def_R; 00248 } 00249 00250 stingy_R = def_R < stingy_R ? def_R : stingy_R; 00251 fastjet::JetDefinition stingy_jet_def(fastjet::cambridge_algorithm, stingy_R); 00252 00253 //FlavourRecombiner recom; 00254 //stingy_jet_def.set_recombiner(&recom); 00255 fastjet::ClusterSequence scs(clusterSeq()->constituents(jet), stingy_jet_def); 00256 std::vector<fastjet::PseudoJet> stingy_jets = sorted_by_pt(scs.inclusive_jets()); 00257 00258 fastjet::PseudoJet reconst_jet(0.0, 0.0, 0.0, 0.0); 00259 00260 for (unsigned isj = 0; isj < std::min(3U, (unsigned int) stingy_jets.size()); ++isj) { 00261 reconst_jet += stingy_jets[isj]; 00262 } 00263 return reconst_jet; 00264 } 00265 00266 00267 } Generated on Tue Mar 24 2015 17:35:26 for The Rivet MC analysis system by ![]() |