FastJet 3.0beta1
|
00001 //---------------------------------------------------------------------- 00002 /// \file 00003 /// \page Example13 13 - boosted top tagging 00004 /// 00005 /// fastjet example program, illustration of carrying out boosted 00006 /// top subjet ID analysis using the Johns Hopkins top tagger as 00007 /// introduced in arXiv:0806.0848 (Kaplan, Rehermann, Schwartz 00008 /// and Tweedie) 00009 /// 00010 /// run it with : ./13-boosted_top < data/boosted_top_event.dat 00011 /// 00012 /// Source code: 13-boosted_top.cc 00013 //---------------------------------------------------------------------- 00014 00015 00016 //STARTHEADER 00017 // $Id: 13-boosted_top.cc 2493 2011-08-03 15:48:16Z salam $ 00018 // 00019 // Copyright (c) 2005-2011, Matteo Cacciari, Gavin Salam and Gregory Soyez 00020 // 00021 //---------------------------------------------------------------------- 00022 // This file is part of FastJet. 00023 // 00024 // FastJet is free software; you can redistribute it and/or modify 00025 // it under the terms of the GNU General Public License as published by 00026 // the Free Software Foundation; either version 2 of the License, or 00027 // (at your option) any later version. 00028 // 00029 // The algorithms that underlie FastJet have required considerable 00030 // development and are described in hep-ph/0512210. If you use 00031 // FastJet as part of work towards a scientific publication, please 00032 // include a citation to the FastJet paper. 00033 // 00034 // FastJet is distributed in the hope that it will be useful, 00035 // but WITHOUT ANY WARRANTY; without even the implied warranty of 00036 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the 00037 // GNU General Public License for more details. 00038 // 00039 // You should have received a copy of the GNU General Public License 00040 // along with FastJet; if not, write to the Free Software 00041 // Foundation, Inc.: 00042 // 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA 00043 //---------------------------------------------------------------------- 00044 //ENDHEADER 00045 00046 #include <iostream> // needed for io 00047 #include <sstream> // needed for internal io 00048 #include <iomanip> 00049 #include <cmath> 00050 00051 #include <fastjet/ClusterSequence.hh> 00052 #include <fastjet/Selector.hh> 00053 #include <fastjet/tools/JHTopTagger.hh> 00054 00055 using namespace std; 00056 using namespace fastjet; 00057 00058 00059 //---------------------------------------------------------------------- 00060 // forward declaration for printing out info about a jet 00061 //---------------------------------------------------------------------- 00062 ostream & operator<<(ostream &, const PseudoJet &); 00063 00064 //---------------------------------------------------------------------- 00065 // core of the program 00066 //---------------------------------------------------------------------- 00067 int main (int argc, char ** argv) { 00068 00069 vector<PseudoJet> particles; 00070 00071 // read in data in format px py pz E b-tag [last of these is optional] 00072 // lines starting with "#" are considered as comments and discarded 00073 //---------------------------------------------------------- 00074 string line; 00075 while (getline(cin,line)) { 00076 if (line.substr(0,1) == "#") {continue;} 00077 istringstream linestream(line); 00078 double px,py,pz,E; 00079 linestream >> px >> py >> pz >> E; 00080 00081 // construct the particle 00082 particles.push_back(PseudoJet(px,py,pz,E)); 00083 } 00084 00085 // compute the parameters to be used through the analysis 00086 // ---------------------------------------------------------- 00087 double Et=0; 00088 for (unsigned int i=0; i<particles.size(); i++) 00089 Et += particles[i].perp(); 00090 00091 double R, delta_p, delta_r; 00092 if (Et>2600){ R=0.4; delta_p=0.05; delta_r=0.19;} 00093 else if (Et>1600){ R=0.6; delta_p=0.05; delta_r=0.19;} 00094 else if (Et>1000){ R=0.8; delta_p=0.10; delta_r=0.19;} 00095 else{ cerr << "Et has to be at least 1 TeV"<< endl; return 1;} 00096 00097 double ptmin = min(500.0, 0.7*Et/2); 00098 00099 // find the jets 00100 // ---------------------------------------------------------- 00101 JetDefinition jet_def(cambridge_algorithm, R); 00102 ClusterSequence cs(particles, jet_def); 00103 vector<PseudoJet> jets = sorted_by_pt(cs.inclusive_jets()); 00104 00105 cout << "Ran: " << jet_def.description() << endl << endl; 00106 cout << "2 Hardest jets: " << jets[0] << endl 00107 << " " << jets[1] << endl << endl; 00108 00109 if (jets[0].perp()<ptmin){ 00110 cout << "No jet above the ptmin threshold" << endl; 00111 return 2; 00112 } 00113 00114 // now do jet tagging using the Johns Hopkins top tagger 00115 // For simplicity, we just apply it to the hardest jet. 00116 // 00117 // In addition to delta_p and delta_r, note that there are two 00118 // further parameters to the JH top tagger that here are implicitly 00119 // set to their defaults: 00120 // 00121 // - cos_theta_W_max (defaults to 0.7) 00122 // - mW (defaults to 80.4). 00123 // 00124 // The value for mW implicitly assumes that momenta are passed in 00125 // GeV. 00126 // ---------------------------------------------------------- 00127 JHTopTagger top_tagger(delta_p, delta_r); 00128 top_tagger.set_top_selector(SelectorMassRange(150,200)); 00129 top_tagger.set_W_selector (SelectorMassRange( 65, 95)); 00130 00131 PseudoJet tagged = top_tagger(jets[0]); 00132 00133 cout << "Ran the following top tagger: " << top_tagger.description() << endl; 00134 00135 if (tagged == 0){ 00136 cout << "No top substructure found" << endl; 00137 return 0; 00138 } 00139 00140 cout << "Found top substructure from the hardest jet:" << endl; 00141 cout << " top candidate: " << tagged << endl; 00142 cout << " |_ W candidate: " << tagged.structure_of<JHTopTagger>().W() << endl; 00143 cout << " | |_ W subjet 1: " << tagged.structure_of<JHTopTagger>().W1() << endl; 00144 cout << " | |_ W subjet 2: " << tagged.structure_of<JHTopTagger>().W2() << endl; 00145 cout << " | cos(theta_W) = " << tagged.structure_of<JHTopTagger>().cos_theta_W() << endl; 00146 cout << " |_ non-W subjet: " << tagged.structure_of<JHTopTagger>().non_W() << endl; 00147 } 00148 00149 00150 //---------------------------------------------------------------------- 00151 // does the actual work for printing out a jet 00152 //---------------------------------------------------------------------- 00153 ostream & operator<<(ostream & ostr, const PseudoJet & jet) { 00154 ostr << "pt, y, phi =" << setprecision(6) 00155 << " " << setw(9) << jet.perp() 00156 << " " << setw(9) << jet.rap() 00157 << " " << setw(9) << jet.phi() 00158 << ", mass = " << setw(9) << jet.m(); 00159 return ostr; 00160 }