FastJet 3.4.2
02-jetdef.cc
Go to the documentation of this file.
1//----------------------------------------------------------------------
2/// \file
3/// \page Example02 02 - changing the jet definition
4///
5/// fastjet basic example program:
6/// illustration of the usage of how to change the jet definition
7/// used for the clustering (see also fastjet::JetDefinition)
8///
9/// run it with : ./02-jetdef < data/single-event.dat
10///
11/// Source code: 02-jetdef.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
45#include "fastjet/ClusterSequence.hh"
46#include <iostream> // needed for io
47#include <cstdio> // needed for io
48
49using namespace std;
50
51/// an example program showing how to use fastjet
52int main(){
53
54 // read in input particles
55 //----------------------------------------------------------
56 vector<fastjet::PseudoJet> input_particles;
57
58 double px, py , pz, E;
59 while (cin >> px >> py >> pz >> E) {
60 // create a fastjet::PseudoJet with these components and put it onto
61 // back of the input_particles vector
62 input_particles.push_back(fastjet::PseudoJet(px,py,pz,E));
63 }
64
65
66 // create a jet definition:
67 // a jet algorithm with a given radius parameter
68 //----------------------------------------------------------
69
70 // select a jet algorithm to use
71 //
72 // this could be one of
73 // {kt_algorithm, cambridge_algorithm, antikt_algorithm,
74 // genkt_algorithm, ee_kt_algorithm, ee_genkt_algorithm}
75 // see example 03-plugin.cc for extra options using plugins
76 // instead of naive algorithms)
78
79 // when appropriate select a radius to use
80 //
81 // this would not be mandatory for e^+ e^- algorithms (see 05-eplus_eminus.cc)
82 double R = 0.6;
83
84 // select an __optional__ recombination scheme
85 //
86 // this could be one of
87 // {E_scheme, pt_scheme, pt2_scheme, Et_scheme, Et2_scheme, BIpt_scheme,
88 // BIpt2_scheme, WTA_pt_scheme, WTA_E_scheme, WTA_modp_scheme,
89 // external_sheme}
90 //
91 // Notes:
92 // - for the usage of a user-defined recombination scheme
93 // (external_scheme), see 11-boosted_higgs.cc
94 // - WTA_E_scheme, WTA_modp_scheme are meant for e+e- clusterings
95 //
96 // By default, the E_scheme is used
98
99 // select an __optional__ strategy
100 //
101 // this could be chosen among
102 // {N2MinHeapTiled, N2Tiled, N2PoorTiled, N2Plain, N3Dumb,
103 // Best,
104 // NlnN, NlnN3pi, NlnN4pi, NlnNCam4pi, NlnNCam2pi2R, NlnNCam}
105 //
106 // By default, the Best strategy is chosen and we advise to keep
107 // that default unless you are targeting a very specific usage. Note
108 // also that the N log (N) strategies for algorithms other than
109 // Cambridge/Aachen need CGAL support.
111
112 // create the JetDefinition from the above information
113 fastjet::JetDefinition jet_def(jet_alg, R, recomb_scheme, strategy);
114
115
116 // run the jet clustering with the above jet definition
117 //----------------------------------------------------------
118 fastjet::ClusterSequence clust_seq(input_particles, jet_def);
119
120
121 // get the resulting jets ordered in pt
122 //----------------------------------------------------------
123 double ptmin = 5.0;
124 vector<fastjet::PseudoJet> inclusive_jets = sorted_by_pt(clust_seq.inclusive_jets(ptmin));
125
126
127 // tell the user what was done
128 // - the description of the algorithm used
129 // - extract the inclusive jets with pt > 5 GeV
130 // show the output as
131 // {index, rap, phi, pt}
132 //----------------------------------------------------------
133 cout << "Ran " << jet_def.description() << endl;
134
135 // label the columns
136 printf("%5s %15s %15s %15s\n","jet #", "rapidity", "phi", "pt");
137
138 // print out the details for each jet
139 for (unsigned int i = 0; i < inclusive_jets.size(); i++) {
140 printf("%5u %15.8f %15.8f %15.8f\n",
141 i, inclusive_jets[i].rap(), inclusive_jets[i].phi(),
142 inclusive_jets[i].perp());
143 }
144
145 return 0;
146}
int main()
an example program showing how to use fastjet
Definition: 02-jetdef.cc:52
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.
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
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
JetAlgorithm
the various families of jet-clustering algorithm
@ kt_algorithm
the longitudinally invariant kt algorithm