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 #include <cstdio>
00048
00049 using namespace std;
00050
00051
00052
00053 namespace fj = fastjet;
00054
00055
00056
00057 void print_jets (const vector<fj::PseudoJet> &);
00058
00059
00060 void print_jet (const fj::PseudoJet & jet);
00061
00062
00063 void print_jets_and_sub (const vector<fj::PseudoJet> & jets,
00064 double dcut);
00065
00066
00067 int main (int argc, char ** argv) {
00068
00069 vector<fj::PseudoJet> input_particles;
00070
00071
00072 double px, py , pz, E;
00073 while (cin >> px >> py >> pz >> E) {
00074
00075
00076 input_particles.push_back(fj::PseudoJet(px,py,pz,E));
00077 }
00078
00079
00080 double R = 1.0;
00081
00082 fj::JetDefinition cam_def(fj::cambridge_algorithm, R);
00083
00084
00085
00086 fj::ClusterSequence cam_seq(input_particles, cam_def);
00087
00088
00089 cout << "Ran " << cam_def.description() << endl;
00090 cout << "Strategy adopted by FastJet was "<<
00091 cam_seq.strategy_string()<<endl<<endl;
00092
00093
00094 double ptmin = 5.0;
00095 vector<fj::PseudoJet> inclusive_jets = cam_seq.inclusive_jets(ptmin);
00096
00097
00098 double smallR = 0.4;
00099 double dcut_cam = pow(smallR/R,2);
00100
00101
00102 cout << "Printing inclusive jets (R = "<<R<<") with pt > "<< ptmin<<" GeV\n";
00103 cout << "and their subjets with smallR = " << smallR << "\n";
00104 cout << "---------------------------------------\n";
00105 print_jets_and_sub(inclusive_jets, dcut_cam);
00106 cout << endl;
00107
00108
00109
00110 vector<fj::PseudoJet> exclusive_jets = cam_seq.exclusive_jets(dcut_cam);
00111 cout << "Printing exclusive jets with dcut = "<< dcut_cam<<" \n";
00112 cout << "--------------------------------------------\n";
00113 print_jets(exclusive_jets);
00114
00115
00116 }
00117
00118
00119
00120
00121 void print_jets (const vector<fj::PseudoJet> & jets) {
00122
00123
00124 vector<fj::PseudoJet> sorted_jets = sorted_by_pt(jets);
00125
00126
00127 printf("%5s %15s %15s %15s %15s\n","jet #", "rapidity",
00128 "phi", "pt", "n constituents");
00129
00130
00131 for (unsigned int i = 0; i < sorted_jets.size(); i++) {
00132 printf("%5u ",i);
00133 print_jet(sorted_jets[i]);
00134 }
00135
00136 }
00137
00138
00139
00140
00141 void print_jets_and_sub (const vector<fj::PseudoJet> & jets,
00142 double dcut) {
00143
00144
00145 vector<fj::PseudoJet> sorted_jets = sorted_by_pt(jets);
00146
00147
00148 printf("%5s %15s %15s %15s %15s\n","jet #", "rapidity",
00149 "phi", "pt", "n constituents");
00150
00151
00152 for (unsigned int i = 0; i < sorted_jets.size(); i++) {
00153 printf("%5u ",i);
00154 print_jet(sorted_jets[i]);
00155 vector<fj::PseudoJet> subjets = sorted_by_pt(sorted_jets[i].exclusive_subjets(dcut));
00156 for (unsigned int j = 0; j < subjets.size(); j++) {
00157 printf(" -sub-%02u ",j);
00158 print_jet(subjets[j]);
00159 }
00160 }
00161
00162 }
00163
00164
00165
00166
00167 void print_jet (const fj::PseudoJet & jet) {
00168 int n_constituents = jet.constituents().size();
00169 printf("%15.8f %15.8f %15.8f %8u\n",
00170 jet.rap(), jet.phi(), jet.perp(), n_constituents);
00171 }