FastJet 3.0beta1
|
00001 //---------------------------------------------------------------------- 00002 /// \file 00003 /// \page Example10 10 - extracting subjets 00004 /// 00005 /// fastjet example program to show how to access subjets; 00006 /// 00007 /// See also 12-boosted_higgs.cc to see the use of subjets for 00008 /// identifying boosted higgs (and other objects) 00009 /// 00010 /// run it with : ./10-subjets < data/single-event.dat 00011 /// 00012 /// Source code: 10-subjets.cc 00013 //---------------------------------------------------------------------- 00014 00015 //STARTHEADER 00016 // $Id: 10-subjets.cc 2470 2011-07-26 09:40:41Z cacciari $ 00017 // 00018 // Copyright (c) 2005-2011, Matteo Cacciari, Gavin Salam and Gregory Soyez 00019 // 00020 //---------------------------------------------------------------------- 00021 // This file is part of FastJet. 00022 // 00023 // FastJet is free software; you can redistribute it and/or modify 00024 // it under the terms of the GNU General Public License as published by 00025 // the Free Software Foundation; either version 2 of the License, or 00026 // (at your option) any later version. 00027 // 00028 // The algorithms that underlie FastJet have required considerable 00029 // development and are described in hep-ph/0512210. If you use 00030 // FastJet as part of work towards a scientific publication, please 00031 // include a citation to the FastJet paper. 00032 // 00033 // FastJet is distributed in the hope that it will be useful, 00034 // but WITHOUT ANY WARRANTY; without even the implied warranty of 00035 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the 00036 // GNU General Public License for more details. 00037 // 00038 // You should have received a copy of the GNU General Public License 00039 // along with FastJet; if not, write to the Free Software 00040 // Foundation, Inc.: 00041 // 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA 00042 //---------------------------------------------------------------------- 00043 //ENDHEADER 00044 00045 #include "fastjet/ClusterSequence.hh" 00046 #include <iostream> // needed for io 00047 #include <cstdio> // needed for io 00048 00049 using namespace std; 00050 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 // for subjet studies, Cambridge/Aachen is the natural algorithm 00067 //---------------------------------------------------------- 00068 double R = 1.0; 00069 fastjet::JetDefinition jet_def(fastjet::cambridge_algorithm, R); 00070 00071 00072 // run the jet clustering with the above jet definition 00073 // and get the jets above 5 GeV 00074 //---------------------------------------------------------- 00075 fastjet::ClusterSequence clust_seq(input_particles, jet_def); 00076 double ptmin = 6.0; 00077 vector<fastjet::PseudoJet> inclusive_jets = sorted_by_pt(clust_seq.inclusive_jets(ptmin)); 00078 00079 // extract the subjets at a smaller angular scale (Rsub=0.5) 00080 // 00081 // This is done by using ClusterSequence::exclusive_subjets(dcut): 00082 // for the Cambridge/Aachen algorithm, running with R and then 00083 // asking for exclusive subjets with dcut should give the same 00084 // subjets as rerunning the algorithm with R'=R*sqrt(dcut) on the 00085 // jet's constituents. 00086 // 00087 // At the same time we output a summary of what has been done and the 00088 // resulting subjets 00089 //---------------------------------------------------------- 00090 double Rsub = 0.5; 00091 double dcut = pow(Rsub/R,2); 00092 00093 // a "header" for the output 00094 cout << "Ran " << jet_def.description() << endl; 00095 cout << "Showing the jets above " << ptmin << " GeV" << endl; 00096 cout << "And their subjets for Rsub = " << Rsub << endl; 00097 printf("%10s %15s %15s %15s %15s\n","jet #", "rapidity", "phi", "pt", "n constituents"); 00098 00099 // show the jets and their subjets 00100 for (unsigned int i = 0; i < inclusive_jets.size(); i++) { 00101 // get the subjets 00102 vector<fastjet::PseudoJet> subjets = sorted_by_pt(inclusive_jets[i].exclusive_subjets(dcut)); 00103 00104 cout << endl; 00105 // print the jet and its subjets 00106 printf("%5u %15.8f %15.8f %15.8f %8u\n", i, 00107 inclusive_jets[i].rap(), inclusive_jets[i].phi(), 00108 inclusive_jets[i].perp(), inclusive_jets[i].constituents().size()); 00109 00110 for (unsigned int j=0; j<subjets.size(); j++) 00111 printf(" sub%4u %15.8f %15.8f %15.8f %8u\n", j, 00112 subjets[j].rap(), subjets[j].phi(), 00113 subjets[j].perp(), subjets[j].constituents().size()); 00114 } 00115 00116 return 0; 00117 }