#include "fastjet/PseudoJet.hh"
#include "fastjet/ClusterSequenceArea.hh"
#include "fastjet/ClusterSequencePassiveArea.hh"
#include "fastjet/config.h"
#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::ClusterSequenceAreaBase &clust_seq, const vector< fastjet::PseudoJet > &unsorted_jets) |
a function that pretty prints a list of jets | |
int | main (int argc, char **argv) |
an example program showing how to use fastjet |
int main | ( | int | argc, | |
char ** | argv | |||
) |
an example program showing how to use fastjet
Definition at line 67 of file fastjet_areas.cc.
References fastjet::Best, fastjet::AreaDefinition::description(), fastjet::JetDefinition::description(), fastjet::ClusterSequence::inclusive_jets(), fastjet::kt_algorithm, fastjet::passive_area, print_jets(), fastjet::ClusterSequence::strategy_string(), and fastjet::ClusterSequence::unclustered_particles().
00067 { 00068 00069 vector<fastjet::PseudoJet> input_particles; 00070 00071 // read in input particles 00072 double px, py , pz, E; 00073 while (cin >> px >> py >> pz >> E) { 00074 // create a fastjet::PseudoJet with these components and put it onto 00075 // back of the input_particles vector 00076 input_particles.push_back(fastjet::PseudoJet(px,py,pz,E)); 00077 } 00078 00079 // create an object that represents your choice of jet algorithm, and 00080 // the associated parameters 00081 double Rparam = 1.0; 00082 fastjet::Strategy strategy = fastjet::Best; 00083 fastjet::JetDefinition jet_def(fastjet::kt_algorithm, Rparam, strategy); 00084 //fastjet::JetDefinition jet_def(fastjet::cambridge_algorithm, Rparam, strategy); 00085 //fastjet::JetDefinition jet_def(fastjet::antikt_algorithm, Rparam, strategy); 00086 //fastjet::JetDefinition jet_def(new fastjet::SISConePlugin(1.0)); 00087 00088 // create an object that specifies how we to define the area 00089 fastjet::AreaDefinition area_def; 00090 bool use_voronoi = false; 00091 if (!use_voronoi) { 00092 double ghost_etamax = 7.0; 00093 int active_area_repeats = 1; 00094 //double ghost_area = 0.01; 00095 //int active_area_repeats = 100; 00096 double ghost_area = 0.05; 00097 //area_def = fastjet::GhostedAreaSpec(ghost_etamax, active_area_repeats, 00098 fastjet::GhostedAreaSpec ghost_spec(ghost_etamax, active_area_repeats, 00099 ghost_area); 00100 area_def = fastjet::AreaDefinition(fastjet::passive_area,ghost_spec); 00101 } else { 00102 double effective_Rfact = 1.0; 00103 area_def = fastjet::VoronoiAreaSpec(effective_Rfact); 00104 } 00105 00106 // run the jet clustering with the above jet definition 00107 fastjet::ClusterSequenceArea clust_seq(input_particles, 00108 jet_def, area_def); 00109 // run the jet clustering with the above jet definition 00110 //fastjet::ClusterSequencePassiveArea clust_seq(input_particles, jet_def, 00111 // area_def.active_spec()); 00112 00113 //cout << clust_seq.empty_area(4.0) << endl; 00114 //cout << clust_seq.n_empty_jets(4.0) << endl; 00115 00116 // tell the user what was done 00117 cout << "Jet definition was: " << jet_def.description() << endl; 00118 cout << "Area definition was: " << area_def.description() << endl; 00119 cout << "Strategy adopted by FastJet was "<< 00120 clust_seq.strategy_string()<<endl<<endl; 00121 00122 // extract the inclusive jets with pt > 5 GeV, sorted by pt 00123 double ptmin = 5.0; 00124 vector<fastjet::PseudoJet> inclusive_jets = clust_seq.inclusive_jets(ptmin); 00125 00126 // print them out 00127 cout << "Printing inclusive jets with pt > "<< ptmin<<" GeV\n"; 00128 cout << "---------------------------------------\n"; 00129 print_jets(clust_seq, inclusive_jets); 00130 cout << endl; 00131 00132 00133 cout << "Number of unclustered particles: " 00134 << clust_seq.unclustered_particles().size() << endl; 00135 00136 00137 }
void print_jets | ( | const fastjet::ClusterSequenceAreaBase & | clust_seq, | |
const vector< fastjet::PseudoJet > & | unsorted_jets | |||
) |
a function that pretty prints a list of jets
Definition at line 142 of file fastjet_areas.cc.
References fastjet::ClusterSequenceAreaBase::area(), fastjet::ClusterSequenceAreaBase::area_error(), and fastjet::sorted_by_pt().
Referenced by main(), and run_jet_finder().
00143 { 00144 00145 // sort jets into increasing pt 00146 vector<fastjet::PseudoJet> jets = sorted_by_pt(unsorted_jets); 00147 00148 printf(" ijet rap phi Pt area +- err\n"); 00149 for (size_t j = 0; j < jets.size(); j++) { 00150 00151 double area = clust_seq.area(jets[j]); 00152 double area_error = clust_seq.area_error(jets[j]); 00153 00154 printf("%5u %9.5f %8.5f %10.3f %8.3f +- %6.3f\n",j,jets[j].rap(), 00155 jets[j].phi(),jets[j].perp(), area, area_error); 00156 } 00157 00158 00159 }