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