FastJet  3.1.0-beta.1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
fastjet_example.cc
Go to the documentation of this file.
1 /// \file
2 /// \page Examples FastJet examples
3 ///
4 /// The FastJet examples have been organised by order of complexity,
5 /// starting by the simplest case and introducing features one after
6 /// another.
7 /// - \subpage Example01
8 /// - \subpage Example02
9 /// - \subpage Example03
10 /// - \subpage Example04
11 /// - \subpage Example05
12 /// - \subpage Example06
13 /// - \subpage Example07 (\subpage Example07old "old version")
14 /// - \subpage Example08
15 /// - \subpage Example09
16 /// - \subpage Example10
17 /// - \subpage Example11
18 /// - \subpage Example12 (\subpage Example12old "old version")
19 /// - \subpage Example13
20 /// - \subpage Example14
21 
22 //STARTHEADER
23 // $Id: fastjet_example.cc 2684 2011-11-14 07:41:44Z soyez $
24 //
25 // Copyright (c) 2005-2011, 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 
51 //----------------------------------------------------------------------
52 // fastjet example program.
53 //
54 // Compile it with: make fastjet_example
55 // run it with : ./fastjet_example < data/single-event.dat
56 //
57 // People who are familiar with the ktjet package are encouraged to
58 // compare this file to the ktjet_example.cc program which does the
59 // same thing in the ktjet framework.
60 //----------------------------------------------------------------------
61 #include "fastjet/PseudoJet.hh"
62 #include "fastjet/ClusterSequence.hh"
63 #include<iostream> // needed for io
64 #include<sstream> // needed for internal io
65 #include<vector>
66 #include <cstdio>
67 
68 using namespace std;
69 
70 // a declaration of a function that pretty prints a list of jets
71 void print_jets (const vector<fastjet::PseudoJet> &);
72 
73 /// an example program showing how to use fastjet
74 int main () {
75 
76  vector<fastjet::PseudoJet> input_particles;
77 
78  // Read in input particles
79  double px, py , pz, E;
80  while (cin >> px >> py >> pz >> E) {
81  // create a fastjet::PseudoJet with these components and put it onto
82  // back of the input_particles vector
83  input_particles.push_back(fastjet::PseudoJet(px,py,pz,E));
84  }
85 
86  // create an object that represents your choice of jet algorithm and
87  // the associated parameters
88  double Rparam = 1.0;
91  fastjet::JetDefinition jet_def(fastjet::kt_algorithm, Rparam, recomb_scheme, strategy);
92 
93  // run the jet clustering with the above jet definition
94  fastjet::ClusterSequence clust_seq(input_particles, jet_def);
95 
96  // tell the user what was done
97  cout << "Ran " << jet_def.description() << endl;
98  cout << "Strategy adopted by FastJet was "<<
99  clust_seq.strategy_string()<<endl<<endl;
100 
101  // extract the inclusive jets with pt > 5 GeV
102  double ptmin = 5.0;
103  vector<fastjet::PseudoJet> inclusive_jets = clust_seq.inclusive_jets(ptmin);
104 
105  // print them out
106  cout << "Printing inclusive jets with pt > "<< ptmin<<" GeV\n";
107  cout << "---------------------------------------\n";
108  print_jets(inclusive_jets);
109  cout << endl;
110 
111  // extract the exclusive jets with dcut = 25 GeV^2
112  double dcut = 25.0;
113  vector<fastjet::PseudoJet> exclusive_jets = clust_seq.exclusive_jets(dcut);
114 
115  // print them out
116  cout << "Printing exclusive jets with dcut = "<< dcut<<" GeV^2\n";
117  cout << "--------------------------------------------\n";
118  print_jets(exclusive_jets);
119 
120 
121 }
122 
123 
124 //----------------------------------------------------------------------
125 /// a function that pretty prints a list of jets
126 void print_jets (const vector<fastjet::PseudoJet> & jets) {
127 
128  // sort jets into increasing pt
129  vector<fastjet::PseudoJet> sorted_jets = sorted_by_pt(jets);
130 
131  // label the columns
132  printf("%5s %15s %15s %15s %15s\n","jet #", "rapidity",
133  "phi", "pt", "n constituents");
134 
135  // print out the details for each jet
136  for (unsigned int i = 0; i < sorted_jets.size(); i++) {
137  // the following is not super efficient since it creates an
138  // intermediate constituents vector
139  int n_constituents = sorted_jets[i].constituents().size();
140  printf("%5u %15.8f %15.8f %15.8f %8u\n",
141  i, sorted_jets[i].rap(), sorted_jets[i].phi(),
142  sorted_jets[i].perp(), n_constituents);
143  }
144 
145 }