fastjet_example.cc

Go to the documentation of this file.
00001 
00002 
00003 
00004 
00005 
00006 
00007 
00008 
00009 
00010 
00011 
00012 
00013 
00014 
00015 
00016 
00017 
00018 
00019 //STARTHEADER
00020 // $Id: fastjet_example.cc 1881 2011-01-27 14:25:24Z soyez $
00021 //
00022 // Copyright (c) 2005-2006, Matteo Cacciari and Gavin Salam
00023 //
00024 //----------------------------------------------------------------------
00025 // This file is part of FastJet.
00026 //
00027 //  FastJet is free software; you can redistribute it and/or modify
00028 //  it under the terms of the GNU General Public License as published by
00029 //  the Free Software Foundation; either version 2 of the License, or
00030 //  (at your option) any later version.
00031 //
00032 //  The algorithms that underlie FastJet have required considerable
00033 //  development and are described in hep-ph/0512210. If you use
00034 //  FastJet as part of work towards a scientific publication, please
00035 //  include a citation to the FastJet paper.
00036 //
00037 //  FastJet is distributed in the hope that it will be useful,
00038 //  but WITHOUT ANY WARRANTY; without even the implied warranty of
00039 //  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
00040 //  GNU General Public License for more details.
00041 //
00042 //  You should have received a copy of the GNU General Public License
00043 //  along with FastJet; if not, write to the Free Software
00044 //  Foundation, Inc.:
00045 //      59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
00046 //----------------------------------------------------------------------
00047 //ENDHEADER
00048 
00049 
00050 //----------------------------------------------------------------------
00051 // fastjet example program. 
00052 //
00053 // Compile it with: make fastjet_example
00054 // run it with    : ./fastjet_example < data/single-event.dat
00055 //
00056 // People who are familiar with the ktjet package are encouraged to
00057 // compare this file to the ktjet_example.cc program which does the
00058 // same thing in the ktjet framework.
00059 //----------------------------------------------------------------------
00060 #include "fastjet/PseudoJet.hh"
00061 #include "fastjet/ClusterSequence.hh"
00062 #include<iostream> // needed for io
00063 #include<sstream>  // needed for internal io
00064 #include<vector> 
00065 #include <cstdio>
00066 
00067 using namespace std;
00068 
00069 // a declaration of a function that pretty prints a list of jets
00070 void print_jets (const vector<fastjet::PseudoJet> &);
00071 
00073 int main (int argc, char ** argv) {
00074   
00075   vector<fastjet::PseudoJet> input_particles;
00076   
00077   // Read in input particles
00078   double px, py , pz, E;
00079   while (cin >> px >> py >> pz >> E) {
00080     // create a fastjet::PseudoJet with these components and put it onto
00081     // back of the input_particles vector
00082     input_particles.push_back(fastjet::PseudoJet(px,py,pz,E)); 
00083   }
00084   
00085   // create an object that represents your choice of jet algorithm and 
00086   // the associated parameters
00087   double Rparam = 1.0;
00088   fastjet::Strategy strategy = fastjet::Best;
00089   fastjet::RecombinationScheme recomb_scheme = fastjet::E_scheme;
00090   fastjet::JetDefinition jet_def(fastjet::kt_algorithm, Rparam, recomb_scheme, strategy);
00091 
00092   // run the jet clustering with the above jet definition
00093   fastjet::ClusterSequence clust_seq(input_particles, jet_def);
00094 
00095   // tell the user what was done
00096   cout << "Ran " << jet_def.description() << endl;
00097   cout << "Strategy adopted by FastJet was "<<
00098        clust_seq.strategy_string()<<endl<<endl;
00099 
00100   // extract the inclusive jets with pt > 5 GeV
00101   double ptmin = 5.0;
00102   vector<fastjet::PseudoJet> inclusive_jets = clust_seq.inclusive_jets(ptmin);
00103 
00104   // print them out
00105   cout << "Printing inclusive jets with pt > "<< ptmin<<" GeV\n";
00106   cout << "---------------------------------------\n";
00107   print_jets(inclusive_jets);
00108   cout << endl;
00109 
00110   // extract the exclusive jets with dcut = 25 GeV^2 
00111   double dcut = 25.0;
00112   vector<fastjet::PseudoJet> exclusive_jets = clust_seq.exclusive_jets(dcut);
00113 
00114   // print them out
00115   cout << "Printing exclusive jets with dcut = "<< dcut<<" GeV^2\n";
00116   cout << "--------------------------------------------\n";
00117   print_jets(exclusive_jets);
00118 
00119 
00120 }
00121 
00122 
00123 //----------------------------------------------------------------------
00125 void print_jets (const vector<fastjet::PseudoJet> & jets) {
00126 
00127   // sort jets into increasing pt
00128   vector<fastjet::PseudoJet> sorted_jets = sorted_by_pt(jets);  
00129 
00130   // label the columns
00131   printf("%5s %15s %15s %15s %15s\n","jet #", "rapidity", 
00132          "phi", "pt", "n constituents");
00133   
00134   // print out the details for each jet
00135   for (unsigned int i = 0; i < sorted_jets.size(); i++) {
00136     // the following is not super efficient since it creates an
00137     // intermediate constituents vector
00138     int n_constituents = sorted_jets[i].constituents().size();
00139     printf("%5u %15.8f %15.8f %15.8f %8u\n",
00140            i, sorted_jets[i].rap(), sorted_jets[i].phi(),
00141            sorted_jets[i].perp(), n_constituents);
00142   }
00143 
00144 }