00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042 #include "fastjet/PseudoJet.hh"
00043 #include "fastjet/ClusterSequence.hh"
00044 #include<iostream>
00045 #include<sstream>
00046 #include<vector>
00047
00048 using namespace std;
00049
00050
00051 void print_jets (const fastjet::ClusterSequence &,
00052 const vector<fastjet::PseudoJet> &);
00053
00055 int main (int argc, char ** argv) {
00056
00057 vector<fastjet::PseudoJet> input_particles;
00058
00059
00060 double px, py , pz, E;
00061 while (cin >> px >> py >> pz >> E) {
00062
00063
00064 input_particles.push_back(fastjet::PseudoJet(px,py,pz,E));
00065 }
00066
00067
00068
00069 double Rparam = 1.0;
00070 fastjet::Strategy strategy = fastjet::Best;
00071 fastjet::RecombinationScheme recomb_scheme = fastjet::E_scheme;
00072 fastjet::JetDefinition jet_def(fastjet::kt_algorithm, Rparam, recomb_scheme, strategy);
00073
00074
00075 fastjet::ClusterSequence clust_seq(input_particles, jet_def);
00076
00077
00078 cout << "Ran " << jet_def.description() << endl;
00079 cout << "Strategy adopted by FastJet was "<<
00080 clust_seq.strategy_string()<<endl<<endl;
00081
00082
00083 double ptmin = 5.0;
00084 vector<fastjet::PseudoJet> inclusive_jets = clust_seq.inclusive_jets(ptmin);
00085
00086
00087 cout << "Printing inclusive jets with pt > "<< ptmin<<" GeV\n";
00088 cout << "---------------------------------------\n";
00089 print_jets(clust_seq, inclusive_jets);
00090 cout << endl;
00091
00092
00093 double dcut = 25.0;
00094 vector<fastjet::PseudoJet> exclusive_jets = clust_seq.exclusive_jets(dcut);
00095
00096
00097 cout << "Printing exclusive jets with dcut = "<< dcut<<" GeV^2\n";
00098 cout << "--------------------------------------------\n";
00099 print_jets(clust_seq, exclusive_jets);
00100
00101
00102 }
00103
00104
00105
00107 void print_jets (const fastjet::ClusterSequence & clust_seq,
00108 const vector<fastjet::PseudoJet> & jets) {
00109
00110
00111 vector<fastjet::PseudoJet> sorted_jets = sorted_by_pt(jets);
00112
00113
00114 printf("%5s %15s %15s %15s %15s\n","jet #", "rapidity",
00115 "phi", "pt", "n constituents");
00116
00117
00118 for (unsigned int i = 0; i < sorted_jets.size(); i++) {
00119 int n_constituents = clust_seq.constituents(sorted_jets[i]).size();
00120 printf("%5u %15.8f %15.8f %15.8f %8u\n",
00121 i, sorted_jets[i].rap(), sorted_jets[i].phi(),
00122 sorted_jets[i].perp(), n_constituents);
00123 }
00124
00125 }