FastJet 3.4.2
05-eplus_eminus.cc
Go to the documentation of this file.
1//----------------------------------------------------------------------
2/// \file
3/// \page Example05 05 - using e+e- algorithms
4///
5/// illustrate the use of e^+ e^- algorithms
6///
7/// They mostly differ from the pp algorithm by the fact that rather
8/// than using a radius parameter and inclusive jets, they use
9/// exclusive jets in one of the following ways:
10/// - a fixed number of them
11/// - with a dcut
12/// - with a ycut
13///
14/// Note that natively, FastJet includes the kt (ee_kt_algorithm) and
15/// genkt (ee_genkt_algorithm) algorithms. Others (like Cambridge for
16/// e+ e-, Jade or SISCone in spherical coordinates) are available as
17/// plugins (see 03-plugin.cc)
18///
19/// run it with : ./05-eplus_eminus < data/single-ee-event.dat
20///
21/// Source code: 05-eplus_eminus.cc
22//----------------------------------------------------------------------
23
24//FJSTARTHEADER
25// $Id$
26//
27// Copyright (c) 2005-2023, Matteo Cacciari, Gavin P. Salam and Gregory Soyez
28//
29//----------------------------------------------------------------------
30// This file is part of FastJet.
31//
32// FastJet is free software; you can redistribute it and/or modify
33// it under the terms of the GNU General Public License as published by
34// the Free Software Foundation; either version 2 of the License, or
35// (at your option) any later version.
36//
37// The algorithms that underlie FastJet have required considerable
38// development. They are described in the original FastJet paper,
39// hep-ph/0512210 and in the manual, arXiv:1111.6097. If you use
40// FastJet as part of work towards a scientific publication, please
41// quote the version you use and include a citation to the manual and
42// optionally also to hep-ph/0512210.
43//
44// FastJet is distributed in the hope that it will be useful,
45// but WITHOUT ANY WARRANTY; without even the implied warranty of
46// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
47// GNU General Public License for more details.
48//
49// You should have received a copy of the GNU General Public License
50// along with FastJet. If not, see <http://www.gnu.org/licenses/>.
51//----------------------------------------------------------------------
52//FJENDHEADER
53
54#include "fastjet/ClusterSequence.hh"
55#include <iostream> // needed for io
56#include <cstdio> // needed for io
57
58using namespace std;
59
60/// an example program showing how to use fastjet
61int main(){
62
63 // read in input particles
64 //----------------------------------------------------------
65 vector<fastjet::PseudoJet> input_particles;
66
67 double px, py , pz, E;
68 while (cin >> px >> py >> pz >> E) {
69 // create a fastjet::PseudoJet with these components and put it onto
70 // back of the input_particles vector
71 input_particles.push_back(fastjet::PseudoJet(px,py,pz,E));
72 }
73
74
75 // create a jet definition for the kt algorithm (note that one
76 // should not specify an R value here)
77 //----------------------------------------------------------
79
80 // run the jet clustering with the above jet definition
81 //----------------------------------------------------------
82 fastjet::ClusterSequence clust_seq(input_particles, jet_def);
83
84 // get 3 exclusive jets
85 //----------------------------------------------------------
86 int n = 3;
87 vector<fastjet::PseudoJet> exclusive_jets = clust_seq.exclusive_jets(n);
88
89
90 // tell the user what was done
91 // - the description of the algorithm used
92 // - extract the inclusive jets with pt > 5 GeV
93 // show the output as
94 // {index, rap, phi, pt, number of constituents}
95 //----------------------------------------------------------
96 cout << "Ran " << jet_def.description() << endl;
97
98 // label the columns
99 printf("%5s %15s\n","jet #", "E");
100
101 // print out the details for each jet
102 for (unsigned int i = 0; i < exclusive_jets.size(); i++) {
103 printf("%5u %15.8f\n",
104 i, exclusive_jets[i].perp());
105 }
106
107 return 0;
108}
int main()
an example program showing how to use fastjet
deals with clustering
std::vector< PseudoJet > exclusive_jets(const double dcut) const
return a vector of all jets (in the sense of the exclusive algorithm) that would be obtained when run...
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
@ ee_kt_algorithm
the e+e- kt algorithm