|
fastjet 2.4.3
|
#include "fastjet/PseudoJet.hh"#include "fastjet/ClusterSequence.hh"#include <iostream>#include <sstream>#include <vector>#include <cstdio>
Go to the source code of this file.
Functions | |
| void | print_jets (const fj::ClusterSequence &clust_seq, const vector< fj::PseudoJet > &jets) |
| a function that pretty prints a list of jets | |
| void | print_jet (const fj::ClusterSequence &clust_seq, const fj::PseudoJet &jet) |
| print a single jet | |
| void | print_jets_and_sub (const fj::ClusterSequence &clust_seq, const vector< fj::PseudoJet > &jets, double dcut) |
| 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 70 of file fastjet_subjets.cc.
References cambridge_algorithm, print_jets(), and print_jets_and_sub().
{
vector<fj::PseudoJet> input_particles;
// read in input particles
double px, py , pz, E;
while (cin >> px >> py >> pz >> E) {
// create a fj::PseudoJet with these components and put it onto
// back of the input_particles vector
input_particles.push_back(fj::PseudoJet(px,py,pz,E));
}
// We'll do two subjet analyses, physically rather different
double R = 1.0;
//fj::JetDefinition kt_def(fj::kt_algorithm, R);
fj::JetDefinition cam_def(fj::cambridge_algorithm, R);
// run the jet clustering with the above jet definition
//fj::ClusterSequence kt_seq(input_particles, kt_def);
fj::ClusterSequence cam_seq(input_particles, cam_def);
// tell the user what was done
cout << "Ran " << cam_def.description() << endl;
cout << "Strategy adopted by FastJet was "<<
cam_seq.strategy_string()<<endl<<endl;
// extract the inclusive jets with pt > 5 GeV
double ptmin = 5.0;
vector<fj::PseudoJet> inclusive_jets = cam_seq.inclusive_jets(ptmin);
// for the subjets
double smallR = 0.4;
double dcut_cam = pow(smallR/R,2);
// print them out
cout << "Printing inclusive jets (R = "<<R<<") with pt > "<< ptmin<<" GeV\n";
cout << "and their subjets with smallR = " << smallR << "\n";
cout << "---------------------------------------\n";
print_jets_and_sub(cam_seq, inclusive_jets, dcut_cam);
cout << endl;
// print them out
vector<fj::PseudoJet> exclusive_jets = cam_seq.exclusive_jets(dcut_cam);
cout << "Printing exclusive jets with dcut = "<< dcut_cam<<" \n";
cout << "--------------------------------------------\n";
print_jets(cam_seq, exclusive_jets);
}
| void print_jet | ( | const fj::ClusterSequence & | clust_seq, |
| const fj::PseudoJet & | jet | ||
| ) |
print a single jet
Definition at line 173 of file fastjet_subjets.cc.
References fastjet::PseudoJet::perp(), fastjet::PseudoJet::phi(), and fastjet::PseudoJet::rap().
Referenced by print_jets(), and print_jets_and_sub().
| void print_jets | ( | const fj::ClusterSequence & | clust_seq, |
| const vector< fj::PseudoJet > & | jets | ||
| ) |
a function that pretty prints a list of jets
Definition at line 124 of file fastjet_subjets.cc.
References print_jet(), and sorted_by_pt().
{
// sort jets into increasing pt
vector<fj::PseudoJet> sorted_jets = sorted_by_pt(jets);
// label the columns
printf("%5s %15s %15s %15s %15s\n","jet #", "rapidity",
"phi", "pt", "n constituents");
// print out the details for each jet
for (unsigned int i = 0; i < sorted_jets.size(); i++) {
printf("%5u ",i);
print_jet(clust_seq, sorted_jets[i]);
}
}
| void print_jets_and_sub | ( | const fj::ClusterSequence & | clust_seq, |
| const vector< fj::PseudoJet > & | jets, | ||
| double | dcut | ||
| ) |
a function that pretty prints a list of jets
Definition at line 145 of file fastjet_subjets.cc.
References print_jet(), and sorted_by_pt().
Referenced by main().
{
// sort jets into increasing pt
vector<fj::PseudoJet> sorted_jets = sorted_by_pt(jets);
// label the columns
printf("%5s %15s %15s %15s %15s\n","jet #", "rapidity",
"phi", "pt", "n constituents");
// print out the details for each jet
for (unsigned int i = 0; i < sorted_jets.size(); i++) {
printf("%5u ",i);
print_jet(clust_seq, sorted_jets[i]);
vector<fj::PseudoJet> subjets =
sorted_by_pt(clust_seq.exclusive_subjets(sorted_jets[i], dcut));
for (unsigned int j = 0; j < subjets.size(); j++) {
printf(" -sub-%02u ",j);
print_jet(clust_seq, subjets[j]);
}
}
}
1.7.3