FastJet  3.1.3
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Modules Pages
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 //STARTHEADER
19 // $Id: 04-constituents.cc 2704 2011-11-16 11:11:10Z soyez $
20 //
21 // Copyright (c) 2005-2011, 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 and are described in hep-ph/0512210. If you use
33 // FastJet as part of work towards a scientific publication, please
34 // include a citation to the FastJet paper.
35 //
36 // FastJet is distributed in the hope that it will be useful,
37 // but WITHOUT ANY WARRANTY; without even the implied warranty of
38 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
39 // GNU General Public License for more details.
40 //
41 // You should have received a copy of the GNU General Public License
42 // along with FastJet. If not, see <http://www.gnu.org/licenses/>.
43 //----------------------------------------------------------------------
44 //ENDHEADER
45 
46 #include "fastjet/ClusterSequence.hh"
47 #include <iostream> // needed for io
48 #include <cstdio> // needed for io
49 
50 using namespace std;
51 
52 /// an example program showing how to use fastjet
53 int main(){
54 
55  // read in input particles
56  //----------------------------------------------------------
57  vector<fastjet::PseudoJet> input_particles;
58 
59  valarray<double> fourvec(4);
60  int index=0;
61  while (cin >> fourvec[0] >> fourvec[1] >> fourvec[2] >> fourvec[3]) {
62  // create a particle with the approprite 4-momentum and
63  // set its user index to keep track of its index.
64  // you can construct a PseudoJet from any object that allows subscripts
65  // from [0] .. [3] (the last one must be the energy)
66  fastjet::PseudoJet particle(fourvec);
67 
68  particle.set_user_index(index);
69  input_particles.push_back(particle);
70 
71  index++;
72  }
73 
74 
75  // create a jet definition:
76  // a jet algorithm with a given radius parameter
77  //----------------------------------------------------------
78  double R = 0.6;
80 
81 
82  // run the jet clustering with the above jet definition
83  //----------------------------------------------------------
84  fastjet::ClusterSequence clust_seq(input_particles, jet_def);
85 
86 
87  // get the resulting jets ordered in pt
88  //----------------------------------------------------------
89  double ptmin = 5.0;
90  vector<fastjet::PseudoJet> inclusive_jets = sorted_by_pt(clust_seq.inclusive_jets(ptmin));
91 
92 
93  // tell the user what was done
94  // - the description of the algorithm used
95  // - extract the inclusive jets with pt > 5 GeV
96  // show the output as
97  // {index, rap, phi, pt, number of constituents}
98  //----------------------------------------------------------
99  cout << "Ran " << jet_def.description() << endl << endl;
100 
101  // label the columns
102  printf("%5s %15s %15s %15s %15s\n","jet #", "rapidity", "phi", "pt", "n constituents");
103  printf(" indices of constituents\n\n");
104 
105  // print out the details for each jet
106  for (unsigned int i = 0; i < inclusive_jets.size(); i++) {
107  // get the constituents of the jet
108  vector<fastjet::PseudoJet> constituents = inclusive_jets[i].constituents();
109 
110  printf("%5u %15.8f %15.8f %15.8f %8u\n",
111  i, inclusive_jets[i].rap(), inclusive_jets[i].phi(),
112  inclusive_jets[i].perp(), (unsigned int) constituents.size());
113 
114  printf(" ");
115  for (unsigned int j=0; j<constituents.size(); j++){
116  printf("%4u ", constituents[j].user_index());
117  if (j%10==9) printf("\n ");
118  }
119  printf("\n\n");
120  }
121 
122  return 0;
123 }
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:334
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.
the longitudinally invariant kt algorithm
std::string description() const
return a textual description of the current jet definition
int main()
an example program showing how to use fastjet
Class to contain pseudojets, including minimal information of use to jet-clustering routines...
Definition: PseudoJet.hh:67
class that is intended to hold a full definition of the jet clusterer