FastJet  3.3.0
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 //STARTHEADER
23 // $Id: 11-filter.cc 3252 2013-10-28 15:22:36Z salam $
24 //
25 // Copyright (c) 2005-2011, 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 and are described in hep-ph/0512210. If you use
37 // FastJet as part of work towards a scientific publication, please
38 // include a citation to the FastJet paper.
39 //
40 // FastJet is distributed in the hope that it will be useful,
41 // but WITHOUT ANY WARRANTY; without even the implied warranty of
42 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
43 // GNU General Public License for more details.
44 //
45 // You should have received a copy of the GNU General Public License
46 // along with FastJet. If not, see <http://www.gnu.org/licenses/>.
47 //----------------------------------------------------------------------
48 //ENDHEADER
49 
50 #include <fastjet/PseudoJet.hh>
51 #include <fastjet/ClusterSequence.hh>
52 #include <fastjet/Selector.hh>
53 #include <iostream>
54 #include "fastjet/tools/Filter.hh"
55 // the following includes are only needed when combining filtering with subtraction
56 #include "fastjet/tools/GridMedianBackgroundEstimator.hh"
57 #include "fastjet/ClusterSequenceArea.hh"
58 #include "fastjet/tools/Subtractor.hh"
59 
60 #include <cstdio> // needed for io
61 
62 using namespace fastjet;
63 using namespace std;
64 
65 // a function returning
66 // min(Rmax, deltaR_factor * deltaR(j1,j2))
67 // where j1 and j2 are the 2 subjets of j
68 // if the jet does not have 2 exactly pieces, Rmax is used.
69 class DynamicRfilt : public FunctionOfPseudoJet<double>{
70 public:
71  // default ctor
72  DynamicRfilt(double Rmax, double deltaR_factor) : _Rmax(Rmax), _deltaR_factor(deltaR_factor){}
73 
74  // action of the function
75  double result(const PseudoJet &j) const{
76  if (! j.has_pieces()) return _Rmax;
77 
78  vector<PseudoJet> pieces = j.pieces();
79  if (pieces.size() != 2) return _Rmax;
80 
81  double deltaR = pieces[0].delta_R(pieces[1]);
82  return min(_Rmax, _deltaR_factor * deltaR);
83  }
84 
85 private:
86  double _Rmax, _deltaR_factor;
87 };
88 
89 /// an example program showing how to use Filter in FastJet
90 int main(){
91  // read in input particles
92  //----------------------------------------------------------
93  vector<PseudoJet> input_particles;
94 
95  double px, py , pz, E;
96  while (cin >> px >> py >> pz >> E) {
97  // create a PseudoJet with these components and put it onto
98  // back of the input_particles vector
99  input_particles.push_back(PseudoJet(px,py,pz,E));
100  }
101 
102  // get the resulting jets ordered in pt
103  //----------------------------------------------------------
104  JetDefinition jet_def(cambridge_algorithm, 1.2);
105  // the use of a ClusterSequenceArea (instead of a plain ClusterSequence)
106  // is only needed because we will later combine filtering with area-based
107  // subtraction
108  ClusterSequenceArea clust_seq(input_particles, jet_def,
109  AreaDefinition(active_area_explicit_ghosts));
110  vector<PseudoJet> inclusive_jets = sorted_by_pt(clust_seq.inclusive_jets(5.0));
111 
112  // label the columns
113  printf("%5s %15s %15s %15s\n","jet #", "rapidity", "phi", "pt");
114 
115  // print out the details for each jet
116  for (unsigned int i = 0; i < inclusive_jets.size(); i++) {
117  printf("%5u %15.8f %15.8f %15.8f\n",
118  i, inclusive_jets[i].rap(), inclusive_jets[i].phi(),
119  inclusive_jets[i].perp());
120  }
121 
122  // simple test to avoid that the example below crashes:
123  // make sure there is at least 3 jets above our 5 GeV
124  if (inclusive_jets.size()<3){
125  cout << "Please provide an event with at least 3 jets above 5 GeV" << endl;
126  return 1;
127  }
128 
129  // the sample PseudoJet that we will filter
130  // - the hardest jet of the event
131  // - the composition of the second and third hardest jets
132  /// (this shows that the Filter can also be applied to a composite jet)
133  //----------------------------------------------------------
134  vector<PseudoJet> candidates;
135  candidates.push_back(inclusive_jets[0]);
136  candidates.push_back(join(inclusive_jets[1],inclusive_jets[2]));
137 
138 
139  // create 5 filters
140  //----------------------------------------------------------
141  vector<Filter> filters;
142 
143  // 1.
144  // the Cambridge/Aachen filter with Rfilt=0.3 (simpliefied version of arXiv:0802.2470)
145  filters.push_back(Filter(JetDefinition(cambridge_algorithm, 0.3), SelectorNHardest(3)));
146 
147  // 2.
148  // the Cambridge/Aachen filter with Rfilt=min(0.3, 0.5*Rbb) as in arXiv:0802.2470
149  SharedPtr<DynamicRfilt> dynamic_Rfilt(new DynamicRfilt(0.3, 0.5));
150  filters.push_back(Filter(dynamic_Rfilt.get(), SelectorNHardest(3)));
151 
152  // 3.
153  // Filtering with a pt cut as for trimming (arXiv:0912.1342)
154  filters.push_back(Filter(JetDefinition(kt_algorithm, 0.2), SelectorPtFractionMin(0.03)));
155 
156  // 4.
157  // First example of filtering with subtraction of the background: provide rho
158  // First, estimate the background for the given event
159  GridMedianBackgroundEstimator bkgd(4.5, 0.55); // uses particles up to |y|=4.5
160  bkgd.set_particles(input_particles);
161  double rho = bkgd.rho();
162  // Then, define the filter
163  filters.push_back(Filter(JetDefinition(cambridge_algorithm, 0.3), SelectorNHardest(3), rho));
164 
165  // 5.
166  // Second example of filtering with subtraction of the background: set a subtractor
167  // First, define a subtractor from a background estimator
168  Subtractor subtractor(&bkgd);
169  // Then, define the filter
171  // Finally, tell the filter about the subtractor
172  filt.set_subtractor(&subtractor);
173  filters.push_back(filt);
174 
175 
176  // apply the various filters to the test PseudoJet
177  // and show the result
178  //----------------------------------------------------------
179 
180  // print out original jet candidates
181  cout << "\nOriginal jets that will be filtered: " << endl;
182  for (vector<PseudoJet>::iterator jit=candidates.begin(); jit!=candidates.end(); jit++){
183  const PseudoJet & c = *jit;
184  cout << " rap = " << c.rap() << ", phi = " << c.phi() << ", pt = " << c.perp()
185  << " [" << c.description() << "]" << endl;
186  }
187 
188  // loop on filters
189  for (vector<Filter>::iterator it=filters.begin(); it!=filters.end(); it++){
190  const Filter & f = *it;
191  cout << "\nUsing filter: " << f.description() << endl;
192 
193  // loop on jet candidates
194  for (vector<PseudoJet>::iterator jit=candidates.begin(); jit!=candidates.end(); jit++){
195  const PseudoJet & c = *jit;
196 
197  // apply filter f to jet c
198  PseudoJet j = f(c);
199 
200  // access properties specific to the Filter
201  //
202  // We first make sure that the jet indeed has a structure
203  // compatible with the result of a Filter (using
204  // has_structure_of()), and then retrieve the pieces rejected by the
205  // filter (using structure_of())
206  assert(j.has_structure_of<Filter>());
207  const Filter::StructureType & fj_struct = j.structure_of<Filter>();
208 
209  // write out result
210  cout << " rap = " << j.rap() << ", phi = " << j.phi() << ", pt = " << j.perp()
211  << " [kept: " << j.pieces().size() << ", rejected: "
212  << fj_struct.rejected().size() << "]" << endl;
213  }
214  }
215 
216  return 0;
217 }
Class that helps perform jet background subtraction.
Definition: Subtractor.hh:62
Selector SelectorPtFractionMin(double fraction)
select objects that carry at least a fraction "fraction" of the reference jet.
Definition: Selector.cc:1365
vector< PseudoJet > sorted_by_pt(const vector< PseudoJet > &jets)
return a vector of jets sorted into decreasing kt2
Definition: PseudoJet.cc:774
General class for user to obtain ClusterSequence with additional area information.
T * get() const
get the stored pointer
Definition: SharedPtr.hh:250
Selector SelectorNHardest(unsigned int n)
select the n hardest objects
Definition: Selector.cc:1074
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
Class that helps perform filtering (Butterworth, Davison, Rubin and Salam, arXiv:0802.2470) and trimming (Krohn, Thaler and Wang, arXiv:0912.1342) on jets, optionally in conjunction with subtraction (Cacciari and Salam, arXiv:0707.1378).
Definition: Filter.hh:97
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
double rho() const
returns rho, the median background density per unit area
class that holds a generic area definition
virtual bool has_pieces() const
returns true if a jet has pieces
Definition: PseudoJet.cc:675
double perp() const
returns the scalar transverse momentum
Definition: PseudoJet.hh:145
Background Estimator based on the median pt/area of a set of grid cells.
the FastJet namespace
void set_particles(const std::vector< PseudoJet > &particles)
tell the background estimator that it has a new event, composed of the specified particles.
virtual std::vector< PseudoJet > pieces() const
retrieve the pieces that make up the jet.
Definition: PseudoJet.cc:684
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
an implementation of C++0x shared pointers (or boost&#39;s)
Definition: SharedPtr.hh:121
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:1012
double phi() const
returns phi (in the range 0..2pi)
Definition: PseudoJet.hh:110
double rap() const
returns the rapidity or some large value when the rapidity is infinite
Definition: PseudoJet.hh:125
base class providing interface for a generic function of a PseudoJet
int main()
an example program showing how to use Filter in FastJet
Definition: 11-filter.cc:90
the longitudinally invariant variant of the cambridge algorithm (aka Aachen algoithm).
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:1002
std::string description() const
return a string describing what kind of PseudoJet we are dealing with
Definition: PseudoJet.cc:413
Class to contain pseudojets, including minimal information of use to jet-clustering routines...
Definition: PseudoJet.hh:67
virtual std::string description() const
class description
Definition: Filter.cc:50
class that is intended to hold a full definition of the jet clusterer