FastJet 3.4.3
Loading...
Searching...
No Matches
11-filter.cc
Go to the documentation of this file.
1//----------------------------------------------------------------------
2/// \file
3/// \page Example11 11 - use of filtering
4///
5/// fastjet example program to illustrate the use of the
6/// fastjet::Filter class
7///
8/// We apply different filter examples to either the hardest jet of
9/// the given event, or to the composition of the two hardest jets:
10///
11/// - two examples of a filter keeping a fixed number of subjets (as
12/// in arXiv:0802.2470)
13/// - a "trimmer" i.e. a filter keeping subjets carrying at least a given
14/// fraction of the pt of the jet (arXiv:0912.1342).
15/// - two examples of filter in combination with background subtraction
16///
17/// run it with : ./11-filter < data/single-event.dat
18///
19/// Source code: 11-filter.cc
20//----------------------------------------------------------------------
21
22//FJSTARTHEADER
23// $Id$
24//
25// Copyright (c) 2005-2024, Matteo Cacciari, Gavin P. Salam and Gregory Soyez
26//
27//----------------------------------------------------------------------
28// This file is part of FastJet.
29//
30// FastJet is free software; you can redistribute it and/or modify
31// it under the terms of the GNU General Public License as published by
32// the Free Software Foundation; either version 2 of the License, or
33// (at your option) any later version.
34//
35// The algorithms that underlie FastJet have required considerable
36// development. They are described in the original FastJet paper,
37// hep-ph/0512210 and in the manual, arXiv:1111.6097. If you use
38// FastJet as part of work towards a scientific publication, please
39// quote the version you use and include a citation to the manual and
40// optionally also to hep-ph/0512210.
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//FJENDHEADER
51
52#include <fastjet/PseudoJet.hh>
53#include <fastjet/ClusterSequence.hh>
54#include <fastjet/Selector.hh>
55#include <iostream>
56#include "fastjet/tools/Filter.hh"
57// the following includes are only needed when combining filtering with subtraction
58#include "fastjet/tools/GridMedianBackgroundEstimator.hh"
59#include "fastjet/ClusterSequenceArea.hh"
60#include "fastjet/tools/Subtractor.hh"
61
62#include <cstdio> // needed for io
63
64using namespace fastjet;
65using namespace std;
66
67// a function returning
68// min(Rmax, deltaR_factor * deltaR(j1,j2))
69// where j1 and j2 are the 2 subjets of j
70// if the jet does not have 2 exactly pieces, Rmax is used.
71class DynamicRfilt : public FunctionOfPseudoJet<double>{
72public:
73 // default ctor
74 DynamicRfilt(double Rmax, double deltaR_factor) : _Rmax(Rmax), _deltaR_factor(deltaR_factor){}
75
76 // action of the function
77 double result(const PseudoJet &j) const{
78 if (! j.has_pieces()) return _Rmax;
79
80 vector<PseudoJet> pieces = j.pieces();
81 if (pieces.size() != 2) return _Rmax;
82
83 double deltaR = pieces[0].delta_R(pieces[1]);
84 return min(_Rmax, _deltaR_factor * deltaR);
85 }
86
87private:
88 double _Rmax, _deltaR_factor;
89};
90
91/// an example program showing how to use Filter in FastJet
92int main(){
93 // read in input particles
94 //----------------------------------------------------------
95 vector<PseudoJet> input_particles;
96
97 double px, py , pz, E;
98 while (cin >> px >> py >> pz >> E) {
99 // create a PseudoJet with these components and put it onto
100 // back of the input_particles vector
101 input_particles.push_back(PseudoJet(px,py,pz,E));
102 }
103
104 // get the resulting jets ordered in pt
105 //----------------------------------------------------------
107 // the use of a ClusterSequenceArea (instead of a plain ClusterSequence)
108 // is only needed because we will later combine filtering with area-based
109 // subtraction
110 ClusterSequenceArea clust_seq(input_particles, jet_def,
111 AreaDefinition(active_area_explicit_ghosts));
112 vector<PseudoJet> inclusive_jets = sorted_by_pt(clust_seq.inclusive_jets(5.0));
113
114 // label the columns
115 printf("%5s %15s %15s %15s\n","jet #", "rapidity", "phi", "pt");
116
117 // print out the details for each jet
118 for (unsigned int i = 0; i < inclusive_jets.size(); i++) {
119 printf("%5u %15.8f %15.8f %15.8f\n",
120 i, inclusive_jets[i].rap(), inclusive_jets[i].phi(),
121 inclusive_jets[i].perp());
122 }
123
124 // simple test to avoid that the example below crashes:
125 // make sure there is at least 3 jets above our 5 GeV
126 if (inclusive_jets.size()<3){
127 cout << "Please provide an event with at least 3 jets above 5 GeV" << endl;
128 return 1;
129 }
130
131 // the sample PseudoJet that we will filter
132 // - the hardest jet of the event
133 // - the composition of the second and third hardest jets
134 /// (this shows that the Filter can also be applied to a composite jet)
135 //----------------------------------------------------------
136 vector<PseudoJet> candidates;
137 candidates.push_back(inclusive_jets[0]);
138 candidates.push_back(join(inclusive_jets[1],inclusive_jets[2]));
139
140
141 // create 5 filters
142 //----------------------------------------------------------
143 vector<Filter> filters;
144
145 // 1.
146 // the Cambridge/Aachen filter with Rfilt=0.3 (simpliefied version of arXiv:0802.2470)
147 filters.push_back(Filter(JetDefinition(cambridge_algorithm, 0.3), SelectorNHardest(3)));
148
149 // 2.
150 // the Cambridge/Aachen filter with Rfilt=min(0.3, 0.5*Rbb) as in arXiv:0802.2470
151 SharedPtr<DynamicRfilt> dynamic_Rfilt(new DynamicRfilt(0.3, 0.5));
152 filters.push_back(Filter(dynamic_Rfilt.get(), SelectorNHardest(3)));
153
154 // 3.
155 // Filtering with a pt cut as for trimming (arXiv:0912.1342)
156 filters.push_back(Filter(JetDefinition(kt_algorithm, 0.2), SelectorPtFractionMin(0.03)));
157
158 // 4.
159 // First example of filtering with subtraction of the background: provide rho
160 // First, estimate the background for the given event
161 GridMedianBackgroundEstimator bkgd(4.5, 0.55); // uses particles up to |y|=4.5
162 bkgd.set_particles(input_particles);
163 double rho = bkgd.rho();
164 // Then, define the filter
165 filters.push_back(Filter(JetDefinition(cambridge_algorithm, 0.3), SelectorNHardest(3), rho));
166
167 // 5.
168 // Second example of filtering with subtraction of the background: set a subtractor
169 // First, define a subtractor from a background estimator
170 Subtractor subtractor(&bkgd);
171 // Then, define the filter
173 // Finally, tell the filter about the subtractor
174 filt.set_subtractor(&subtractor);
175 filters.push_back(filt);
176
177
178 // apply the various filters to the test PseudoJet
179 // and show the result
180 //----------------------------------------------------------
181
182 // print out original jet candidates
183 cout << "\nOriginal jets that will be filtered: " << endl;
184 for (vector<PseudoJet>::iterator jit=candidates.begin(); jit!=candidates.end(); jit++){
185 const PseudoJet & c = *jit;
186 cout << " rap = " << c.rap() << ", phi = " << c.phi() << ", pt = " << c.perp()
187 << " [" << c.description() << "]" << endl;
188 }
189
190 // loop on filters
191 for (vector<Filter>::iterator it=filters.begin(); it!=filters.end(); it++){
192 const Filter & f = *it;
193 cout << "\nUsing filter: " << f.description() << endl;
194
195 // loop on jet candidates
196 for (vector<PseudoJet>::iterator jit=candidates.begin(); jit!=candidates.end(); jit++){
197 const PseudoJet & c = *jit;
198
199 // apply filter f to jet c
200 PseudoJet j = f(c);
201
202 // access properties specific to the Filter
203 //
204 // We first make sure that the jet indeed has a structure
205 // compatible with the result of a Filter (using
206 // has_structure_of()), and then retrieve the pieces rejected by the
207 // filter (using structure_of())
208 assert(j.has_structure_of<Filter>());
209 const Filter::StructureType & fj_struct = j.structure_of<Filter>();
210
211 // write out result
212 cout << " rap = " << j.rap() << ", phi = " << j.phi() << ", pt = " << j.perp()
213 << " [kept: " << j.pieces().size() << ", rejected: "
214 << fj_struct.rejected().size() << "]" << endl;
215 }
216 }
217
218 return 0;
219}
int main()
an example program showing how to use Filter in FastJet
Definition 11-filter.cc:92
class that holds a generic area definition
General class for user to obtain ClusterSequence with additional area information.
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
void set_subtractor(const FunctionOfPseudoJet< PseudoJet > *subtractor_in)
Set a subtractor that is applied to all individual subjets before deciding which ones to keep.
Definition Filter.hh:144
virtual std::string description() const
class description
Definition Filter.cc:50
base class providing interface for a generic function of a PseudoJet
virtual TOut result(const PseudoJet &pj) const =0
the action of the function this has to be overloaded in derived classes
Background Estimator based on the median pt/area of a set of grid cells.
void set_particles(const std::vector< PseudoJet > &particles) override
tell the background estimator that it has a new event, composed of the specified particles.
double rho() const override
returns rho, the median background density per unit area
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
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...
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...
virtual bool has_pieces() const
returns true if a jet has pieces
Definition PseudoJet.cc:772
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
An implementation of shared pointers that is broadly similar to C++11 shared_ptr (https://en....
Definition SharedPtr.hh:341
T * get() const
get the stored pointer
Definition SharedPtr.hh:473
Class that helps perform jet background subtraction.
Definition Subtractor.hh:62
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