FastJet 3.4.2
fastjet_areas.cc
1
2//FJSTARTHEADER
3// $Id$
4//
5// Copyright (c) 2005-2023, Matteo Cacciari, Gavin P. Salam and Gregory Soyez
6//
7//----------------------------------------------------------------------
8// This file is part of FastJet.
9//
10// FastJet is free software; you can redistribute it and/or modify
11// it under the terms of the GNU General Public License as published by
12// the Free Software Foundation; either version 2 of the License, or
13// (at your option) any later version.
14//
15// The algorithms that underlie FastJet have required considerable
16// development. They are described in the original FastJet paper,
17// hep-ph/0512210 and in the manual, arXiv:1111.6097. If you use
18// FastJet as part of work towards a scientific publication, please
19// quote the version you use and include a citation to the manual and
20// optionally also to hep-ph/0512210.
21//
22// FastJet is distributed in the hope that it will be useful,
23// but WITHOUT ANY WARRANTY; without even the implied warranty of
24// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
25// GNU General Public License for more details.
26//
27// You should have received a copy of the GNU General Public License
28// along with FastJet. If not, see <http://www.gnu.org/licenses/>.
29//----------------------------------------------------------------------
30//FJENDHEADER
31
32
33
34
35
36
37
38//----------------------------------------------------------------------
39// fastjet areas example program.
40//
41// Compile it with: make fastjet_areas
42// run it with : ./fastjet_areas < data/single-event.dat
43//
44//----------------------------------------------------------------------
45#include "fastjet/PseudoJet.hh"
46#include "fastjet/ClusterSequenceArea.hh"
47#include "fastjet/ClusterSequencePassiveArea.hh"
48
49// get info on how fastjet was configured
50#include "fastjet/config.h"
51
52#ifdef ENABLE_PLUGIN_SISCONE
53#include "fastjet/SISConePlugin.hh"
54#endif
55
56#include<iostream> // needed for io
57#include<sstream> // needed for internal io
58#include<vector>
59
60using namespace std;
61
62// a declaration of a function that pretty prints a list of jets
63void print_jets (const vector<fastjet::PseudoJet> &);
64
65/// an example program showing how to use fastjet
66int main () {
67
68 vector<fastjet::PseudoJet> input_particles;
69
70 // read in input particles
71 double px, py , pz, E;
72 while (cin >> px >> py >> pz >> E) {
73 // create a fastjet::PseudoJet with these components and put it onto
74 // back of the input_particles vector
75 input_particles.push_back(fastjet::PseudoJet(px,py,pz,E));
76 }
77
78 // create an object that represents your choice of jet algorithm, and
79 // the associated parameters
80 double Rparam = 1.0;
83 //fastjet::JetDefinition jet_def(fastjet::cambridge_algorithm, Rparam, strategy);
84 //fastjet::JetDefinition jet_def(fastjet::antikt_algorithm, Rparam, strategy);
85 //fastjet::JetDefinition jet_def(new fastjet::SISConePlugin(Rparam,0.75));
86
87 // create an object that specifies how we to define the area
89 bool use_voronoi = false;
90 if (!use_voronoi) {
91 double ghost_etamax = 6.0;
92 double ghost_area = 0.01;
93 int active_area_repeats = 1;
94
95 // now create the object that holds info about ghosts, and from that
96 // get an area definition
97 fastjet::GhostedAreaSpec ghost_spec(ghost_etamax, active_area_repeats,
98 ghost_area);
99 area_def = fastjet::AreaDefinition(fastjet::active_area,ghost_spec);
100 //area_def = fastjet::AreaDefinition(fastjet::passive_area,ghost_spec);
101 } else {
102 double effective_Rfact = 1.0;
103 area_def = fastjet::VoronoiAreaSpec(effective_Rfact);
104 }
105
106 // run the jet clustering with the above jet definition
107 fastjet::ClusterSequenceArea clust_seq(input_particles,
108 jet_def, area_def);
109 // you can also run the individual area classes directly
110 //fastjet::ClusterSequencePassiveArea clust_seq(input_particles, jet_def,
111 // area_def.ghost_spec());
112
113 // you may want to find out how much area in a given range (|y|<range)
114 // is empty of real jets (or corresponds to pure "ghost" jets).
115 //double range = 4.0;
116 //cout << clust_seq.empty_area(range) << endl;
117 //cout << clust_seq.n_empty_jets(range) << endl;
118
119 // tell the user what was done
120 cout << "Jet definition was: " << jet_def.description() << endl;
121 cout << "Area definition was: " << area_def.description() << endl;
122 cout << "Strategy adopted by FastJet was "<<
123 clust_seq.strategy_string()<<endl<<endl;
124
125 // extract the inclusive jets with pt > 5 GeV, sorted by pt
126 double ptmin = 5.0;
127 vector<fastjet::PseudoJet> inclusive_jets = clust_seq.inclusive_jets(ptmin);
128
129 // print them out
130 cout << "Printing inclusive jets with pt > "<< ptmin<<" GeV\n";
131 cout << "---------------------------------------\n";
132 print_jets(inclusive_jets);
133 cout << endl;
134
135
136 cout << "Number of unclustered particles: "
137 << clust_seq.unclustered_particles().size() << endl;
138
139
140}
141
142
143//----------------------------------------------------------------------
144/// a function that pretty prints a list of jets
145void print_jets (const vector<fastjet::PseudoJet> & unsorted_jets) {
146
147 // sort jets into increasing pt
148 vector<fastjet::PseudoJet> jets = sorted_by_pt(unsorted_jets);
149
150 printf(" ijet rap phi Pt area +- err\n");
151 for (unsigned int j = 0; j < jets.size(); j++) {
152
153 double area = jets[j].area();
154 double area_error = jets[j].area_error();
155
156 printf("%5u %9.5f %8.5f %10.3f %8.3f +- %6.3f\n",j,jets[j].rap(),
157 jets[j].phi(),jets[j].perp(), area, area_error);
158 }
159
160
161}
162
163
164
165
166
167
168
int main()
an example program showing how to use fastjet
Definition: 01-basic.cc:52
class that holds a generic area definition
std::string description() const
return a description of the current area definition
General class for user to obtain ClusterSequence with additional area information.
Parameters to configure the computation of jet areas using ghosts.
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
Specification for the computation of the Voronoi jet area.
void print_jets(const vector< fastjet::PseudoJet > &)
a function that pretty prints a list of jets
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....
@ E_scheme
summing the 4-momenta
@ kt_algorithm
the longitudinally invariant kt algorithm
vector< PseudoJet > sorted_by_pt(const vector< PseudoJet > &jets)
return a vector of jets sorted into decreasing kt2
Definition: PseudoJet.cc:871