FastJet 3.0beta1
|
00001 //---------------------------------------------------------------------- 00002 /// \file 00003 /// \page Example04 04 - accessing clustering information in a PseudoJet 00004 /// 00005 /// illustrate how a jet can carry information about its clustering 00006 /// 00007 /// We do it by associating a user index to each of the input particles 00008 /// and show what particles are in each jets (above 5 GeV) 00009 /// 00010 /// We also illustrate a few other features about how a fastjet::PseudoJet 00011 /// can access its underlying structure. 00012 /// 00013 /// run it with : ./04-constituents < data/single-event.dat 00014 /// 00015 /// Source code: 04-constituents.cc 00016 //---------------------------------------------------------------------- 00017 00018 //STARTHEADER 00019 // $Id: 04-constituents.cc 2003 2011-03-11 16:10:08Z soyez $ 00020 // 00021 // Copyright (c) 2005-2011, Matteo Cacciari, Gavin Salam and Gregory Soyez 00022 // 00023 //---------------------------------------------------------------------- 00024 // This file is part of FastJet. 00025 // 00026 // FastJet is free software; you can redistribute it and/or modify 00027 // it under the terms of the GNU General Public License as published by 00028 // the Free Software Foundation; either version 2 of the License, or 00029 // (at your option) any later version. 00030 // 00031 // The algorithms that underlie FastJet have required considerable 00032 // development and are described in hep-ph/0512210. If you use 00033 // FastJet as part of work towards a scientific publication, please 00034 // include a citation to the FastJet paper. 00035 // 00036 // FastJet is distributed in the hope that it will be useful, 00037 // but WITHOUT ANY WARRANTY; without even the implied warranty of 00038 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the 00039 // GNU General Public License for more details. 00040 // 00041 // You should have received a copy of the GNU General Public License 00042 // along with FastJet; if not, write to the Free Software 00043 // Foundation, Inc.: 00044 // 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA 00045 //---------------------------------------------------------------------- 00046 //ENDHEADER 00047 00048 #include "fastjet/ClusterSequence.hh" 00049 #include <iostream> // needed for io 00050 #include <cstdio> // needed for io 00051 00052 using namespace std; 00053 00054 /// an example program showing how to use fastjet 00055 int main (int argc, char ** argv) { 00056 00057 // read in input particles 00058 //---------------------------------------------------------- 00059 vector<fastjet::PseudoJet> input_particles; 00060 00061 valarray<double> fourvec(4); 00062 int index=0; 00063 while (cin >> fourvec[0] >> fourvec[1] >> fourvec[2] >> fourvec[3]) { 00064 // create a particle with the approprite 4-momentum and 00065 // set its user index to keep track of its index. 00066 // you can construct a PseudoJet from any object that allows subscripts 00067 // from [0] .. [3] (the last one must be the energy) 00068 fastjet::PseudoJet particle(fourvec); 00069 00070 particle.set_user_index(index); 00071 input_particles.push_back(particle); 00072 00073 index++; 00074 } 00075 00076 00077 // create a jet definition: 00078 // a jet algorithm with a given radius parameter 00079 //---------------------------------------------------------- 00080 double R = 0.6; 00081 fastjet::JetDefinition jet_def(fastjet::kt_algorithm, R); 00082 00083 00084 // run the jet clustering with the above jet definition 00085 //---------------------------------------------------------- 00086 fastjet::ClusterSequence clust_seq(input_particles, jet_def); 00087 00088 00089 // get the resulting jets ordered in pt 00090 //---------------------------------------------------------- 00091 double ptmin = 5.0; 00092 vector<fastjet::PseudoJet> inclusive_jets = sorted_by_pt(clust_seq.inclusive_jets(ptmin)); 00093 00094 00095 // tell the user what was done 00096 // - the description of the algorithm used 00097 // - extract the inclusive jets with pt > 5 GeV 00098 // show the output as 00099 // {index, rap, phi, pt, number of constituents} 00100 //---------------------------------------------------------- 00101 cout << "Ran " << jet_def.description() << endl << endl; 00102 00103 // label the columns 00104 printf("%5s %15s %15s %15s %15s\n","jet #", "rapidity", "phi", "pt", "n constituents"); 00105 printf(" indices of constituents\n\n"); 00106 00107 // print out the details for each jet 00108 for (unsigned int i = 0; i < inclusive_jets.size(); i++) { 00109 // get the constituents of the jet 00110 vector<fastjet::PseudoJet> constituents = inclusive_jets[i].constituents(); 00111 00112 printf("%5u %15.8f %15.8f %15.8f %8u\n", 00113 i, inclusive_jets[i].rap(), inclusive_jets[i].phi(), 00114 inclusive_jets[i].perp(), constituents.size()); 00115 00116 printf(" "); 00117 for (unsigned int j=0; j<constituents.size(); j++){ 00118 printf("%4u ", constituents[j].user_index()); 00119 if (j%10==9) printf("\n "); 00120 } 00121 printf("\n\n"); 00122 } 00123 00124 return 0; 00125 }