FastJet  3.1.3
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Modules Pages
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 //STARTHEADER
25 // $Id: 05-eplus_eminus.cc 2684 2011-11-14 07:41:44Z soyez $
26 //
27 // Copyright (c) 2005-2011, 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 and are described in hep-ph/0512210. If you use
39 // FastJet as part of work towards a scientific publication, please
40 // include a citation to the FastJet paper.
41 //
42 // FastJet is distributed in the hope that it will be useful,
43 // but WITHOUT ANY WARRANTY; without even the implied warranty of
44 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
45 // GNU General Public License for more details.
46 //
47 // You should have received a copy of the GNU General Public License
48 // along with FastJet. If not, see <http://www.gnu.org/licenses/>.
49 //----------------------------------------------------------------------
50 //ENDHEADER
51 
52 #include "fastjet/ClusterSequence.hh"
53 #include <iostream> // needed for io
54 #include <cstdio> // needed for io
55 
56 using namespace std;
57 
58 /// an example program showing how to use fastjet
59 int main(){
60 
61  // read in input particles
62  //----------------------------------------------------------
63  vector<fastjet::PseudoJet> input_particles;
64 
65  double px, py , pz, E;
66  while (cin >> px >> py >> pz >> E) {
67  // create a fastjet::PseudoJet with these components and put it onto
68  // back of the input_particles vector
69  input_particles.push_back(fastjet::PseudoJet(px,py,pz,E));
70  }
71 
72 
73  // create a jet definition for the kt algorithm (note that one
74  // should not specify an R value here)
75  //----------------------------------------------------------
77 
78  // run the jet clustering with the above jet definition
79  //----------------------------------------------------------
80  fastjet::ClusterSequence clust_seq(input_particles, jet_def);
81 
82  // get 3 exclusive jets
83  //----------------------------------------------------------
84  int n = 3;
85  vector<fastjet::PseudoJet> exclusive_jets = clust_seq.exclusive_jets(n);
86 
87 
88  // tell the user what was done
89  // - the description of the algorithm used
90  // - extract the inclusive jets with pt > 5 GeV
91  // show the output as
92  // {index, rap, phi, pt, number of constituents}
93  //----------------------------------------------------------
94  cout << "Ran " << jet_def.description() << endl;
95 
96  // label the columns
97  printf("%5s %15s\n","jet #", "E");
98 
99  // print out the details for each jet
100  for (unsigned int i = 0; i < exclusive_jets.size(); i++) {
101  printf("%5u %15.8f\n",
102  i, exclusive_jets[i].perp());
103  }
104 
105  return 0;
106 }
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...
deals with clustering
std::string description() const
return a textual description of the current jet definition
int main()
an example program showing how to use fastjet
the e+e- kt algorithm
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