FastJet 3.4.1
03-plugin.cc
Go to the documentation of this file.
1//----------------------------------------------------------------------
2/// \file
3/// \page Example03 03 - using plugins
4///
5/// fastjet plugins example program:
6/// we illustrate the plugin usage
7/// here, we use the SISCone plugin though different choices are possible
8/// see the output of 'fastjet-config --list-plugins' for more details
9///
10/// Note that when using plugins, the code needs to be linked against
11/// the libfastjetplugins library (with the default monolithic
12/// build. For non-monolithic build, individual libraries have to be
13/// used for each plugin).
14/// This is ensured in practice by calling
15/// fastjet-config --libs --plugins
16///
17/// run it with : ./03-plugin < data/single-event.dat
18///
19/// Source code: 03-plugin.cc
20//----------------------------------------------------------------------
21
22//STARTHEADER
23// $Id$
24//
25// Copyright (c) 2005-2018, Matteo Cacciari, Gavin P. Salam and Gregory Soyez
26//
27//----------------------------------------------------------------------
28// This file is part of FastJet.
29//
30// FastJet is free software; you can redistribute it and/or modify
31// it under the terms of the GNU General Public License as published by
32// the Free Software Foundation; either version 2 of the License, or
33// (at your option) any later version.
34//
35// The algorithms that underlie FastJet have required considerable
36// development and are described in hep-ph/0512210. If you use
37// FastJet as part of work towards a scientific publication, please
38// include a citation to the FastJet paper.
39//
40// FastJet is distributed in the hope that it will be useful,
41// but WITHOUT ANY WARRANTY; without even the implied warranty of
42// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
43// GNU General Public License for more details.
44//
45// You should have received a copy of the GNU General Public License
46// along with FastJet. If not, see <http://www.gnu.org/licenses/>.
47//----------------------------------------------------------------------
48//ENDHEADER
49
50#include "fastjet/ClusterSequence.hh"
51#include <iostream> // needed for io
52#include <cstdio> // needed for io
53
54// include the SISCone plugin header if enabled
55#include "fastjet/config.h"
56#ifdef FASTJET_ENABLE_PLUGIN_SISCONE
57#include "fastjet/SISConePlugin.hh"
58#else
59#warning "SISCone plugin not enabled. Skipping the example"
60#endif // FASTJET_ENABLE_PLUGIN_SISCONE
61
62
63using namespace std;
64
65int main(){
66
67#ifdef FASTJET_ENABLE_PLUGIN_SISCONE
68
69 // read in input particles
70 //----------------------------------------------------------
71 vector<fastjet::PseudoJet> input_particles;
72
73 double px, py , pz, E;
74 while (cin >> px >> py >> pz >> E) {
75 // create a fastjet::PseudoJet with these components and put it onto
76 // back of the input_particles vector
77 input_particles.push_back(fastjet::PseudoJet(px,py,pz,E));
78 }
79
80
81 // create a jet definition fron a plugin.
82 // It basically requires declaring a JetDefinition from a pointer to
83 // the plugin
84 //
85 // we will use the SISCone plugin here. Its (mandatory) parameters
86 // are a cone radius and an overlap threshold, plus other optional
87 // parameters
88 //
89 // for other plugin, see individual documentations for a description
90 // of their parameters
91 //
92 // the list of available plugins for a given build of FastJet can be
93 // obtained using
94 // fastjet-config --list-plugins
95 // from the command line.
96 //----------------------------------------------------------
97 double cone_radius = 0.7;
98 double overlap_threshold = 0.75;
99 fastjet::SISConePlugin siscone(cone_radius, overlap_threshold);
100 fastjet::JetDefinition jet_def(& siscone);
101
102
103 // run the jet clustering with the above jet definition
104 //----------------------------------------------------------
105 fastjet::ClusterSequence clust_seq(input_particles, jet_def);
106
107
108 // get the resulting jets ordered in pt
109 //----------------------------------------------------------
110 double ptmin = 5.0;
111 vector<fastjet::PseudoJet> inclusive_jets = sorted_by_pt(clust_seq.inclusive_jets(ptmin));
112
113
114 // tell the user what was done
115 // - the description of the algorithm used
116 // - extract the inclusive jets with pt > 5 GeV
117 // show the output as
118 // {index, rap, phi, pt}
119 //----------------------------------------------------------
120 cout << "Ran " << jet_def.description() << endl;
121
122 // label the columns
123 printf("%5s %15s %15s %15s\n","jet #", "rapidity", "phi", "pt");
124
125 // print out the details for each jet
126 for (unsigned int i = 0; i < inclusive_jets.size(); i++) {
127 printf("%5u %15.8f %15.8f %15.8f\n",
128 i, inclusive_jets[i].rap(), inclusive_jets[i].phi(),
129 inclusive_jets[i].perp());
130 }
131
132#endif
133
134 return 0;
135
136}
int main()
an example program showing how to use fastjet
Definition: 01-basic.cc:50
deals with clustering
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
Implementation of the SISCone algorithm (plugin for fastjet v2.1 upwards)
vector< PseudoJet > sorted_by_pt(const vector< PseudoJet > &jets)
return a vector of jets sorted into decreasing kt2
Definition: PseudoJet.cc:871