FastJet 3.0beta1
|
00001 //---------------------------------------------------------------------- 00002 /// \file 00003 /// \page Example06 06 - using jet areas 00004 /// 00005 /// fastjet example program for jet areas 00006 /// It mostly illustrates the usage of the 00007 /// fastjet::AreaDefinition and fastjet::ClusterSequenceArea classes 00008 /// 00009 /// run it with : ./06-area < data/single-event.dat 00010 /// 00011 /// Source code: 06-area.cc 00012 //---------------------------------------------------------------------- 00013 00014 //STARTHEADER 00015 // $Id: 06-area.cc 2003 2011-03-11 16:10:08Z soyez $ 00016 // 00017 // Copyright (c) 2005-2011, Matteo Cacciari, Gavin Salam and Gregory Soyez 00018 // 00019 //---------------------------------------------------------------------- 00020 // This file is part of FastJet. 00021 // 00022 // FastJet is free software; you can redistribute it and/or modify 00023 // it under the terms of the GNU General Public License as published by 00024 // the Free Software Foundation; either version 2 of the License, or 00025 // (at your option) any later version. 00026 // 00027 // The algorithms that underlie FastJet have required considerable 00028 // development and are described in hep-ph/0512210. If you use 00029 // FastJet as part of work towards a scientific publication, please 00030 // include a citation to the FastJet paper. 00031 // 00032 // FastJet is distributed in the hope that it will be useful, 00033 // but WITHOUT ANY WARRANTY; without even the implied warranty of 00034 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the 00035 // GNU General Public License for more details. 00036 // 00037 // You should have received a copy of the GNU General Public License 00038 // along with FastJet; if not, write to the Free Software 00039 // Foundation, Inc.: 00040 // 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA 00041 //---------------------------------------------------------------------- 00042 //ENDHEADER 00043 00044 #include "fastjet/ClusterSequenceArea.hh" // use this instead of the "usual" ClusterSequence to get area support 00045 #include <iostream> // needed for io 00046 #include <cstdio> // needed for io 00047 00048 using namespace std; 00049 00050 /// an example program showing how to use fastjet 00051 int main (int argc, char ** argv) { 00052 00053 // read in input particles 00054 //---------------------------------------------------------- 00055 vector<fastjet::PseudoJet> input_particles; 00056 00057 double px, py , pz, E; 00058 while (cin >> px >> py >> pz >> E) { 00059 // create a fastjet::PseudoJet with these components and put it onto 00060 // back of the input_particles vector 00061 input_particles.push_back(fastjet::PseudoJet(px,py,pz,E)); 00062 } 00063 00064 00065 // create a jet definition: 00066 // a jet algorithm with a given radius parameter 00067 //---------------------------------------------------------- 00068 double R = 0.6; 00069 fastjet::JetDefinition jet_def(fastjet::kt_algorithm, R); 00070 00071 00072 // Now we also need an AreaDefinition to define the properties of the 00073 // area we want 00074 // 00075 // This is made of 2 building blocks: 00076 // - the area type: 00077 // passive, active, active with explicit ghosts, or Voronoi area 00078 // - the specifications: 00079 // a VoronoiSpec or a GhostedAreaSpec for the 3 ghost-bases ones 00080 // 00081 //---------------------------------------------------------- For 00082 // GhostedAreaSpec (as below), the minimal info you have to provide 00083 // is up to what rapidity ghosts are placed. 00084 // Other commonm parameters (that mostly have an impact on the 00085 // precision on the area) include the number of repetitions 00086 // (i.e. the number of different sets of ghosts that are used) and 00087 // the ghost density (controlled through the ghost_area). 00088 // Other, more exotic, parameters (not shown here) control how ghosts 00089 // are placed. 00090 // 00091 // The ghost rapidity interval should be large enough to cover the 00092 // jets for which you want to calculate. E.g. if you want to 00093 // calculate the area of jets up to |y|=4, you need to put ghosts up 00094 // to at least 4+R (or, optionally, up to the largest particle 00095 // rapidity if this is smaller). 00096 double maxrap = 5.0; 00097 unsigned int n_repeat = 3; // default is 1 00098 double ghost_area = 0.01; // this is the default 00099 fastjet::GhostedAreaSpec area_spec(maxrap, n_repeat, ghost_area); 00100 00101 fastjet::AreaDefinition area_def(fastjet::active_area, area_spec); 00102 00103 // run the jet clustering with the above jet and area definitions 00104 // 00105 // The only change is the usage of a ClusterSequenceArea rather than 00106 //a ClusterSequence 00107 //---------------------------------------------------------- 00108 fastjet::ClusterSequenceArea clust_seq(input_particles, jet_def, area_def); 00109 00110 00111 // get the resulting jets ordered in pt 00112 //---------------------------------------------------------- 00113 double ptmin = 5.0; 00114 vector<fastjet::PseudoJet> inclusive_jets = sorted_by_pt(clust_seq.inclusive_jets(ptmin)); 00115 00116 00117 // tell the user what was done 00118 // - the description of the algorithm and area used 00119 // - extract the inclusive jets with pt > 5 GeV 00120 // show the output as 00121 // {index, rap, phi, pt, number of constituents} 00122 //---------------------------------------------------------- 00123 cout << endl; 00124 cout << "Ran " << jet_def.description() << endl; 00125 cout << "Area: " << area_def.description() << endl << endl; 00126 00127 // label the columns 00128 printf("%5s %15s %15s %15s %15s %15s\n","jet #", "rapidity", "phi", "pt", "area", "area error"); 00129 00130 // print out the details for each jet 00131 for (unsigned int i = 0; i < inclusive_jets.size(); i++) { 00132 printf("%5u %15.8f %15.8f %15.8f %15.8f %15.8f\n", i, 00133 inclusive_jets[i].rap(), inclusive_jets[i].phi(), inclusive_jets[i].perp(), 00134 inclusive_jets[i].area(), inclusive_jets[i].area_error()); 00135 } 00136 00137 return 0; 00138 }