00001 //---------------------------------------------------------------------- 00002 /// \file 00003 /// \page Example08 08 - using the Selector tool 00004 /// 00005 /// fastjet sample program to illustrate the use of fastjet::Selector 00006 /// 00007 /// Note that to use fastjet::Selector (and fastjet tools), you need the 00008 /// fastjettools library. It is included by default by 00009 /// fastjet-config --libs 00010 /// 00011 /// run it with : ./08-selector < data/single-event.dat 00012 /// 00013 /// Source code: 08-selector.cc 00014 //---------------------------------------------------------------------- 00015 00016 //STARTHEADER 00017 // $Id: 08-selector.cc 2003 2011-03-11 16:10:08Z soyez $ 00018 // 00019 // Copyright (c) 2005-2011, Matteo Cacciari, Gavin Salam and Gregory Soyez 00020 // 00021 //---------------------------------------------------------------------- 00022 // This file is part of FastJet. 00023 // 00024 // FastJet is free software; you can redistribute it and/or modify 00025 // it under the terms of the GNU General Public License as published by 00026 // the Free Software Foundation; either version 2 of the License, or 00027 // (at your option) any later version. 00028 // 00029 // The algorithms that underlie FastJet have required considerable 00030 // development and are described in hep-ph/0512210. If you use 00031 // FastJet as part of work towards a scientific publication, please 00032 // include a citation to the FastJet paper. 00033 // 00034 // FastJet is distributed in the hope that it will be useful, 00035 // but WITHOUT ANY WARRANTY; without even the implied warranty of 00036 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the 00037 // GNU General Public License for more details. 00038 // 00039 // You should have received a copy of the GNU General Public License 00040 // along with FastJet; if not, write to the Free Software 00041 // Foundation, Inc.: 00042 // 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA 00043 //---------------------------------------------------------------------- 00044 //ENDHEADER 00045 00046 #include "fastjet/ClusterSequence.hh" 00047 #include "fastjet/Selector.hh" 00048 #include <iostream> // needed for io 00049 #include <cstdio> // needed for io 00050 00051 using namespace std; 00052 00053 int main (int argc, char ** argv) { 00054 00055 // read in input particles 00056 //---------------------------------------------------------- 00057 vector<fastjet::PseudoJet> input_particles; 00058 00059 double px, py , pz, E; 00060 while (cin >> px >> py >> pz >> E) { 00061 // create a fastjet::PseudoJet with these components and put it onto 00062 // back of the input_particles vector 00063 input_particles.push_back(fastjet::PseudoJet(px,py,pz,E)); 00064 } 00065 00066 // Selector application #1: keep particles within a given acceptance 00067 // e.g. all particles with 1<|y|<2.5, all particles with pt>1 for |y|<1 00068 //---------------------------------------------------------- 00069 fastjet::Selector particle_selector = fastjet::SelectorAbsRapRange(1.0,2.5) 00070 || (fastjet::SelectorAbsRapMax(1.0) && fastjet::SelectorPtMin(1.0)); 00071 cout << input_particles.size() << " particles before selector" << endl; 00072 input_particles = particle_selector(input_particles); 00073 cout << input_particles.size() << " particles after selector" << endl; 00074 00075 // create a jet definition: 00076 // a jet algorithm with a given radius parameter 00077 //---------------------------------------------------------- 00078 double R = 0.6; 00079 fastjet::JetDefinition jet_def(fastjet::kt_algorithm, R); 00080 00081 00082 // run the jet clustering with the above jet definition 00083 //---------------------------------------------------------- 00084 fastjet::ClusterSequence clust_seq(input_particles, jet_def); 00085 00086 00087 // get the 5 hardest jets within |y|<2 00088 // 00089 // Note that this nicely illustrates that you should watch out that 00090 // Selectors do not necessarily commute. 00091 // 00092 // The && operator behaves like a logical and i.e. keeps objects 00093 // that satisfy both criteria (independently). It does commute. 00094 // 00095 // The * operator applies Selectors successively (starting from the 00096 // rightmost as in a usual operator product). Here, order may matter. 00097 //---------------------------------------------------------- 00098 fastjet::Selector jet_selector = fastjet::SelectorNHardest(5) * fastjet::SelectorAbsRapMax(2.0); 00099 vector<fastjet::PseudoJet> inclusive_jets = sorted_by_pt(jet_selector(clust_seq.inclusive_jets())); 00100 00101 00102 // tell the user what was done 00103 // - the description of the algorithm used 00104 // - the description of teh selectors used 00105 // - the jets 00106 // show the output as 00107 // {index, rap, phi, pt} 00108 //---------------------------------------------------------- 00109 double rapmin, rapmax; 00110 00111 cout << "Ran " << jet_def.description() << endl; 00112 00113 cout << "Selected particles: " << particle_selector.description() << endl; 00114 particle_selector.get_rapidity_extent(rapmin, rapmax); 00115 cout << " with a total rapidity range of [" << rapmin << ", " << rapmax << "]" << endl; 00116 00117 cout << "Selected jets: " << jet_selector.description() << endl; 00118 jet_selector.get_rapidity_extent(rapmin, rapmax); 00119 cout << " with a total rapidity range of [" << rapmin << ", " << rapmax << "]" << endl; 00120 00121 // label the columns 00122 printf("%5s %15s %15s %15s\n","jet #", "rapidity", "phi", "pt"); 00123 00124 // print out the details for each jet 00125 for (unsigned int i = 0; i < inclusive_jets.size(); i++) { 00126 printf("%5u %15.8f %15.8f %15.8f\n", 00127 i, inclusive_jets[i].rap(), inclusive_jets[i].phi(), 00128 inclusive_jets[i].perp()); 00129 } 00130 00131 return 0; 00132 }