FastJet 3.4.2
14-groomers.cc
Go to the documentation of this file.
1//----------------------------------------------------------------------
2/// \file
3/// \page Example14 14 - unified use of transformers
4///
5/// fastjet example program to illustrate the use of the
6/// fastjet::Filter and fastjet::Pruner classes in a unified way
7/// through their derivation from fastjet::Transformer
8///
9/// The two hardest jets in a boosted top event, clustered with an
10/// (abnormally) large R, are then groomed using different tools. One
11/// notes the reduction in the mass of the jets after grooming.
12///
13/// run it with : ./14-groomers < data/boosted_top_event.dat
14///
15/// Source code: 14-groomers.cc
16//----------------------------------------------------------------------
17
18//FJSTARTHEADER
19// $Id$
20//
21// Copyright (c) 2005-2023, Matteo Cacciari, Gavin P. Salam and Gregory Soyez
22//
23//----------------------------------------------------------------------
24// This file is part of FastJet.
25//
26// FastJet is free software; you can redistribute it and/or modify
27// it under the terms of the GNU General Public License as published by
28// the Free Software Foundation; either version 2 of the License, or
29// (at your option) any later version.
30//
31// The algorithms that underlie FastJet have required considerable
32// development. They are described in the original FastJet paper,
33// hep-ph/0512210 and in the manual, arXiv:1111.6097. If you use
34// FastJet as part of work towards a scientific publication, please
35// quote the version you use and include a citation to the manual and
36// optionally also to hep-ph/0512210.
37//
38// FastJet is distributed in the hope that it will be useful,
39// but WITHOUT ANY WARRANTY; without even the implied warranty of
40// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
41// GNU General Public License for more details.
42//
43// You should have received a copy of the GNU General Public License
44// along with FastJet. If not, see <http://www.gnu.org/licenses/>.
45//----------------------------------------------------------------------
46//FJENDHEADER
47
48#include <fastjet/PseudoJet.hh>
49#include <fastjet/ClusterSequence.hh>
50#include <fastjet/Selector.hh>
51#include <iostream>
52#include "fastjet/tools/Filter.hh"
53#include "fastjet/tools/Pruner.hh"
54
55#include <cstdio> // needed for io
56
57using namespace fastjet;
58using namespace std;
59
60/// an example program showing how to use Filter and Pruner in FastJet
61int main(){
62 // read in input particles
63 //----------------------------------------------------------
64 vector<PseudoJet> input_particles;
65
66 double px, py , pz, E;
67 while (cin >> px >> py >> pz >> E) {
68 // create a PseudoJet with these components and put it onto
69 // back of the input_particles vector
70 input_particles.push_back(PseudoJet(px,py,pz,E));
71 }
72
73 // get the resulting jets ordered in pt
74 //----------------------------------------------------------
76 ClusterSequence clust_seq(input_particles, jet_def);
77 vector<PseudoJet> inclusive_jets =
78 sorted_by_pt(clust_seq.inclusive_jets(5.0));
79
80 // label the columns
81 printf("%5s %15s %15s %15s %15s\n","jet #", "rapidity", "phi", "pt", "mass");
82
83 // print out the details for each jet
84 for (unsigned int i = 0; i < inclusive_jets.size(); i++) {
85 printf("%5u %15.8f %15.8f %15.8f %15.8f\n",
86 i, inclusive_jets[i].rap(), inclusive_jets[i].phi(),
87 inclusive_jets[i].perp(),inclusive_jets[i].m());
88 }
89
90 // simple test to avoid that the example below crashes:
91 // make sure there is at least 3 jets above our 5 GeV
92 if (inclusive_jets.size()<2){
93 cout << "Please provide an event with at least 2 jets above 5 GeV" << endl;
94 return 1;
95 }
96
97 // We will groom the two hardest jets of the event
98 //----------------------------------------------------------
99 vector<PseudoJet> candidates = SelectorNHardest(2)(inclusive_jets);
100
101 // create 3 groomers
102 //----------------------------------------------------------
103 vector<Transformer *> groomers;
104
105 // 1.
106 // the Cambridge/Aachen filter with Rfilt=0.3
107 // (simplified version of arXiv:0802.2470)
108 double Rfilt = 0.3;
109 unsigned int nfilt = 3;
110 groomers.push_back(new Filter(JetDefinition(cambridge_algorithm, Rfilt),
111 SelectorNHardest(nfilt) ) );
112
113 // 2.
114 // Filtering with a pt cut as for trimming (arXiv:0912.1342)
115 double Rtrim = 0.2;
116 double ptfrac = 0.03;
117 groomers.push_back(new Filter(JetDefinition(kt_algorithm, Rtrim),
118 SelectorPtFractionMin(ptfrac) ) );
119
120 // 3.
121 // Pruning (arXiv:0903.5081)
122 double zcut = 0.1;
123 double rcut_factor = 0.5;
124 groomers.push_back(new Pruner(cambridge_algorithm, zcut, rcut_factor));
125
126 // apply the various groomers to the test PseudoJet's
127 // and show the result
128 //----------------------------------------------------------
129
130 // print out original jet candidates
131 cout << "\nOriginal jets that will be grooomed: " << endl;
132 for (vector<PseudoJet>::iterator jit=candidates.begin(); jit!=candidates.end(); jit++){
133 const PseudoJet & c = *jit;
134 cout << " rap = " << c.rap() << ", phi = " << c.phi() << ", pt = " << c.perp()
135 << ", mass = " << c.m()
136 << " [" << c.description() << "]" << endl;
137 }
138
139 // loop on groomers
140 for (unsigned int i=0; i < groomers.size(); i++){
141 const Transformer & f = *groomers[i];
142 cout << "\nUsing groomer: " << f.description() << endl;
143
144 // loop on jet candidates
145 for (vector<PseudoJet>::iterator jit=candidates.begin(); jit!=candidates.end(); jit++){
146 const PseudoJet & c = *jit;
147
148 // apply groomer f to jet c
149 PseudoJet j = f(c);
150
151 // access properties specific to the given transformer
152 //
153 // We first make sure that the jet indeed has a structure
154 // compatible with the result of a Filter or Pruner (using
155 // has_structure_of()), and then retrieve the pieces rejected by the
156 // groomer (using structure_of())
157 int n_rejected;
158 if (j.has_structure_of<Filter>()) {
159 const Filter::StructureType & fj_struct = j.structure_of<Filter>();
160 n_rejected = fj_struct.rejected().size();
161 }else {
162 assert(j.has_structure_of<Pruner>()); // make sure
163 const Pruner::StructureType & fj_struct = j.structure_of<Pruner>();
164 n_rejected = fj_struct.rejected().size();
165 }
166
167 // write out result
168 cout << " rap = " << j.rap() << ", phi = " << j.phi() << ", pt = " << j.perp()
169 << " mass = " << j.m()
170 << " [kept: " << j.pieces().size()
171 << ", rejected: " << n_rejected;
172 if (j.has_structure_of<Pruner>()) {
173 cout << ", Rcut: " << j.structure_of<Pruner>().Rcut();
174 }
175 cout << "]" << endl;
176 }
177 }
178
179 // a bit of memory cleaning
180 for (unsigned int i=0; i < groomers.size(); i++) delete groomers[i];
181
182 return 0;
183}
int main()
an example program showing how to use Filter and Pruner in FastJet
Definition: 14-groomers.cc:61
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.
Class to contain structure information for a filtered jet.
Definition: Filter.hh:203
const std::vector< PseudoJet > & rejected() const
returns the subjets that were not kept during the filtering procedure (subtracted if the filter reque...
Definition: Filter.hh:228
Class that helps perform filtering (Butterworth, Davison, Rubin and Salam, arXiv:0802....
Definition: Filter.hh:97
class that is intended to hold a full definition of the jet clusterer
The structure associated with a PseudoJet thas has gone through a Pruner transformer.
Definition: Pruner.hh:189
std::vector< PseudoJet > rejected() const
return the constituents that have been rejected
Definition: Pruner.hh:200
Transformer that prunes a jet.
Definition: Pruner.hh:107
Class to contain pseudojets, including minimal information of use to jet-clustering routines.
Definition: PseudoJet.hh:68
double rap() const
returns the rapidity or some large value when the rapidity is infinite
Definition: PseudoJet.hh:138
const TransformerType::StructureType & structure_of() const
this is a helper to access any structure created by a Transformer (that is, of type Transformer::Stru...
Definition: PseudoJet.hh:1146
double phi() const
returns phi (in the range 0..2pi)
Definition: PseudoJet.hh:123
double perp() const
returns the scalar transverse momentum
Definition: PseudoJet.hh:158
bool has_structure_of() const
check if the PseudoJet has the structure resulting from a Transformer (that is, its structure is comp...
Definition: PseudoJet.hh:1136
double m() const
returns the invariant mass (If m2() is negative then -sqrt(-m2()) is returned, as in CLHEP)
Definition: PseudoJet.hh:1064
std::string description() const
return a string describing what kind of PseudoJet we are dealing with
Definition: PseudoJet.cc:506
virtual std::vector< PseudoJet > pieces() const
retrieve the pieces that make up the jet.
Definition: PseudoJet.cc:781
Base (abstract) class for a jet transformer.
Definition: Transformer.hh:71
virtual std::string description() const =0
This should be overloaded to return a description of the Transformer.
Selector SelectorNHardest(unsigned int n)
select the n hardest objects
Definition: Selector.cc:1074
Selector SelectorPtFractionMin(double fraction)
select objects that carry at least a fraction "fraction" of the reference jet.
Definition: Selector.cc:1365
the FastJet namespace
@ cambridge_algorithm
the longitudinally invariant variant of the cambridge algorithm (aka Aachen algoithm).
@ kt_algorithm
the longitudinally invariant kt algorithm
vector< PseudoJet > sorted_by_pt(const vector< PseudoJet > &jets)
return a vector of jets sorted into decreasing kt2
Definition: PseudoJet.cc:871