Main Page | Namespace List | Class Hierarchy | Class List | Directories | File List | Namespace Members | Class Members | File Members

fastjet_areas.cc

Go to the documentation of this file.
00001 
00002 //STARTHEADER
00003 // $Id: fastjet_areas.cc 432 2007-01-23 12:12:59Z salam $
00004 //
00005 // Copyright (c) 2005-2006, Matteo Cacciari and Gavin Salam
00006 //
00007 //----------------------------------------------------------------------
00008 // This file is part of FastJet.
00009 //
00010 //  FastJet is free software; you can redistribute it and/or modify
00011 //  it under the terms of the GNU General Public License as published by
00012 //  the Free Software Foundation; either version 2 of the License, or
00013 //  (at your option) any later version.
00014 //
00015 //  The algorithms that underlie FastJet have required considerable
00016 //  development and are described in hep-ph/0512210. If you use
00017 //  FastJet as part of work towards a scientific publication, please
00018 //  include a citation to the FastJet paper.
00019 //
00020 //  FastJet is distributed in the hope that it will be useful,
00021 //  but WITHOUT ANY WARRANTY; without even the implied warranty of
00022 //  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
00023 //  GNU General Public License for more details.
00024 //
00025 //  You should have received a copy of the GNU General Public License
00026 //  along with FastJet; if not, write to the Free Software
00027 //  Foundation, Inc.:
00028 //      59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
00029 //----------------------------------------------------------------------
00030 //ENDHEADER
00031 
00032 
00033 
00034 
00035 
00036 
00037 
00038 //----------------------------------------------------------------------
00039 // fastjet areas example program. 
00040 //
00041 // Compile it with: make fastjet_areas
00042 // run it with    : ./fastjet_areas < data/single-event.dat
00043 //
00044 //----------------------------------------------------------------------
00045 #include "fastjet/PseudoJet.hh"
00046 #include "fastjet/ClusterSequenceActiveArea.hh"
00047 #include<iostream> // needed for io
00048 #include<sstream>  // needed for internal io
00049 #include<vector> 
00050 
00051 using namespace std;
00052 
00053 // a declaration of a function that pretty prints a list of jets
00054 void print_jets (const fastjet::ClusterSequenceWithArea &, 
00055                  const vector<fastjet::PseudoJet> &);
00056 
00058 int main (int argc, char ** argv) {
00059   
00060   vector<fastjet::PseudoJet> input_particles;
00061   
00062   // read in input particles
00063   double px, py , pz, E;
00064   while (cin >> px >> py >> pz >> E) {
00065     // create a fastjet::PseudoJet with these components and put it onto
00066     // back of the input_particles vector
00067     input_particles.push_back(fastjet::PseudoJet(px,py,pz,E)); 
00068   }
00069   
00070   // create an object that represents your choice of jet finder and 
00071   // the associated parameters
00072   double Rparam = 1.0;
00073   fastjet::Strategy strategy = fastjet::Best;
00074   fastjet::JetDefinition jet_def(fastjet::kt_algorithm, Rparam, strategy);
00075 
00076   // create an object that specifies how we to define the (active) area
00077   double ghost_etamax = 6.0;
00078   int    active_area_repeats = 3;
00079   double ghost_area    = 0.01;
00080   fastjet::ActiveAreaSpec area_spec(ghost_etamax, active_area_repeats, 
00081                                     ghost_area);
00082 
00083   // run the jet clustering with the above jet definition
00084   fastjet::ClusterSequenceActiveArea clust_seq(input_particles, 
00085                                                jet_def, area_spec);
00086 
00087   // tell the user what was done
00088   cout << "Strategy adopted by FastJet was "<<
00089        clust_seq.strategy_string()<<endl<<endl;
00090 
00091   // extract the inclusive jets with pt > 5 GeV, sorted by pt
00092   double ptmin = 5.0;
00093   vector<fastjet::PseudoJet> inclusive_jets = clust_seq.inclusive_jets(ptmin);
00094 
00095   // print them out
00096   cout << "Printing inclusive jets with pt > "<< ptmin<<" GeV\n";
00097   cout << "---------------------------------------\n";
00098   print_jets(clust_seq, inclusive_jets);
00099   cout << endl;
00100 
00101 
00102 }
00103 
00104 
00105 //----------------------------------------------------------------------
00107 void print_jets (const fastjet::ClusterSequenceWithArea & clust_seq, 
00108                  const vector<fastjet::PseudoJet> & unsorted_jets) {
00109 
00110   // sort jets into increasing pt
00111   vector<fastjet::PseudoJet> jets = sorted_by_pt(unsorted_jets);  
00112 
00113   printf(" ijet   rap      phi        Pt         area  +-   err\n");
00114   for (size_t j = 0; j < jets.size(); j++) {
00115 
00116     double area     = clust_seq.area(jets[j]);
00117     double area_error = clust_seq.area_error(jets[j]);
00118 
00119     printf("%5u %9.5f %8.5f %10.3f %8.3f +- %6.3f\n",j,jets[j].rap(),
00120            jets[j].phi(),jets[j].perp(), area, area_error);
00121   }
00122 
00123 
00124 }
00125 
00126 
00127 
00128 
00129 
00130 
00131 

Generated on Mon Apr 2 20:57:48 2007 for fastjet by  doxygen 1.4.2