#include "FjPseudoJet.hh"#include "FjClusterSequence.hh"#include <iostream>#include <sstream>#include <vector>Include dependency graph for fastjet_example_v1_interface.cc:

Go to the source code of this file.
Functions | |
| void | print_jets (const FjClusterSequence &, const vector< FjPseudoJet > &) |
| 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 54 of file fastjet_example_v1_interface.cc. References Best, fastjet::ClusterSequence::inclusive_jets(), print_jets(), and fastjet::ClusterSequence::strategy_string(). 00054 {
00055
00056 vector<FjPseudoJet> input_particles;
00057
00058 // read in input particles
00059 double px, py , pz, E;
00060 while (cin >> px >> py >> pz >> E) {
00061 // create a FjPseudoJet with these components and put it onto
00062 // back of the input_particles vector
00063 input_particles.push_back(FjPseudoJet(px,py,pz,E));
00064 }
00065
00066 // run the jet clustering with option R=1.0 and strategy=Best
00067 double Rparam = 1.0;
00068 FjClusterSequence clust_seq(input_particles, Rparam, Best);
00069
00070 // tell the user what was done
00071 cout << "Strategy adopted by FastJet was "<<
00072 clust_seq.strategy_string()<<endl<<endl;
00073
00074 // extract the inclusive jets with pt > 5 GeV, sorted by pt
00075 double ptmin = 5.0;
00076 vector<FjPseudoJet> inclusive_jets = clust_seq.inclusive_jets(ptmin);
00077
00078 // print them out
00079 cout << "Printing inclusive jets with pt > "<< ptmin<<" GeV\n";
00080 cout << "---------------------------------------\n";
00081 print_jets(clust_seq, inclusive_jets);
00082 cout << endl;
00083
00084 // extract the exclusive jets with dcut = 25 GeV^2
00085 double dcut = 25.0;
00086 vector<FjPseudoJet> exclusive_jets = clust_seq.exclusive_jets(dcut);
00087
00088 // print them out
00089 cout << "Printing exclusive jets with dcut = "<< dcut<<" GeV^2\n";
00090 cout << "--------------------------------------------\n";
00091 print_jets(clust_seq, exclusive_jets);
00092
00093
00094 }
|
|
||||||||||||
|
a function that pretty prints a list of jets
Definition at line 99 of file fastjet_example_v1_interface.cc. References fastjet::ClusterSequence::constituents(), and fastjet::sorted_by_pt(). 00100 {
00101
00102 // sort jets into increasing pt
00103 vector<FjPseudoJet> sorted_jets = sorted_by_pt(jets);
00104
00105 // label the columns
00106 printf("%5s %15s %15s %15s %15s\n","jet #", "rapidity",
00107 "phi", "pt", "n constituents");
00108
00109 // print out the details for each jet
00110 for (unsigned int i = 0; i < sorted_jets.size(); i++) {
00111 int n_constituents = clust_seq.constituents(sorted_jets[i]).size();
00112 printf("%5u %15.8f %15.8f %15.8f %8u\n",
00113 i, sorted_jets[i].rap(), sorted_jets[i].phi(),
00114 sorted_jets[i].perp(), n_constituents);
00115 }
00116
00117 }
|
1.4.2