FastJet 3.4.2
07-subtraction-old.cc
Go to the documentation of this file.
1//----------------------------------------------------------------------
2/// \file
3/// \page Example07old 07 - subtracting jet background contamination (old version)
4///
5/// fastjet subtraction example program.
6///
7/// Note that this example is deprecated --- see 07-subtraction.cc
8/// for the newest version --- so it is not built by default
9///
10/// run it with : ./07-subtraction-old < data/Pythia-Zp2jets-lhc-pileup-1ev.dat
11///
12/// Source code: 07-subtraction-old.cc
13//----------------------------------------------------------------------
14
15//FJSTARTHEADER
16// $Id$
17//
18// Copyright (c) 2005-2023, Matteo Cacciari, Gavin P. Salam and Gregory Soyez
19//
20//----------------------------------------------------------------------
21// This file is part of FastJet.
22//
23// FastJet is free software; you can redistribute it and/or modify
24// it under the terms of the GNU General Public License as published by
25// the Free Software Foundation; either version 2 of the License, or
26// (at your option) any later version.
27//
28// The algorithms that underlie FastJet have required considerable
29// development. They are described in the original FastJet paper,
30// hep-ph/0512210 and in the manual, arXiv:1111.6097. If you use
31// FastJet as part of work towards a scientific publication, please
32// quote the version you use and include a citation to the manual and
33// optionally also to hep-ph/0512210.
34//
35// FastJet is distributed in the hope that it will be useful,
36// but WITHOUT ANY WARRANTY; without even the implied warranty of
37// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
38// GNU General Public License for more details.
39//
40// You should have received a copy of the GNU General Public License
41// along with FastJet. If not, see <http://www.gnu.org/licenses/>.
42//----------------------------------------------------------------------
43//FJENDHEADER
44
45#include "fastjet/PseudoJet.hh"
46#include "fastjet/ClusterSequenceArea.hh"
47#include <iostream> // needed for io
48
49using namespace std;
50
51int main (int argc, char ** argv) {
52
53 // read in input particles
54 //
55 // since we use here simulated data we can split the hard event
56 // from the full (i.e. with pileup added) one
57 //----------------------------------------------------------
58
59 vector<fastjet::PseudoJet> hard_event, full_event;
60
61 // read in input particles. Keep the hard event generated by PYTHIA
62 // separated from the full event, so as to be able to gauge the
63 // "goodness" of the subtraction from the full event, which also
64 // includes pileup
65 double particle_maxrap = 5.0;
66
67 string line;
68 int nsub = 0; // counter to keep track of which sub-event we're reading
69 while (getline(cin, line)) {
70 istringstream linestream(line);
71 // take substrings to avoid problems when there are extra "pollution"
72 // characters (e.g. line-feed).
73 if (line.substr(0,4) == "#END") {break;}
74 if (line.substr(0,9) == "#SUBSTART") {
75 // if more sub events follow, make copy of first one (the hard one) here
76 if (nsub == 1) hard_event = full_event;
77 nsub += 1;
78 }
79 if (line.substr(0,1) == "#") {continue;}
80 double px,py,pz,E;
81 linestream >> px >> py >> pz >> E;
82 // you can construct
83 fastjet::PseudoJet particle(px,py,pz,E);
84
85 // push event onto back of full_event vector
86 if (abs(particle.rap()) <= particle_maxrap) full_event.push_back(particle);
87 }
88
89 // if we have read in only one event, copy it across here...
90 if (nsub == 1) hard_event = full_event;
91
92 // if there was nothing in the event
93 if (nsub == 0) {
94 cerr << "Error: read empty event\n";
95 exit(-1);
96 }
97
98
99 // create a jet definition for the clustering
100 // We use the anti-kt algorithm with a radius of 0.5
101 //----------------------------------------------------------
102 double R = 0.5;
104
105 // create an area definition for the clustering
106 //----------------------------------------------------------
107 // ghosts should go up to the acceptance of the detector or
108 // (with infinite acceptance) at least 2R beyond the region
109 // where you plan to investigate jets.
110 double ghost_maxrap = 6.0;
111 fastjet::GhostedAreaSpec area_spec(ghost_maxrap);
112 fastjet::AreaDefinition area_def(fastjet::active_area, area_spec);
113
114 // run the jet clustering with the above jet and area definitions
115 // for both the hard and full event
116 //
117 // We retrieve the jets above 7 GeV in both case (note that the
118 // 7-GeV cut we be applied again later on after we subtract the jets
119 // from the full event)
120 // ----------------------------------------------------------
121 fastjet::ClusterSequenceArea clust_seq_hard(hard_event, jet_def, area_def);
122 fastjet::ClusterSequenceArea clust_seq_full(full_event, jet_def, area_def);
123
124 double ptmin = 7.0;
125 vector<fastjet::PseudoJet> hard_jets = sorted_by_pt(clust_seq_hard.inclusive_jets(ptmin));
126 vector<fastjet::PseudoJet> full_jets = sorted_by_pt(clust_seq_full.inclusive_jets(ptmin));
127
128 // Now turn to the estimation of the background (for the full event)
129 //
130 // This also requires a ClusterSequenceArea.
131 // In general, this ClusterSequenceArea does not need to be the same
132 // as the one used (above) to cluster and extract the jets from the
133 // event:
134 // - We strongly recommend using the kt or Cambridge/Aachen algorithm
135 // (a warning will be issued otherwise)
136 // - The choice of the radius is a bit more subtle. R=0.4 has been
137 // chosen to limit the impact of hard jets; in samples of
138 // dominantly sparse events it may cause the UE/pileup to be
139 // underestimated a little, a slightly larger value (0.5 or 0.6)
140 // may be better.
141 // - For the area definition, we recommend the use of explicit
142 // ghosts (i.e. active_area_explicit_ghosts)
143 // As mentionned in the area example (06-area.cc), ghosts should
144 // extend sufficiently far in rapidity to cover the jets used in
145 // the computation of the background (see also the comment below)
146 //
147 // ----------------------------------------------------------
149 fastjet::GhostedAreaSpec area_spec_bkgd(ghost_maxrap);
150 fastjet::AreaDefinition area_def_bkgd(fastjet::active_area_explicit_ghosts, area_spec_bkgd);
151 fastjet::ClusterSequenceArea clust_seq_bkgd(full_event, jet_def_bkgd, area_def_bkgd);
152
153 // Once you have the ClusterSequenceArea, you can compute the
154 // background. This is estimated over a given range
155 // (RangeDefinition) i.e. only jets within that range will be used
156 // to estimate the background. You shold thus make sure the ghosts
157 // extend far enough in rapidity to cover the range, a warning will
158 // be issued otherwise.
159 //
160 // The simplest way to define a RangeDefinition is through its
161 // maximal |y| extent but other options are possible e.g. through a
162 // minimal and maximal rapidity and minimal and maximal azimuthal
163 // angle. If needed, you can even define your own ranges (a few are
164 // provided with FastJet)
165 //
166 // Finally, the estimation of the background properties rho (the
167 // average density per unit area) and sigma (the average
168 // fluctuations per unit area) is done using
169 // ClusterSequenceArea::get_median_rho_and_sigma(). This takes
170 // 2 main parameters: the range discussed above and a boolean
171 // controlling the use of 4-vector or scalar areas (we suggest using
172 // 4-vector areas)
173 //
174 // ----------------------------------------------------------
175 double range_maxrap = 4.5; // we have a ghost_maxrap of 6.0, particles up to 5
176 fastjet::RangeDefinition range(range_maxrap);
177
178 bool use_4vector_area = true;
179
180 double rho, sigma;
181 clust_seq_bkgd.get_median_rho_and_sigma(range, use_4vector_area, rho, sigma);
182
183 // show a summary of what was done so far
184 // - the description of the algorithms, areas and ranges used
185 // - the background properties
186 // - the jets in the hard event
187 //----------------------------------------------------------
188 cout << "Main clustering:" << endl;
189 cout << " Ran: " << jet_def.description() << endl;
190 cout << " Area: " << area_def.description() << endl;
191 cout << " Particles up to |y|=" << particle_maxrap << endl;
192 cout << endl;
193
194 cout << "Background estimation:" << endl;
195 cout << " Ran " << jet_def_bkgd.description() << endl;
196 cout << " Area: " << area_def_bkgd.description() << endl;
197 cout << " Range: " << range.description() << endl;
198 cout << " Giving, for the full event" << endl;
199 cout << " rho = " << rho << endl;
200 cout << " sigma = " << sigma << endl;
201 cout << endl;
202
203 cout << "Jets above " << ptmin << " GeV in the hard event (" << hard_event.size() << " particles)" << endl;
204 cout << "---------------------------------------\n";
205 printf("%5s %15s %15s %15s %15s\n","jet #", "rapidity", "phi", "pt", "area");
206 for (unsigned int i = 0; i < hard_jets.size(); i++) {
207 printf("%5u %15.8f %15.8f %15.8f %15.8f\n", i,
208 hard_jets[i].rap(), hard_jets[i].phi(), hard_jets[i].perp(),
209 hard_jets[i].area());
210 }
211 cout << endl;
212
213 // Once the background properties have been computed, subtraction
214 // can be applied on the jets
215 //
216 // This uses ClusterSequenceArea::subtracted_jet(jet, rho), with the
217 // ClusterSequence used to cluster the jet and the background
218 // density we have just computed
219 //
220 // (Note that when using scalar areas, subtracted_pt should be used
221 // instead of subtracted_jet)
222 //
223 // We output the jets before and after subtraction
224 // ----------------------------------------------------------
225 cout << "Jets above " << ptmin << " GeV in the full event (" << full_event.size() << " particles)" << endl;
226 cout << "---------------------------------------\n";
227 printf("%5s %15s %15s %15s %15s %15s %15s %15s\n","jet #", "rapidity", "phi", "pt", "area", "rap_sub", "phi_sub", "pt_sub");
228 unsigned int idx=0;
229 for (unsigned int i=0; i<full_jets.size(); i++){
230 // get the subtracted jet
231 fastjet::PseudoJet subtracted_jet = clust_seq_full.subtracted_jet(full_jets[i], rho);
232
233 // re-apply the pt cut
234 if (subtracted_jet.perp2() >= ptmin*ptmin){
235 printf("%5u %15.8f %15.8f %15.8f %15.8f %15.8f %15.8f %15.8f\n", idx,
236 full_jets[i].rap(), full_jets[i].phi(), full_jets[i].perp(),
237 full_jets[i].area(),
238 subtracted_jet.rap(), subtracted_jet.phi(),
239 subtracted_jet.perp());
240 idx++;
241 }
242 }
243
244 return 0;
245}
int main()
an example program showing how to use fastjet
Definition: 01-basic.cc:52
class that holds a generic area definition
General class for user to obtain ClusterSequence with additional area information.
Parameters to configure the computation of jet areas using ghosts.
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
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
double perp2() const
returns the squared transverse momentum
Definition: PseudoJet.hh:156
class for holding a range definition specification, given by limits on rapidity and azimuth.
@ antikt_algorithm
like the k_t but with distance measures dij = min(1/kti^2,1/ktj^2) Delta R_{ij}^2 / R^2 diB = 1/kti^2
@ 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