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