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