FastJet 3.4.2
fastjet_example.cc
Go to the documentation of this file.
1/// \file
2/// \page Examples FastJet examples
3///
4/// The FastJet examples have been organised by order of complexity,
5/// starting by the simplest case and introducing features one after
6/// another.
7/// - \subpage Example01
8/// - \subpage Example02
9/// - \subpage Example03
10/// - \subpage Example04
11/// - \subpage Example05
12/// - \subpage Example06
13/// - \subpage Example07 (\subpage Example07old "old version")
14/// - \subpage Example08
15/// - \subpage Example09
16/// - \subpage Example10
17/// - \subpage Example11
18/// - \subpage Example12 (\subpage Example12old "old version")
19/// - \subpage Example13
20/// - \subpage Example14
21
22//FJSTARTHEADER
23// $Id$
24//
25// Copyright (c) 2005-2023, Matteo Cacciari, Gavin P. Salam and Gregory Soyez
26//
27//----------------------------------------------------------------------
28// This file is part of FastJet.
29//
30// FastJet is free software; you can redistribute it and/or modify
31// it under the terms of the GNU General Public License as published by
32// the Free Software Foundation; either version 2 of the License, or
33// (at your option) any later version.
34//
35// The algorithms that underlie FastJet have required considerable
36// development. They are described in the original FastJet paper,
37// hep-ph/0512210 and in the manual, arXiv:1111.6097. If you use
38// FastJet as part of work towards a scientific publication, please
39// quote the version you use and include a citation to the manual and
40// optionally also to hep-ph/0512210.
41//
42// FastJet is distributed in the hope that it will be useful,
43// but WITHOUT ANY WARRANTY; without even the implied warranty of
44// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
45// GNU General Public License for more details.
46//
47// You should have received a copy of the GNU General Public License
48// along with FastJet. If not, see <http://www.gnu.org/licenses/>.
49//----------------------------------------------------------------------
50//FJENDHEADER
51
52
53//----------------------------------------------------------------------
54// fastjet example program.
55//
56// Compile it with: make fastjet_example
57// run it with : ./fastjet_example < data/single-event.dat
58//
59// People who are familiar with the ktjet package are encouraged to
60// compare this file to the ktjet_example.cc program which does the
61// same thing in the ktjet framework.
62//----------------------------------------------------------------------
63#include "fastjet/PseudoJet.hh"
64#include "fastjet/ClusterSequence.hh"
65#include<iostream> // needed for io
66#include<sstream> // needed for internal io
67#include<vector>
68#include <cstdio>
69
70using namespace std;
71
72// a declaration of a function that pretty prints a list of jets
73void print_jets (const vector<fastjet::PseudoJet> &);
74
75/// an example program showing how to use fastjet
76int main () {
77
78 vector<fastjet::PseudoJet> input_particles;
79
80 // Read in input particles
81 double px, py , pz, E;
82 while (cin >> px >> py >> pz >> E) {
83 // create a fastjet::PseudoJet with these components and put it onto
84 // back of the input_particles vector
85 input_particles.push_back(fastjet::PseudoJet(px,py,pz,E));
86 }
87
88 // create an object that represents your choice of jet algorithm and
89 // the associated parameters
90 double Rparam = 1.0;
93 fastjet::JetDefinition jet_def(fastjet::kt_algorithm, Rparam, recomb_scheme, strategy);
94
95 // run the jet clustering with the above jet definition
96 fastjet::ClusterSequence clust_seq(input_particles, jet_def);
97
98 // tell the user what was done
99 cout << "Ran " << jet_def.description() << endl;
100 cout << "Strategy adopted by FastJet was "<<
101 clust_seq.strategy_string()<<endl<<endl;
102
103 // extract the inclusive jets with pt > 5 GeV
104 double ptmin = 5.0;
105 vector<fastjet::PseudoJet> inclusive_jets = clust_seq.inclusive_jets(ptmin);
106
107 // print them out
108 cout << "Printing inclusive jets with pt > "<< ptmin<<" GeV\n";
109 cout << "---------------------------------------\n";
110 print_jets(inclusive_jets);
111 cout << endl;
112
113 // extract the exclusive jets with dcut = 25 GeV^2
114 double dcut = 25.0;
115 vector<fastjet::PseudoJet> exclusive_jets = clust_seq.exclusive_jets(dcut);
116
117 // print them out
118 cout << "Printing exclusive jets with dcut = "<< dcut<<" GeV^2\n";
119 cout << "--------------------------------------------\n";
120 print_jets(exclusive_jets);
121
122
123}
124
125
126//----------------------------------------------------------------------
127/// a function that pretty prints a list of jets
128void print_jets (const vector<fastjet::PseudoJet> & jets) {
129
130 // sort jets into increasing pt
131 vector<fastjet::PseudoJet> sorted_jets = sorted_by_pt(jets);
132
133 // label the columns
134 printf("%5s %15s %15s %15s %15s\n","jet #", "rapidity",
135 "phi", "pt", "n constituents");
136
137 // print out the details for each jet
138 for (unsigned int i = 0; i < sorted_jets.size(); i++) {
139 // the following is not super efficient since it creates an
140 // intermediate constituents vector
141 int n_constituents = sorted_jets[i].constituents().size();
142 printf("%5u %15.8f %15.8f %15.8f %8u\n",
143 i, sorted_jets[i].rap(), sorted_jets[i].phi(),
144 sorted_jets[i].perp(), n_constituents);
145 }
146
147}
deals with clustering
std::vector< PseudoJet > inclusive_jets(const double ptmin=0.0) const
return a vector of all jets (in the sense of the inclusive algorithm) with pt >= ptmin.
std::string strategy_string() const
return the name of the strategy used to cluster the event
std::vector< PseudoJet > exclusive_jets(const double dcut) const
return a vector of all jets (in the sense of the exclusive algorithm) that would be obtained when run...
class that is intended to hold a full definition of the jet clusterer
std::string description() const
return a textual description of the current jet definition
Class to contain pseudojets, including minimal information of use to jet-clustering routines.
Definition: PseudoJet.hh:68
void print_jets(const vector< fastjet::PseudoJet > &)
a function that pretty prints a list of jets
int main()
an example program showing how to use fastjet
Strategy
the various options for the algorithmic strategy to adopt in clustering events with kt and cambridge ...
@ Best
automatic selection of the best (based on N), including the LazyTiled strategies that are new to FJ3....
RecombinationScheme
The various recombination schemes.
@ E_scheme
summing the 4-momenta
@ kt_algorithm
the longitudinally invariant kt algorithm