FastJet  3.1.3
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Modules Pages
08-selector.cc
Go to the documentation of this file.
1 //----------------------------------------------------------------------
2 /// \file
3 /// \page Example08 08 - using the Selector tool
4 ///
5 /// fastjet sample program to illustrate the use of fastjet::Selector
6 ///
7 /// run it with : ./08-selector < data/single-event.dat
8 ///
9 /// Source code: 08-selector.cc
10 //----------------------------------------------------------------------
11 
12 //STARTHEADER
13 // $Id: 08-selector.cc 2684 2011-11-14 07:41:44Z soyez $
14 //
15 // Copyright (c) 2005-2011, Matteo Cacciari, Gavin P. Salam and Gregory Soyez
16 //
17 //----------------------------------------------------------------------
18 // This file is part of FastJet.
19 //
20 // FastJet is free software; you can redistribute it and/or modify
21 // it under the terms of the GNU General Public License as published by
22 // the Free Software Foundation; either version 2 of the License, or
23 // (at your option) any later version.
24 //
25 // The algorithms that underlie FastJet have required considerable
26 // development and are described in hep-ph/0512210. If you use
27 // FastJet as part of work towards a scientific publication, please
28 // include a citation to the FastJet paper.
29 //
30 // FastJet is distributed in the hope that it will be useful,
31 // but WITHOUT ANY WARRANTY; without even the implied warranty of
32 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
33 // GNU General Public License for more details.
34 //
35 // You should have received a copy of the GNU General Public License
36 // along with FastJet. If not, see <http://www.gnu.org/licenses/>.
37 //----------------------------------------------------------------------
38 //ENDHEADER
39 
40 #include "fastjet/ClusterSequence.hh"
41 #include "fastjet/Selector.hh"
42 #include <iostream> // needed for io
43 #include <cstdio> // needed for io
44 
45 using namespace std;
46 using namespace fastjet;
47 
48 int main(){
49 
50  // read in input particles
51  //----------------------------------------------------------
52  vector<PseudoJet> input_particles;
53 
54  double px, py , pz, E;
55  while (cin >> px >> py >> pz >> E) {
56  // create a PseudoJet with these components and put it onto
57  // the back of the input_particles vector
58  input_particles.push_back(PseudoJet(px,py,pz,E));
59  }
60 
61  // Selector application #1: keep particles within a given acceptance
62  // e.g. all particles with 1<|y|<2.5, and all particles with pt>1 for |y|<1
63  // (we include the (redundant) fastjet:: prefix just as a reminder that this
64  // is the namespace where Selector is to be found).
65  //----------------------------------------------------------
66  fastjet::Selector particle_selector = fastjet::SelectorAbsRapRange(1.0,2.5)
68  cout << input_particles.size() << " particles before selector" << endl;
69  input_particles = particle_selector(input_particles);
70  cout << input_particles.size() << " particles after selector" << endl;
71 
72  // create a jet definition:
73  // a jet algorithm with a given radius parameter
74  //----------------------------------------------------------
75  double R = 0.6;
76  JetDefinition jet_def(kt_algorithm, R);
77 
78 
79  // run the jet clustering with the above jet definition
80  //----------------------------------------------------------
81  ClusterSequence clust_seq(input_particles, jet_def);
82 
83 
84  // get the 5 hardest jets within |y|<2
85  //
86  // Note that this nicely illustrates that you should watch out that
87  // Selectors do not necessarily commute.
88  //
89  // The && operator behaves like a logical and, i.e. it keeps objects
90  // that satisfy both criteria (independently). It does commute.
91  //
92  // The * operator applies Selectors successively (starting from the
93  // rightmost one as in a usual operator product). Here, order may matter.
94  //----------------------------------------------------------
95  Selector jet_selector = SelectorNHardest(5) * SelectorAbsRapMax(2.0);
96  vector<PseudoJet> inclusive_jets = sorted_by_pt(jet_selector(clust_seq.inclusive_jets()));
97 
98 
99  // tell the user what was done
100  // - the description of the algorithm used
101  // - the description of teh selectors used
102  // - the jets
103  // show the output as
104  // {index, rap, phi, pt}
105  //----------------------------------------------------------
106  double rapmin, rapmax;
107 
108  cout << "Ran " << jet_def.description() << endl;
109 
110  cout << "Selected particles: " << particle_selector.description() << endl;
111  particle_selector.get_rapidity_extent(rapmin, rapmax);
112  cout << " with a total rapidity range of [" << rapmin << ", " << rapmax << "]" << endl;
113 
114  cout << "Selected jets: " << jet_selector.description() << endl;
115  jet_selector.get_rapidity_extent(rapmin, rapmax);
116  cout << " with a total rapidity range of [" << rapmin << ", " << rapmax << "]" << endl;
117 
118  // label the columns
119  printf("%5s %15s %15s %15s\n","jet #", "rapidity", "phi", "pt");
120 
121  // print out the details for each jet
122  for (unsigned int i = 0; i < inclusive_jets.size(); i++) {
123  printf("%5u %15.8f %15.8f %15.8f\n",
124  i, inclusive_jets[i].rap(), inclusive_jets[i].phi(),
125  inclusive_jets[i].perp());
126  }
127 
128  return 0;
129 }
vector< PseudoJet > sorted_by_pt(const vector< PseudoJet > &jets)
return a vector of jets sorted into decreasing kt2
Definition: PseudoJet.cc:799
void get_rapidity_extent(double &rapmin, double &rapmax) const
returns the rapidity range for which it may return "true"
Definition: Selector.hh:232
deals with clustering
Selector SelectorNHardest(unsigned int n)
select the n hardest objects
Definition: Selector.cc:1074
Selector SelectorAbsRapMax(double absrapmax)
select objects with |rap| <= absrapmax
Definition: Selector.cc:880
Selector SelectorPtMin(double ptmin)
select objects with pt >= ptmin
Definition: Selector.cc:681
Selector SelectorAbsRapRange(double rapmin, double rapmax)
select objects with absrapmin <= |rap| <= absrapmax
Definition: Selector.cc:885
int main()
an example program showing how to use fastjet
Definition: 01-basic.cc:50
the FastJet namespace
std::string description() const
returns a textual description of the selector
Definition: Selector.hh:237
Class that encodes information about cuts and other selection criteria that can be applied to PseudoJ...
Definition: Selector.hh:149
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