#include "fastjet/PseudoJet.hh"
#include "fastjet/ClusterSequenceActiveArea.hh"
#include <iostream>
#include <sstream>
#include <vector>
Include dependency graph for fastjet_areas.cc:
Go to the source code of this file.
Functions | |
void | print_jets (const fastjet::ClusterSequenceWithArea &, const vector< fastjet::PseudoJet > &) |
a function that pretty prints a list of jets | |
int | main (int argc, char **argv) |
an example program showing how to use fastjet |
|
an example program showing how to use fastjet
Definition at line 58 of file fastjet_areas.cc. References fastjet::Best, fastjet::ClusterSequence::inclusive_jets(), fastjet::kt_algorithm, print_jets(), and fastjet::ClusterSequence::strategy_string(). 00058 { 00059 00060 vector<fastjet::PseudoJet> input_particles; 00061 00062 // read in input particles 00063 double px, py , pz, E; 00064 while (cin >> px >> py >> pz >> E) { 00065 // create a fastjet::PseudoJet with these components and put it onto 00066 // back of the input_particles vector 00067 input_particles.push_back(fastjet::PseudoJet(px,py,pz,E)); 00068 } 00069 00070 // create an object that represents your choice of jet finder and 00071 // the associated parameters 00072 double Rparam = 1.0; 00073 fastjet::Strategy strategy = fastjet::Best; 00074 fastjet::JetDefinition jet_def(fastjet::kt_algorithm, Rparam, strategy); 00075 00076 // create an object that specifies how we to define the (active) area 00077 double ghost_etamax = 6.0; 00078 int active_area_repeats = 3; 00079 double ghost_area = 0.01; 00080 fastjet::ActiveAreaSpec area_spec(ghost_etamax, active_area_repeats, 00081 ghost_area); 00082 00083 // run the jet clustering with the above jet definition 00084 fastjet::ClusterSequenceActiveArea clust_seq(input_particles, 00085 jet_def, area_spec); 00086 00087 // tell the user what was done 00088 cout << "Strategy adopted by FastJet was "<< 00089 clust_seq.strategy_string()<<endl<<endl; 00090 00091 // extract the inclusive jets with pt > 5 GeV, sorted by pt 00092 double ptmin = 5.0; 00093 vector<fastjet::PseudoJet> inclusive_jets = clust_seq.inclusive_jets(ptmin); 00094 00095 // print them out 00096 cout << "Printing inclusive jets with pt > "<< ptmin<<" GeV\n"; 00097 cout << "---------------------------------------\n"; 00098 print_jets(clust_seq, inclusive_jets); 00099 cout << endl; 00100 00101 00102 }
|
|
a function that pretty prints a list of jets
Definition at line 107 of file fastjet_areas.cc. References fastjet::ClusterSequenceWithArea::area(), fastjet::ClusterSequenceWithArea::area_error(), and fastjet::sorted_by_pt(). Referenced by main(), and run_jet_finder(). 00108 { 00109 00110 // sort jets into increasing pt 00111 vector<fastjet::PseudoJet> jets = sorted_by_pt(unsorted_jets); 00112 00113 printf(" ijet rap phi Pt area +- err\n"); 00114 for (size_t j = 0; j < jets.size(); j++) { 00115 00116 double area = clust_seq.area(jets[j]); 00117 double area_error = clust_seq.area_error(jets[j]); 00118 00119 printf("%5u %9.5f %8.5f %10.3f %8.3f +- %6.3f\n",j,jets[j].rap(), 00120 jets[j].phi(),jets[j].perp(), area, area_error); 00121 } 00122 00123 00124 }
|