FastJet 3.0beta1
01-basic.cc
Go to the documentation of this file.
00001 //----------------------------------------------------------------------
00002 /// \file
00003 /// \page Example01 01 - basic usage example
00004 ///
00005 /// fastjet basic example program:
00006 ///   simplest illustration of the usage of the basic classes:
00007 ///   fastjet::PseudoJet, fastjet::JetDefinition and 
00008 ///   fastjet::ClusterSequence
00009 ///
00010 /// run it with    : ./01-basic < data/single-event.dat
00011 ///
00012 /// Source code: 01-basic.cc
00013 //----------------------------------------------------------------------
00014 
00015 //STARTHEADER
00016 // $Id: 01-basic.cc 2003 2011-03-11 16:10:08Z soyez $
00017 //
00018 // Copyright (c) 2005-2011, Matteo Cacciari, Gavin Salam and Gregory Soyez
00019 //
00020 //----------------------------------------------------------------------
00021 // This file is part of FastJet.
00022 //
00023 //  FastJet is free software; you can redistribute it and/or modify
00024 //  it under the terms of the GNU General Public License as published by
00025 //  the Free Software Foundation; either version 2 of the License, or
00026 //  (at your option) any later version.
00027 //
00028 //  The algorithms that underlie FastJet have required considerable
00029 //  development and are described in hep-ph/0512210. If you use
00030 //  FastJet as part of work towards a scientific publication, please
00031 //  include a citation to the FastJet paper.
00032 //
00033 //  FastJet is distributed in the hope that it will be useful,
00034 //  but WITHOUT ANY WARRANTY; without even the implied warranty of
00035 //  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
00036 //  GNU General Public License for more details.
00037 //
00038 //  You should have received a copy of the GNU General Public License
00039 //  along with FastJet; if not, write to the Free Software
00040 //  Foundation, Inc.:
00041 //      59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
00042 //----------------------------------------------------------------------
00043 //ENDHEADER
00044 
00045 #include "fastjet/ClusterSequence.hh"
00046 #include <iostream> // needed for io
00047 #include <cstdio>   // needed for io
00048 
00049 using namespace std;
00050 
00051 /// an example program showing how to use fastjet
00052 int main (int argc, char ** argv) {
00053   
00054   // read in input particles
00055   //----------------------------------------------------------
00056   vector<fastjet::PseudoJet> input_particles;
00057   
00058   double px, py , pz, E;
00059   while (cin >> px >> py >> pz >> E) {
00060     // create a fastjet::PseudoJet with these components and put it onto
00061     // back of the input_particles vector
00062     input_particles.push_back(fastjet::PseudoJet(px,py,pz,E)); 
00063   }
00064   
00065 
00066   // create a jet definition: 
00067   // a jet algorithm with a given radius parameter
00068   //----------------------------------------------------------
00069   double R = 0.6;
00070   fastjet::JetDefinition jet_def(fastjet::antikt_algorithm, R);
00071 
00072 
00073   // run the jet clustering with the above jet definition
00074   //----------------------------------------------------------
00075   fastjet::ClusterSequence clust_seq(input_particles, jet_def);
00076 
00077 
00078   // get the resulting jets ordered in pt
00079   //----------------------------------------------------------
00080   double ptmin = 5.0;
00081   vector<fastjet::PseudoJet> inclusive_jets = sorted_by_pt(clust_seq.inclusive_jets(ptmin));
00082 
00083 
00084   // tell the user what was done
00085   //  - the description of the algorithm used
00086   //  - extract the inclusive jets with pt > 5 GeV
00087   //    show the output as 
00088   //      {index, rap, phi, pt}
00089   //----------------------------------------------------------
00090   cout << "Ran " << jet_def.description() << endl;
00091 
00092   // label the columns
00093   printf("%5s %15s %15s %15s\n","jet #", "rapidity", "phi", "pt");
00094  
00095   // print out the details for each jet
00096   for (unsigned int i = 0; i < inclusive_jets.size(); i++) {
00097     printf("%5u %15.8f %15.8f %15.8f\n",
00098            i, inclusive_jets[i].rap(), inclusive_jets[i].phi(),
00099            inclusive_jets[i].perp());
00100   }
00101 
00102   return 0;
00103 }
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends