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