FastJet  3.1.3
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Modules Pages
02-jetdef.cc
Go to the documentation of this file.
1 //----------------------------------------------------------------------
2 /// \file
3 /// \page Example02 02 - changing the jet definition
4 ///
5 /// fastjet basic example program:
6 /// illustration of the usage of how to change the jet definition
7 /// used for the clustering (see also fastjet::JetDefinition)
8 ///
9 /// run it with : ./02-jetdef < data/single-event.dat
10 ///
11 /// Source code: 02-jetdef.cc
12 //----------------------------------------------------------------------
13 //
14 //STARTHEADER
15 // $Id: 02-jetdef.cc 3663 2014-09-05 07:23:19Z soyez $
16 //
17 // Copyright (c) 2005-2011, Matteo Cacciari, Gavin P. Salam and Gregory Soyez
18 //
19 //----------------------------------------------------------------------
20 // This file is part of FastJet.
21 //
22 // FastJet is free software; you can redistribute it and/or modify
23 // it under the terms of the GNU General Public License as published by
24 // the Free Software Foundation; either version 2 of the License, or
25 // (at your option) any later version.
26 //
27 // The algorithms that underlie FastJet have required considerable
28 // development and are described in hep-ph/0512210. If you use
29 // FastJet as part of work towards a scientific publication, please
30 // include a citation to the FastJet paper.
31 //
32 // FastJet is distributed in the hope that it will be useful,
33 // but WITHOUT ANY WARRANTY; without even the implied warranty of
34 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
35 // GNU General Public License for more details.
36 //
37 // You should have received a copy of the GNU General Public License
38 // along with FastJet. If not, see <http://www.gnu.org/licenses/>.
39 //----------------------------------------------------------------------
40 //ENDHEADER
41 
42 
43 #include "fastjet/ClusterSequence.hh"
44 #include <iostream> // needed for io
45 #include <cstdio> // needed for io
46 
47 using namespace std;
48 
49 /// an example program showing how to use fastjet
50 int main(){
51 
52  // read in input particles
53  //----------------------------------------------------------
54  vector<fastjet::PseudoJet> input_particles;
55 
56  double px, py , pz, E;
57  while (cin >> px >> py >> pz >> E) {
58  // create a fastjet::PseudoJet with these components and put it onto
59  // back of the input_particles vector
60  input_particles.push_back(fastjet::PseudoJet(px,py,pz,E));
61  }
62 
63 
64  // create a jet definition:
65  // a jet algorithm with a given radius parameter
66  //----------------------------------------------------------
67 
68  // select a jet algorithm to use
69  //
70  // this could be one of
71  // {kt_algorithm, cambridge_algorithm, antikt_algorithm,
72  // genkt_algorithm, ee_kt_algorithm, ee_genkt_algorithm}
73  // see example 03-plugin.cc for extra options using plugins
74  // instead of naive algorithms)
76 
77  // when appropriate select a radius to use
78  //
79  // this would not be mandatory for e^+ e^- algorithms (see 05-eplus_eminus.cc)
80  double R = 0.6;
81 
82  // select an __optional__ recombination scheme
83  //
84  // this could be one of
85  // {E_scheme, pt_scheme, pt2_scheme, Et_scheme, Et2_scheme, BIpt_scheme,
86  // BIpt2_scheme, WTA_pt_scheme, WTA_E_scheme, WTA_modp_scheme,
87  // external_sheme}
88  //
89  // Notes:
90  // - for the usage of a user-defined recombination scheme
91  // (external_scheme), see 11-boosted_higgs.cc
92  // - WTA_E_scheme, WTA_modp_scheme are meant for e+e- clusterings
93  //
94  // By default, the E_scheme is used
96 
97  // select an __optional__ strategy
98  //
99  // this could be chosen among
100  // {N2MinHeapTiled, N2Tiled, N2PoorTiled, N2Plain, N3Dumb,
101  // Best,
102  // NlnN, NlnN3pi, NlnN4pi, NlnNCam4pi, NlnNCam2pi2R, NlnNCam}
103  //
104  // By default, the Best strategy is chosen and we advise to keep
105  // that default unless you are targeting a very specific usage. Note
106  // also that the N log (N) strategies for algorithms other than
107  // Cambridge/Aachen need CGAL support.
108  fastjet::Strategy strategy = fastjet::Best;
109 
110  // create the JetDefinition from the above information
111  fastjet::JetDefinition jet_def(jet_alg, R, recomb_scheme, strategy);
112 
113 
114  // run the jet clustering with the above jet definition
115  //----------------------------------------------------------
116  fastjet::ClusterSequence clust_seq(input_particles, jet_def);
117 
118 
119  // get the resulting jets ordered in pt
120  //----------------------------------------------------------
121  double ptmin = 5.0;
122  vector<fastjet::PseudoJet> inclusive_jets = sorted_by_pt(clust_seq.inclusive_jets(ptmin));
123 
124 
125  // tell the user what was done
126  // - the description of the algorithm used
127  // - extract the inclusive jets with pt > 5 GeV
128  // show the output as
129  // {index, rap, phi, pt}
130  //----------------------------------------------------------
131  cout << "Ran " << jet_def.description() << endl;
132 
133  // label the columns
134  printf("%5s %15s %15s %15s\n","jet #", "rapidity", "phi", "pt");
135 
136  // print out the details for each jet
137  for (unsigned int i = 0; i < inclusive_jets.size(); i++) {
138  printf("%5u %15.8f %15.8f %15.8f\n",
139  i, inclusive_jets[i].rap(), inclusive_jets[i].phi(),
140  inclusive_jets[i].perp());
141  }
142 
143  return 0;
144 }
summing the 4-momenta
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
RecombinationScheme
The various recombination schemes.
automatic selection of the best (based on N), including the LazyTiled strategies that are new to FJ3...
JetAlgorithm
the various families of jet-clustering algorithm
Strategy
the various options for the algorithmic strategy to adopt in clustering events with kt and cambridge ...
int main()
an example program showing how to use fastjet
Definition: 02-jetdef.cc:50
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