FastJet 3.4.2
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//FJSTARTHEADER
13// $Id$
14//
15// Copyright (c) 2005-2023, 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. They are described in the original FastJet paper,
27// hep-ph/0512210 and in the manual, arXiv:1111.6097. If you use
28// FastJet as part of work towards a scientific publication, please
29// quote the version you use and include a citation to the manual and
30// optionally also to hep-ph/0512210.
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//FJENDHEADER
41
42#include "fastjet/ClusterSequence.hh"
43#include "fastjet/Selector.hh"
44#include <iostream> // needed for io
45#include <cstdio> // needed for io
46
47using namespace std;
48using namespace fastjet;
49
50int main(){
51
52 // read in input particles
53 //----------------------------------------------------------
54 vector<PseudoJet> input_particles;
55
56 double px, py , pz, E;
57 while (cin >> px >> py >> pz >> E) {
58 // create a PseudoJet with these components and put it onto
59 // the back of the input_particles vector
60 input_particles.push_back(PseudoJet(px,py,pz,E));
61 }
62
63 // Selector application #1: keep particles within a given acceptance
64 // e.g. all particles with 1<|y|<2.5, and all particles with pt>1 for |y|<1
65 // (we include the (redundant) fastjet:: prefix just as a reminder that this
66 // is the namespace where Selector is to be found).
67 //----------------------------------------------------------
68 fastjet::Selector particle_selector = fastjet::SelectorAbsRapRange(1.0,2.5)
70 cout << input_particles.size() << " particles before selector" << endl;
71 input_particles = particle_selector(input_particles);
72 cout << input_particles.size() << " particles after selector" << endl;
73
74 // create a jet definition:
75 // a jet algorithm with a given radius parameter
76 //----------------------------------------------------------
77 double R = 0.6;
78 JetDefinition jet_def(kt_algorithm, R);
79
80
81 // run the jet clustering with the above jet definition
82 //----------------------------------------------------------
83 ClusterSequence clust_seq(input_particles, jet_def);
84
85
86 // get the 5 hardest jets within |y|<2
87 //
88 // Note that this nicely illustrates that you should watch out that
89 // Selectors do not necessarily commute.
90 //
91 // The && operator behaves like a logical and, i.e. it keeps objects
92 // that satisfy both criteria (independently). It does commute.
93 //
94 // The * operator applies Selectors successively (starting from the
95 // rightmost one as in a usual operator product). Here, order may matter.
96 //----------------------------------------------------------
97 Selector jet_selector = SelectorNHardest(5) * SelectorAbsRapMax(2.0);
98 vector<PseudoJet> inclusive_jets = sorted_by_pt(jet_selector(clust_seq.inclusive_jets()));
99
100
101 // tell the user what was done
102 // - the description of the algorithm used
103 // - the description of teh selectors used
104 // - the jets
105 // show the output as
106 // {index, rap, phi, pt}
107 //----------------------------------------------------------
108 double rapmin, rapmax;
109
110 cout << "Ran " << jet_def.description() << endl;
111
112 cout << "Selected particles: " << particle_selector.description() << endl;
113 particle_selector.get_rapidity_extent(rapmin, rapmax);
114 cout << " with a total rapidity range of [" << rapmin << ", " << rapmax << "]" << endl;
115
116 cout << "Selected jets: " << jet_selector.description() << endl;
117 jet_selector.get_rapidity_extent(rapmin, rapmax);
118 cout << " with a total rapidity range of [" << rapmin << ", " << rapmax << "]" << endl;
119
120 // label the columns
121 printf("%5s %15s %15s %15s\n","jet #", "rapidity", "phi", "pt");
122
123 // print out the details for each jet
124 for (unsigned int i = 0; i < inclusive_jets.size(); i++) {
125 printf("%5u %15.8f %15.8f %15.8f\n",
126 i, inclusive_jets[i].rap(), inclusive_jets[i].phi(),
127 inclusive_jets[i].perp());
128 }
129
130 return 0;
131}
int main()
an example program showing how to use fastjet
Definition: 01-basic.cc:52
deals with clustering
class that is intended to hold a full definition of the jet clusterer
Class to contain pseudojets, including minimal information of use to jet-clustering routines.
Definition: PseudoJet.hh:68
Class that encodes information about cuts and other selection criteria that can be applied to PseudoJ...
Definition: Selector.hh:149
std::string description() const
returns a textual description of the selector
Definition: Selector.hh:237
void get_rapidity_extent(double &rapmin, double &rapmax) const
returns the rapidity range for which it may return "true"
Definition: Selector.hh:232
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
the FastJet namespace
vector< PseudoJet > sorted_by_pt(const vector< PseudoJet > &jets)
return a vector of jets sorted into decreasing kt2
Definition: PseudoJet.cc:871