FastJet 3.4.2
fastjet_subjets.cc
1//FJSTARTHEADER
2// $Id$
3//
4// Copyright (c) 2005-2023, Matteo Cacciari, Gavin P. Salam and Gregory Soyez
5//
6//----------------------------------------------------------------------
7// This file is part of FastJet.
8//
9// FastJet is free software; you can redistribute it and/or modify
10// it under the terms of the GNU General Public License as published by
11// the Free Software Foundation; either version 2 of the License, or
12// (at your option) any later version.
13//
14// The algorithms that underlie FastJet have required considerable
15// development. They are described in the original FastJet paper,
16// hep-ph/0512210 and in the manual, arXiv:1111.6097. If you use
17// FastJet as part of work towards a scientific publication, please
18// quote the version you use and include a citation to the manual and
19// optionally also to hep-ph/0512210.
20//
21// FastJet is distributed in the hope that it will be useful,
22// but WITHOUT ANY WARRANTY; without even the implied warranty of
23// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
24// GNU General Public License for more details.
25//
26// You should have received a copy of the GNU General Public License
27// along with FastJet. If not, see <http://www.gnu.org/licenses/>.
28//----------------------------------------------------------------------
29//FJENDHEADER
30
31
32//----------------------------------------------------------------------
33// fastjet example program to show how to access subjets;
34//
35// See also fastjet_higgs_decomp.cc to see the use of subjets for
36// identifying boosted higgs (and other objects)
37//
38// Compile it with: make fastjet_subjets
39// run it with : ./fastjet_subjets < data/single-event.dat
40//
41//----------------------------------------------------------------------
42#include "fastjet/PseudoJet.hh"
43#include "fastjet/ClusterSequence.hh"
44#include<iostream> // needed for io
45#include<sstream> // needed for internal io
46#include<vector>
47#include <cstdio>
48
49using namespace std;
50
51// to avoid excessive typing, define an abbreviation for the
52// fastjet namespace
53namespace fj = fastjet;
54
55
56// a declaration of a function that pretty prints a list of jets
57void print_jets (const vector<fj::PseudoJet> &);
58
59// and this pretty prinst a single jet
60void print_jet (const fj::PseudoJet & jet);
61
62// pretty print the jets and their subjets
63void print_jets_and_sub (const vector<fj::PseudoJet> & jets,
64 double dcut);
65
66/// an example program showing how to use fastjet
67int main (int argc, char ** argv) {
68
69 vector<fj::PseudoJet> input_particles;
70
71 // read in input particles
72 double px, py , pz, E;
73 while (cin >> px >> py >> pz >> E) {
74 // create a fj::PseudoJet with these components and put it onto
75 // back of the input_particles vector
76 input_particles.push_back(fj::PseudoJet(px,py,pz,E));
77 }
78
79 // We'll do two subjet analyses, physically rather different
80 double R = 1.0;
81 //fj::JetDefinition kt_def(fj::kt_algorithm, R);
82 fj::JetDefinition cam_def(fj::cambridge_algorithm, R);
83
84 // run the jet clustering with the above jet definition
85 //fj::ClusterSequence kt_seq(input_particles, kt_def);
86 fj::ClusterSequence cam_seq(input_particles, cam_def);
87
88 // tell the user what was done
89 cout << "Ran " << cam_def.description() << endl;
90 cout << "Strategy adopted by FastJet was "<<
91 cam_seq.strategy_string()<<endl<<endl;
92
93 // extract the inclusive jets with pt > 5 GeV
94 double ptmin = 5.0;
95 vector<fj::PseudoJet> inclusive_jets = cam_seq.inclusive_jets(ptmin);
96
97 // for the subjets
98 double smallR = 0.4;
99 double dcut_cam = pow(smallR/R,2);
100
101 // print them out
102 cout << "Printing inclusive jets (R = "<<R<<") with pt > "<< ptmin<<" GeV\n";
103 cout << "and their subjets with smallR = " << smallR << "\n";
104 cout << "---------------------------------------\n";
105 print_jets_and_sub(inclusive_jets, dcut_cam);
106 cout << endl;
107
108
109 // print them out
110 vector<fj::PseudoJet> exclusive_jets = cam_seq.exclusive_jets(dcut_cam);
111 cout << "Printing exclusive jets with dcut = "<< dcut_cam<<" \n";
112 cout << "--------------------------------------------\n";
113 print_jets(exclusive_jets);
114
115
116}
117
118
119//----------------------------------------------------------------------
120/// a function that pretty prints a list of jets
121void print_jets (const vector<fj::PseudoJet> & jets) {
122
123 // sort jets into increasing pt
124 vector<fj::PseudoJet> sorted_jets = sorted_by_pt(jets);
125
126 // label the columns
127 printf("%5s %15s %15s %15s %15s\n","jet #", "rapidity",
128 "phi", "pt", "n constituents");
129
130 // print out the details for each jet
131 for (unsigned int i = 0; i < sorted_jets.size(); i++) {
132 printf("%5u ",i);
133 print_jet(sorted_jets[i]);
134 }
135
136}
137
138
139//----------------------------------------------------------------------
140/// a function that pretty prints a list of jets
141void print_jets_and_sub (const vector<fj::PseudoJet> & jets,
142 double dcut) {
143
144 // sort jets into increasing pt
145 vector<fj::PseudoJet> sorted_jets = sorted_by_pt(jets);
146
147 // label the columns
148 printf("%5s %15s %15s %15s %15s\n","jet #", "rapidity",
149 "phi", "pt", "n constituents");
150
151 // print out the details for each jet
152 for (unsigned int i = 0; i < sorted_jets.size(); i++) {
153 printf("%5u ",i);
154 print_jet(sorted_jets[i]);
155 vector<fj::PseudoJet> subjets = sorted_by_pt(sorted_jets[i].exclusive_subjets(dcut));
156 for (unsigned int j = 0; j < subjets.size(); j++) {
157 printf(" -sub-%02u ",j);
158 print_jet(subjets[j]);
159 }
160 }
161
162}
163
164
165//----------------------------------------------------------------------
166/// print a single jet
167void print_jet (const fj::PseudoJet & jet) {
168 int n_constituents = jet.constituents().size();
169 printf("%15.8f %15.8f %15.8f %8u\n",
170 jet.rap(), jet.phi(), jet.perp(), n_constituents);
171}
int main()
an example program showing how to use fastjet
Definition: 01-basic.cc:52
deals with clustering
class that is intended to hold a full definition of the jet clusterer
Class to contain pseudojets, including minimal information of use to jet-clustering routines.
Definition: PseudoJet.hh:68
virtual std::vector< PseudoJet > constituents() const
retrieve the constituents.
Definition: PseudoJet.cc:681
double rap() const
returns the rapidity or some large value when the rapidity is infinite
Definition: PseudoJet.hh:138
double phi() const
returns phi (in the range 0..2pi)
Definition: PseudoJet.hh:123
double perp() const
returns the scalar transverse momentum
Definition: PseudoJet.hh:158
void print_jets(const vector< fastjet::PseudoJet > &)
a function that pretty prints a list of jets
the FastJet namespace
vector< PseudoJet > sorted_by_pt(const vector< PseudoJet > &jets)
return a vector of jets sorted into decreasing kt2
Definition: PseudoJet.cc:871