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 fastjet::ClusterSequence &clust_seq, const vector< fastjet::PseudoJet > &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 56 of file fastjet_example.cc.
References Best, E_scheme, kt_algorithm, and print_jets().
{ vector<fastjet::PseudoJet> input_particles; // read in input particles double px, py , pz, E; while (cin >> px >> py >> pz >> E) { // create a fastjet::PseudoJet with these components and put it onto // back of the input_particles vector input_particles.push_back(fastjet::PseudoJet(px,py,pz,E)); } // create an object that represents your choice of jet algorithm and // the associated parameters double Rparam = 1.0; fastjet::Strategy strategy = fastjet::Best; fastjet::RecombinationScheme recomb_scheme = fastjet::E_scheme; fastjet::JetDefinition jet_def(fastjet::kt_algorithm, Rparam, recomb_scheme, strategy); // run the jet clustering with the above jet definition fastjet::ClusterSequence clust_seq(input_particles, jet_def); // tell the user what was done cout << "Ran " << jet_def.description() << endl; cout << "Strategy adopted by FastJet was "<< clust_seq.strategy_string()<<endl<<endl; // extract the inclusive jets with pt > 5 GeV double ptmin = 5.0; vector<fastjet::PseudoJet> inclusive_jets = clust_seq.inclusive_jets(ptmin); // print them out cout << "Printing inclusive jets with pt > "<< ptmin<<" GeV\n"; cout << "---------------------------------------\n"; print_jets(clust_seq, inclusive_jets); cout << endl; // extract the exclusive jets with dcut = 25 GeV^2 double dcut = 25.0; vector<fastjet::PseudoJet> exclusive_jets = clust_seq.exclusive_jets(dcut); // print them out cout << "Printing exclusive jets with dcut = "<< dcut<<" GeV^2\n"; cout << "--------------------------------------------\n"; print_jets(clust_seq, exclusive_jets); }
void print_jets | ( | const fastjet::ClusterSequence & | clust_seq, |
const vector< fastjet::PseudoJet > & | jets | ||
) |
a function that pretty prints a list of jets
Definition at line 108 of file fastjet_example.cc.
References d0::inline_maths::phi(), and sorted_by_pt().
{ // sort jets into increasing pt vector<fastjet::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++) { int n_constituents = clust_seq.constituents(sorted_jets[i]).size(); printf("%5u %15.8f %15.8f %15.8f %8u\n", i, sorted_jets[i].rap(), sorted_jets[i].phi(), sorted_jets[i].perp(), n_constituents); } }