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