FastJet 3.4.2
07-subtraction.cc
Go to the documentation of this file.
1//----------------------------------------------------------------------
2/// \file
3/// \page Example07 07 - subtracting jet background contamination
4///
5/// fastjet subtraction example program.
6///
7/// run it with : ./07-subtraction < data/Pythia-Zp2jets-lhc-pileup-1ev.dat
8///
9/// Source code: 07-subtraction.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/PseudoJet.hh"
43#include "fastjet/ClusterSequenceArea.hh"
44#include "fastjet/Selector.hh"
45#include "fastjet/tools/JetMedianBackgroundEstimator.hh"
46#include "fastjet/tools/Subtractor.hh"
47#include <iostream> // needed for io
48#include <fastjet/config.h> // for the FASTJET_VERSION_NUMBER preprocessor symbol
49
50using namespace std;
51using namespace fastjet;
52
53int main(){
54
55 // read in input particles
56 //
57 // since we use here simulated data we can split the hard event
58 // from the full (i.e. with pileup added) one
59 //----------------------------------------------------------
60
61 vector<PseudoJet> hard_event, full_event;
62
63 // read in input particles. Keep the hard event generated by PYTHIA
64 // separated from the full event, so as to be able to gauge the
65 // "goodness" of the subtraction from the full event, which also
66 // includes pileup
67 double particle_maxrap = 5.0;
68
69 string line;
70 int nsub = 0; // counter to keep track of which sub-event we're reading
71 while (getline(cin, line)) {
72 istringstream linestream(line);
73 // take substrings to avoid problems when there are extra "pollution"
74 // characters (e.g. line-feed).
75 if (line.substr(0,4) == "#END") {break;}
76 if (line.substr(0,9) == "#SUBSTART") {
77 // if more sub events follow, make copy of first one (the hard one) here
78 if (nsub == 1) hard_event = full_event;
79 nsub += 1;
80 }
81 if (line.substr(0,1) == "#") {continue;}
82 double px,py,pz,E;
83 linestream >> px >> py >> pz >> E;
84 // you can construct
85 PseudoJet particle(px,py,pz,E);
86
87 // push event onto back of full_event vector
88 if (abs(particle.rap()) <= particle_maxrap) full_event.push_back(particle);
89 }
90
91 // if we have read in only one event, copy it across here...
92 if (nsub == 1) hard_event = full_event;
93
94 // if there was nothing in the event
95 if (nsub == 0) {
96 cerr << "Error: read empty event\n";
97 exit(-1);
98 }
99
100
101 // create a jet definition for the clustering
102 // We use the anti-kt algorithm with a radius of 0.5
103 //----------------------------------------------------------
104 double R = 0.5;
105 JetDefinition jet_def(antikt_algorithm, R);
106
107 // create an area definition for the clustering
108 //----------------------------------------------------------
109 // ghosts should go up to the acceptance of the detector or
110 // (with infinite acceptance) at least 2R beyond the region
111 // where you plan to investigate jets.
112 double ghost_maxrap = 6.0;
113 GhostedAreaSpec area_spec(ghost_maxrap);
114 AreaDefinition area_def(active_area, area_spec);
115
116 // run the jet clustering with the above jet and area definitions
117 // for both the hard and full event
118 //
119 // We retrieve the jets above 7 GeV in both case (note that the
120 // 7-GeV cut we be applied again later on after we subtract the jets
121 // from the full event)
122 // ----------------------------------------------------------
123 ClusterSequenceArea clust_seq_hard(hard_event, jet_def, area_def);
124 ClusterSequenceArea clust_seq_full(full_event, jet_def, area_def);
125
126 double ptmin = 7.0;
127 vector<PseudoJet> hard_jets = sorted_by_pt(clust_seq_hard.inclusive_jets(ptmin));
128 vector<PseudoJet> full_jets = sorted_by_pt(clust_seq_full.inclusive_jets(ptmin));
129
130 // Now turn to the estimation of the background (for the full event)
131 //
132 // There are different ways to do that. In general, this also
133 // requires clustering the particles that will be handled internally
134 // in FastJet.
135 //
136 // The suggested way to proceed is to use a BackgroundEstimator
137 // constructed from the following 3 arguments:
138 // - a jet definition used to cluster the particles.
139 // . We strongly recommend using the kt or Cambridge/Aachen
140 // algorithm (a warning will be issued otherwise)
141 // . The choice of the radius is a bit more subtle. R=0.4 has
142 // been chosen to limit the impact of hard jets; in samples of
143 // dominantly sparse events it may cause the UE/pileup to be
144 // underestimated a little, a slightly larger value (0.5 or
145 // 0.6) may be better.
146 // - An area definition for which we recommend the use of explicit
147 // ghosts (i.e. active_area_explicit_ghosts)
148 // As mentionned in the area example (06-area.cc), ghosts should
149 // extend sufficiently far in rapidity to cover the jets used in
150 // the computation of the background (see also the comment below)
151 // - A Selector specifying the range over which we will keep the
152 // jets entering the estimation of the background (you should
153 // thus make sure the ghosts extend far enough in rapidity to
154 // cover the range, a warning will be issued otherwise).
155 // In this particular example, the two hardest jets in the event
156 // are removed from the background estimation
157 // ----------------------------------------------------------
158 JetDefinition jet_def_bkgd(kt_algorithm, 0.4);
159 AreaDefinition area_def_bkgd(active_area_explicit_ghosts,
160 GhostedAreaSpec(ghost_maxrap));
161 Selector selector = SelectorAbsRapMax(4.5) * (!SelectorNHardest(2));
162 JetMedianBackgroundEstimator bkgd_estimator(selector, jet_def_bkgd, area_def_bkgd);
163
164 // To help manipulate the background estimator, we also provide a
165 // transformer that allows to apply directly the background
166 // subtraction on the jets. This will use the background estimator
167 // to compute rho for the jets to be subtracted.
168 // ----------------------------------------------------------
169 Subtractor subtractor(&bkgd_estimator);
170
171 // since FastJet 3.1.0, rho_m is supported natively in background
172 // estimation (both JetMedianBackgroundEstimator and
173 // GridMedianBackgroundEstimator).
174 //
175 // For backward-compatibility reasons its use is by default switched off
176 // (as is the enforcement of m>0 for the subtracted jets). The
177 // following 2 lines of code switch these on. They are strongly
178 // recommended and should become the default in future versions of
179 // FastJet.
180 //
181 // Note that we also illustrate the use of the
182 // FASTJET_VERSION_NUMBER macro
183#if FASTJET_VERSION_NUMBER >= 30100
184 subtractor.set_use_rho_m(true);
185 subtractor.set_safe_mass(true);
186#endif
187
188 // Finally, once we have an event, we can just tell the background
189 // estimator to use that list of particles
190 // This could be done directly when declaring the background
191 // estimator but the usage below can more easily be accomodated to a
192 // loop over a set of events.
193 // ----------------------------------------------------------
194 bkgd_estimator.set_particles(full_event);
195
196 // show a summary of what was done so far
197 // - the description of the algorithms, areas and ranges used
198 // - the background properties
199 // - the jets in the hard event
200 //----------------------------------------------------------
201 cout << "Main clustering:" << endl;
202 cout << " Ran: " << jet_def.description() << endl;
203 cout << " Area: " << area_def.description() << endl;
204 cout << " Particles up to |y|=" << particle_maxrap << endl;
205 cout << endl;
206
207 cout << "Background estimation:" << endl;
208 cout << " " << bkgd_estimator.description() << endl << endl;;
209 cout << " Giving, for the full event" << endl;
210 BackgroundEstimate bkgd_estimate = bkgd_estimator.estimate();
211 cout << " rho = " << bkgd_estimate.rho() << endl;
212 cout << " sigma = " << bkgd_estimate.sigma() << endl;
213 cout << " rho_m = " << bkgd_estimate.rho_m() << endl;
214 cout << " sigma_m = " << bkgd_estimate.sigma_m() << endl;
215 cout << endl;
216
217 cout << "Jets above " << ptmin << " GeV in the hard event (" << hard_event.size() << " particles)" << endl;
218 cout << "---------------------------------------\n";
219 printf("%5s %15s %15s %15s %15s %15s\n","jet #", "rapidity", "phi", "pt", "m", "area");
220 for (unsigned int i = 0; i < hard_jets.size(); i++) {
221 printf("%5u %15.8f %15.8f %15.8f %15.8f %15.8f\n", i,
222 hard_jets[i].rap(), hard_jets[i].phi(), hard_jets[i].pt(), hard_jets[i].m(),
223 hard_jets[i].area());
224 }
225 cout << endl;
226
227 // Once the background properties have been computed, subtraction
228 // can be applied on the jets. Subtraction is performed on the
229 // full 4-vector
230 //
231 // We output the jets before and after subtraction
232 // ----------------------------------------------------------
233 cout << "Jets above " << ptmin << " GeV in the full event (" << full_event.size() << " particles)" << endl;
234 cout << "---------------------------------------\n";
235 printf("%5s %15s %15s %15s %15s %15s %15s %15s %15s %15s\n","jet #", "rapidity", "phi", "pt", "m", "area", "rap_sub", "phi_sub", "pt_sub", "m_sub");
236 unsigned int idx=0;
237
238 // get the subtracted jets
239 vector<PseudoJet> subtracted_jets = subtractor(full_jets);
240
241 for (unsigned int i=0; i<full_jets.size(); i++){
242 // re-apply the pt cut
243 if (subtracted_jets[i].pt2() >= ptmin*ptmin){
244 printf("%5u %15.8f %15.8f %15.8f %15.8f %15.8f %15.8f %15.8f %15.8f %15.8f\n", idx,
245 full_jets[i].rap(), full_jets[i].phi(), full_jets[i].pt(), full_jets[i].m(),
246 full_jets[i].area(),
247 subtracted_jets[i].rap(), subtracted_jets[i].phi(),
248 subtracted_jets[i].pt(),
249 subtracted_jets[i].m());
250 idx++;
251 }
252 }
253
254 return 0;
255}
int main()
an example program showing how to use fastjet
Definition: 01-basic.cc:52
class that holds a generic area definition
/// a class that holds the result of the calculation
double sigma_m() const
fluctuations in the purely longitudinal (particle-mass-induced) component of the background density p...
double rho() const
background density per unit area
double sigma() const
background fluctuations per unit square-root area must be multipled by sqrt(area) to get fluctuations...
double rho_m() const
purely longitudinal (particle-mass-induced) component of the background density per unit area
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 estimate the pt density of the background per unit area, using the median of the distributio...
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
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 SelectorAbsRapMax(double absrapmax)
select objects with |rap| <= absrapmax
Definition: Selector.cc:880
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