FastJet 3.0beta1
|
00001 //STARTHEADER 00002 // $Id: fastjet_subjets.cc 1764 2010-09-18 08:08:15Z salam $ 00003 // 00004 // Copyright (c) 2005-2006, Matteo Cacciari and Gavin Salam 00005 // 00006 //---------------------------------------------------------------------- 00007 // This file is part of FastJet. 00008 // 00009 // FastJet is free software; you can redistribute it and/or modify 00010 // it under the terms of the GNU General Public License as published by 00011 // the Free Software Foundation; either version 2 of the License, or 00012 // (at your option) any later version. 00013 // 00014 // The algorithms that underlie FastJet have required considerable 00015 // development and are described in hep-ph/0512210. If you use 00016 // FastJet as part of work towards a scientific publication, please 00017 // include a citation to the FastJet paper. 00018 // 00019 // FastJet is distributed in the hope that it will be useful, 00020 // but WITHOUT ANY WARRANTY; without even the implied warranty of 00021 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the 00022 // GNU General Public License for more details. 00023 // 00024 // You should have received a copy of the GNU General Public License 00025 // along with FastJet; if not, write to the Free Software 00026 // Foundation, Inc.: 00027 // 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA 00028 //---------------------------------------------------------------------- 00029 //ENDHEADER 00030 00031 00032 //---------------------------------------------------------------------- 00033 // fastjet example program to show how to access subjets; 00034 // 00035 // See also fastjet_higgs_decomp.cc to see the use of subjets for 00036 // identifying boosted higgs (and other objects) 00037 // 00038 // Compile it with: make fastjet_subjets 00039 // run it with : ./fastjet_subjets < data/single-event.dat 00040 // 00041 //---------------------------------------------------------------------- 00042 #include "fastjet/PseudoJet.hh" 00043 #include "fastjet/ClusterSequence.hh" 00044 #include<iostream> // needed for io 00045 #include<sstream> // needed for internal io 00046 #include<vector> 00047 #include <cstdio> 00048 00049 using namespace std; 00050 00051 // to avoid excessive typing, define an abbreviation for the 00052 // fastjet namespace 00053 namespace fj = fastjet; 00054 00055 00056 // a declaration of a function that pretty prints a list of jets 00057 void print_jets (const vector<fj::PseudoJet> &); 00058 00059 // and this pretty prinst a single jet 00060 void print_jet (const fj::PseudoJet & jet); 00061 00062 // pretty print the jets and their subjets 00063 void print_jets_and_sub (const vector<fj::PseudoJet> & jets, 00064 double dcut); 00065 00066 /// an example program showing how to use fastjet 00067 int main (int argc, char ** argv) { 00068 00069 vector<fj::PseudoJet> input_particles; 00070 00071 // read in input particles 00072 double px, py , pz, E; 00073 while (cin >> px >> py >> pz >> E) { 00074 // create a fj::PseudoJet with these components and put it onto 00075 // back of the input_particles vector 00076 input_particles.push_back(fj::PseudoJet(px,py,pz,E)); 00077 } 00078 00079 // We'll do two subjet analyses, physically rather different 00080 double R = 1.0; 00081 //fj::JetDefinition kt_def(fj::kt_algorithm, R); 00082 fj::JetDefinition cam_def(fj::cambridge_algorithm, R); 00083 00084 // run the jet clustering with the above jet definition 00085 //fj::ClusterSequence kt_seq(input_particles, kt_def); 00086 fj::ClusterSequence cam_seq(input_particles, cam_def); 00087 00088 // tell the user what was done 00089 cout << "Ran " << cam_def.description() << endl; 00090 cout << "Strategy adopted by FastJet was "<< 00091 cam_seq.strategy_string()<<endl<<endl; 00092 00093 // extract the inclusive jets with pt > 5 GeV 00094 double ptmin = 5.0; 00095 vector<fj::PseudoJet> inclusive_jets = cam_seq.inclusive_jets(ptmin); 00096 00097 // for the subjets 00098 double smallR = 0.4; 00099 double dcut_cam = pow(smallR/R,2); 00100 00101 // print them out 00102 cout << "Printing inclusive jets (R = "<<R<<") with pt > "<< ptmin<<" GeV\n"; 00103 cout << "and their subjets with smallR = " << smallR << "\n"; 00104 cout << "---------------------------------------\n"; 00105 print_jets_and_sub(inclusive_jets, dcut_cam); 00106 cout << endl; 00107 00108 00109 // print them out 00110 vector<fj::PseudoJet> exclusive_jets = cam_seq.exclusive_jets(dcut_cam); 00111 cout << "Printing exclusive jets with dcut = "<< dcut_cam<<" \n"; 00112 cout << "--------------------------------------------\n"; 00113 print_jets(exclusive_jets); 00114 00115 00116 } 00117 00118 00119 //---------------------------------------------------------------------- 00120 /// a function that pretty prints a list of jets 00121 void print_jets (const vector<fj::PseudoJet> & jets) { 00122 00123 // sort jets into increasing pt 00124 vector<fj::PseudoJet> sorted_jets = sorted_by_pt(jets); 00125 00126 // label the columns 00127 printf("%5s %15s %15s %15s %15s\n","jet #", "rapidity", 00128 "phi", "pt", "n constituents"); 00129 00130 // print out the details for each jet 00131 for (unsigned int i = 0; i < sorted_jets.size(); i++) { 00132 printf("%5u ",i); 00133 print_jet(sorted_jets[i]); 00134 } 00135 00136 } 00137 00138 00139 //---------------------------------------------------------------------- 00140 /// a function that pretty prints a list of jets 00141 void print_jets_and_sub (const vector<fj::PseudoJet> & jets, 00142 double dcut) { 00143 00144 // sort jets into increasing pt 00145 vector<fj::PseudoJet> sorted_jets = sorted_by_pt(jets); 00146 00147 // label the columns 00148 printf("%5s %15s %15s %15s %15s\n","jet #", "rapidity", 00149 "phi", "pt", "n constituents"); 00150 00151 // print out the details for each jet 00152 for (unsigned int i = 0; i < sorted_jets.size(); i++) { 00153 printf("%5u ",i); 00154 print_jet(sorted_jets[i]); 00155 vector<fj::PseudoJet> subjets = sorted_by_pt(sorted_jets[i].exclusive_subjets(dcut)); 00156 for (unsigned int j = 0; j < subjets.size(); j++) { 00157 printf(" -sub-%02u ",j); 00158 print_jet(subjets[j]); 00159 } 00160 } 00161 00162 } 00163 00164 00165 //---------------------------------------------------------------------- 00166 /// print a single jet 00167 void print_jet (const fj::PseudoJet & jet) { 00168 int n_constituents = jet.constituents().size(); 00169 printf("%15.8f %15.8f %15.8f %8u\n", 00170 jet.rap(), jet.phi(), jet.perp(), n_constituents); 00171 }