FastJet 3.4.3
Loading...
Searching...
No Matches
04-constituents.cc
Go to the documentation of this file.
1//----------------------------------------------------------------------
2/// \file
3/// \page Example04 04 - accessing clustering information in a PseudoJet
4///
5/// illustrate how a jet can carry information about its clustering
6///
7/// We do it by associating a user index to each of the input particles
8/// and show what particles are in each jets (above 5 GeV)
9///
10/// We also illustrate a few other features about how a fastjet::PseudoJet
11/// can access its underlying structure.
12///
13/// run it with : ./04-constituents < data/single-event.dat
14///
15/// Source code: 04-constituents.cc
16//----------------------------------------------------------------------
17
18//FJSTARTHEADER
19// $Id$
20//
21// Copyright (c) 2005-2024, Matteo Cacciari, Gavin P. Salam and Gregory Soyez
22//
23//----------------------------------------------------------------------
24// This file is part of FastJet.
25//
26// FastJet is free software; you can redistribute it and/or modify
27// it under the terms of the GNU General Public License as published by
28// the Free Software Foundation; either version 2 of the License, or
29// (at your option) any later version.
30//
31// The algorithms that underlie FastJet have required considerable
32// development. They are described in the original FastJet paper,
33// hep-ph/0512210 and in the manual, arXiv:1111.6097. If you use
34// FastJet as part of work towards a scientific publication, please
35// quote the version you use and include a citation to the manual and
36// optionally also to hep-ph/0512210.
37//
38// FastJet is distributed in the hope that it will be useful,
39// but WITHOUT ANY WARRANTY; without even the implied warranty of
40// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
41// GNU General Public License for more details.
42//
43// You should have received a copy of the GNU General Public License
44// along with FastJet. If not, see <http://www.gnu.org/licenses/>.
45//----------------------------------------------------------------------
46//FJENDHEADER
47
48#include "fastjet/ClusterSequence.hh"
49#include <iostream> // needed for io
50#include <cstdio> // needed for io
51
52using namespace std;
53
54/// an example program showing how to use fastjet
55int main(){
56
57 // read in input particles
58 //----------------------------------------------------------
59 vector<fastjet::PseudoJet> input_particles;
60
61 valarray<double> fourvec(4);
62 int index=0;
63 while (cin >> fourvec[0] >> fourvec[1] >> fourvec[2] >> fourvec[3]) {
64 // create a particle with the approprite 4-momentum and
65 // set its user index to keep track of its index.
66 // you can construct a PseudoJet from any object that allows subscripts
67 // from [0] .. [3] (the last one must be the energy)
68 fastjet::PseudoJet particle(fourvec);
69
70 particle.set_user_index(index);
71 input_particles.push_back(particle);
72
73 index++;
74 }
75
76
77 // create a jet definition:
78 // a jet algorithm with a given radius parameter
79 //----------------------------------------------------------
80 double R = 0.6;
82
83
84 // run the jet clustering with the above jet definition
85 //----------------------------------------------------------
86 fastjet::ClusterSequence clust_seq(input_particles, jet_def);
87
88
89 // get the resulting jets ordered in pt
90 //----------------------------------------------------------
91 double ptmin = 5.0;
92 vector<fastjet::PseudoJet> inclusive_jets = sorted_by_pt(clust_seq.inclusive_jets(ptmin));
93
94
95 // tell the user what was done
96 // - the description of the algorithm used
97 // - extract the inclusive jets with pt > 5 GeV
98 // show the output as
99 // {index, rap, phi, pt, number of constituents}
100 //----------------------------------------------------------
101 cout << "Ran " << jet_def.description() << endl << endl;
102
103 // label the columns
104 printf("%5s %15s %15s %15s %15s\n","jet #", "rapidity", "phi", "pt", "n constituents");
105 printf(" indices of constituents\n\n");
106
107 // print out the details for each jet
108 for (unsigned int i = 0; i < inclusive_jets.size(); i++) {
109 // get the constituents of the jet
110 vector<fastjet::PseudoJet> constituents = inclusive_jets[i].constituents();
111
112 printf("%5u %15.8f %15.8f %15.8f %8u\n",
113 i, inclusive_jets[i].rap(), inclusive_jets[i].phi(),
114 inclusive_jets[i].perp(), (unsigned int) constituents.size());
115
116 printf(" ");
117 for (unsigned int j=0; j<constituents.size(); j++){
118 printf("%4u ", constituents[j].user_index());
119 if (j%10==9) printf("\n ");
120 }
121 printf("\n\n");
122 }
123
124 return 0;
125}
int main()
an example program showing how to use fastjet
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
void set_user_index(const int index)
set the user_index, intended to allow the user to add simple identifying information to a particle/je...
Definition PseudoJet.hh:392
@ kt_algorithm
the longitudinally invariant kt algorithm