FastJet 3.4.2
06-area.cc
Go to the documentation of this file.
1//----------------------------------------------------------------------
2/// \file
3/// \page Example06 06 - using jet areas
4///
5/// fastjet example program for jet areas
6/// It mostly illustrates the usage of the
7/// fastjet::AreaDefinition and fastjet::ClusterSequenceArea classes
8///
9/// run it with : ./06-area < data/single-event.dat
10///
11/// Source code: 06-area.cc
12//----------------------------------------------------------------------
13
14//FJSTARTHEADER
15// $Id$
16//
17// Copyright (c) 2005-2023, Matteo Cacciari, Gavin P. Salam and Gregory Soyez
18//
19//----------------------------------------------------------------------
20// This file is part of FastJet.
21//
22// FastJet is free software; you can redistribute it and/or modify
23// it under the terms of the GNU General Public License as published by
24// the Free Software Foundation; either version 2 of the License, or
25// (at your option) any later version.
26//
27// The algorithms that underlie FastJet have required considerable
28// development. They are described in the original FastJet paper,
29// hep-ph/0512210 and in the manual, arXiv:1111.6097. If you use
30// FastJet as part of work towards a scientific publication, please
31// quote the version you use and include a citation to the manual and
32// optionally also to hep-ph/0512210.
33//
34// FastJet is distributed in the hope that it will be useful,
35// but WITHOUT ANY WARRANTY; without even the implied warranty of
36// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
37// GNU General Public License for more details.
38//
39// You should have received a copy of the GNU General Public License
40// along with FastJet. If not, see <http://www.gnu.org/licenses/>.
41//----------------------------------------------------------------------
42//FJENDHEADER
43
44#include "fastjet/ClusterSequenceArea.hh" // use this instead of the "usual" ClusterSequence to get area support
45#include <iostream> // needed for io
46#include <cstdio> // needed for io
47
48using namespace std;
49
50/// an example program showing how to use fastjet
51int main(){
52
53 // read in input particles
54 //----------------------------------------------------------
55 vector<fastjet::PseudoJet> input_particles;
56
57 double px, py , pz, E;
58 while (cin >> px >> py >> pz >> E) {
59 // create a fastjet::PseudoJet with these components and put it onto
60 // back of the input_particles vector
61 input_particles.push_back(fastjet::PseudoJet(px,py,pz,E));
62 }
63
64
65 // create a jet definition:
66 // a jet algorithm with a given radius parameter
67 //----------------------------------------------------------
68 double R = 0.6;
70
71
72 // Now we also need an AreaDefinition to define the properties of the
73 // area we want
74 //
75 // This is made of 2 building blocks:
76 // - the area type:
77 // passive, active, active with explicit ghosts, or Voronoi area
78 // - the specifications:
79 // a VoronoiSpec or a GhostedAreaSpec for the 3 ghost-bases ones
80 //
81 //---------------------------------------------------------- For
82 // GhostedAreaSpec (as below), the minimal info you have to provide
83 // is up to what rapidity ghosts are placed.
84 // Other commonm parameters (that mostly have an impact on the
85 // precision on the area) include the number of repetitions
86 // (i.e. the number of different sets of ghosts that are used) and
87 // the ghost density (controlled through the ghost_area).
88 // Other, more exotic, parameters (not shown here) control how ghosts
89 // are placed.
90 //
91 // The ghost rapidity interval should be large enough to cover the
92 // jets for which you want to calculate. E.g. if you want to
93 // calculate the area of jets up to |y|=4, you need to put ghosts up
94 // to at least 4+R (or, optionally, up to the largest particle
95 // rapidity if this is smaller).
96 double maxrap = 5.0;
97 unsigned int n_repeat = 3; // default is 1
98 double ghost_area = 0.01; // this is the default
99 fastjet::GhostedAreaSpec area_spec(maxrap, n_repeat, ghost_area);
100
101 fastjet::AreaDefinition area_def(fastjet::active_area, area_spec);
102
103 // run the jet clustering with the above jet and area definitions
104 //
105 // The only change is the usage of a ClusterSequenceArea rather than
106 //a ClusterSequence
107 //----------------------------------------------------------
108 fastjet::ClusterSequenceArea clust_seq(input_particles, jet_def, area_def);
109
110
111 // get the resulting jets ordered in pt
112 //----------------------------------------------------------
113 double ptmin = 5.0;
114 vector<fastjet::PseudoJet> inclusive_jets = sorted_by_pt(clust_seq.inclusive_jets(ptmin));
115
116
117 // tell the user what was done
118 // - the description of the algorithm and area used
119 // - extract the inclusive jets with pt > 5 GeV
120 // show the output as
121 // {index, rap, phi, pt, number of constituents}
122 //----------------------------------------------------------
123 cout << endl;
124 cout << "Ran " << jet_def.description() << endl;
125 cout << "Area: " << area_def.description() << endl << endl;
126
127 // label the columns
128 printf("%5s %15s %15s %15s %15s %15s\n","jet #", "rapidity", "phi", "pt", "area", "area error");
129
130 // print out the details for each jet
131 for (unsigned int i = 0; i < inclusive_jets.size(); i++) {
132 printf("%5u %15.8f %15.8f %15.8f %15.8f %15.8f\n", i,
133 inclusive_jets[i].rap(), inclusive_jets[i].phi(), inclusive_jets[i].perp(),
134 inclusive_jets[i].area(), inclusive_jets[i].area_error());
135 }
136
137 return 0;
138}
int main()
an example program showing how to use fastjet
Definition: 06-area.cc:51
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.
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.
Parameters to configure the computation of jet areas using ghosts.
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
@ kt_algorithm
the longitudinally invariant kt algorithm