FastJet 3.4.2
fastjet_subtraction.cc
1
2//FJSTARTHEADER
3// $Id$
4//
5// Copyright (c) 2005-2023, Matteo Cacciari, Gavin P. Salam and Gregory Soyez
6//
7//----------------------------------------------------------------------
8// This file is part of FastJet.
9//
10// FastJet is free software; you can redistribute it and/or modify
11// it under the terms of the GNU General Public License as published by
12// the Free Software Foundation; either version 2 of the License, or
13// (at your option) any later version.
14//
15// The algorithms that underlie FastJet have required considerable
16// development. They are described in the original FastJet paper,
17// hep-ph/0512210 and in the manual, arXiv:1111.6097. If you use
18// FastJet as part of work towards a scientific publication, please
19// quote the version you use and include a citation to the manual and
20// optionally also to hep-ph/0512210.
21//
22// FastJet is distributed in the hope that it will be useful,
23// but WITHOUT ANY WARRANTY; without even the implied warranty of
24// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
25// GNU General Public License for more details.
26//
27// You should have received a copy of the GNU General Public License
28// along with FastJet. If not, see <http://www.gnu.org/licenses/>.
29//----------------------------------------------------------------------
30//FJENDHEADER
31
32
33
34
35
36
37
38//----------------------------------------------------------------------
39// fastjet subtraction example program.
40//
41// Compile it with: make fastjet_subtraction
42// run it with : ./fastjet_subtraction < data/Pythia-dijet-ptmin100-lhc-pileup-1ev.dat
43//
44//----------------------------------------------------------------------
45#include "fastjet/PseudoJet.hh"
46#include "fastjet/ClusterSequenceArea.hh"
47#include<iostream> // needed for io
48#include<sstream> // needed for internal io
49#include<vector>
50
51using namespace std;
52
53// A declaration of a function that pretty prints a list of jets
54// The subtraction is also performed inside this function
56 const vector<fastjet::PseudoJet> &);
57
58/// an example program showing how to use fastjet
59int main (int argc, char ** argv) {
60
61 // since we use here simulated data we can split the hard event
62 // from the full (i.e. with pileup added) one
63 vector<fastjet::PseudoJet> hard_event, full_event;
64
65 // maximum rapidity that we accept for the input particles
66 double etamax = 6.0;
67
68 // read in input particles. Keep the hard event generated by PYTHIA
69 // separated from the full event, so as to be able to gauge the
70 // "goodness" of the subtraction from the full event which also
71 // includes pileup
72 string line;
73 int nsub = 0;
74 while (getline(cin, line)) {
75 istringstream linestream(line);
76 // take substrings to avoid problems when there are extra "pollution"
77 // characters (e.g. line-feed).
78 if (line.substr(0,4) == "#END") {break;}
79 if (line.substr(0,9) == "#SUBSTART") {
80 // if more sub events follow, make copy of hard one here
81 if (nsub == 1) hard_event = full_event;
82 nsub += 1;
83 }
84 if (line.substr(0,1) == "#") {continue;}
85 valarray<double> fourvec(4);
86 linestream >> fourvec[0] >> fourvec[1] >> fourvec[2] >> fourvec[3];
87 fastjet::PseudoJet psjet(fourvec);
88 psjet.set_user_index(0);
89 // push event onto back of full_event vector
90 if (abs(psjet.rap()) < etamax) {full_event.push_back(psjet);}
91 }
92
93 // if we have read in only one event, copy it across here...
94 if (nsub == 1) hard_event = full_event;
95
96 // if there was nothing in the event
97 if (nsub == 0) {
98 cerr << "Error: read empty event\n";
99 exit(-1);
100 }
101
102
103 // create an object that represents your choice of jet algorithm and
104 // the associated parameters
105 double R = 0.7;
107 fastjet::JetDefinition jet_def(fastjet::kt_algorithm, R, strategy);
108 //fastjet::JetDefinition jet_def(fastjet::cambridge_algorithm, R, strategy);
109 //fastjet::JetDefinition jet_def(fastjet::antikt_algorithm, R, strategy);
110
111 // create an object that specifies how we to define the (active) area
112 double ghost_etamax = 6.0;
113 int active_area_repeats = 1;
114 double ghost_area = 0.01;
115 fastjet::GhostedAreaSpec ghost_spec(ghost_etamax, active_area_repeats,
116 ghost_area);
117 fastjet::AreaDefinition area_def(ghost_spec);
118 //fastjet::AreaDefinition area_def(fastjet::VoronoiAreaSpec(1.0));
119
120 // run the jet clustering with the above jet definition. hard event first
121 fastjet::ClusterSequenceArea clust_seq(hard_event, jet_def, area_def);
122 //fastjet::ClusterSequencePassiveArea clust_seq(hard_event, jet_def);
123
124
125 // extract the inclusive jets with pt > 5 GeV, sorted by pt
126 double ptmin = 5.0;
127 vector<fastjet::PseudoJet> inclusive_jets = clust_seq.inclusive_jets(ptmin);
128
129 // print them out
130 cout << "" <<endl;
131 cout << "Hard event only"<<endl;
132 cout << "Number of input particles: "<<hard_event.size()<<endl;
133 cout << "Strategy used: "<<clust_seq.strategy_string()<<endl;
134 cout << "Printing inclusive jets with pt > "<< ptmin<<" GeV\n";
135 cout << "---------------------------------------\n";
136 print_jets(clust_seq, inclusive_jets);
137 cout << endl;
138
139
140 // repeat everything on the full event
141
142 // run the jet clustering with the above jet definition
143 fastjet::ClusterSequenceArea clust_seq_full(full_event, jet_def, area_def);
144
145 // extract the inclusive jets with pt > 20 GeV, sorted by pt
146 ptmin = 20.0;
147 inclusive_jets = clust_seq_full.inclusive_jets(ptmin);
148
149 // print them out
150 cout << "Full event, with pileup, and its subtraction"<<endl;
151 cout << "Number of input particles: "<<full_event.size()<<endl;
152 cout << "Strategy used: "<<clust_seq_full.strategy_string()<<endl;
153 cout << "Printing inclusive jets with pt > "<< ptmin<<" GeV (before subtraction)\n";
154 cout << "---------------------------------------\n";
155 print_jets(clust_seq_full, inclusive_jets);
156 cout << endl;
157
158
159}
160
161
162//----------------------------------------------------------------------
163/// a function that pretty prints a list of jets, and performs the subtraction
164/// in two different ways, using a generic ClusterSequenceAreaBase
165/// type object.
166void print_jets (const fastjet::ClusterSequenceAreaBase & clust_seq,
167 const vector<fastjet::PseudoJet> & unsorted_jets ) {
168
169 // sort jets into increasing pt
170 vector<fastjet::PseudoJet> jets = sorted_by_pt(unsorted_jets);
171
172 // the corrected jets will go in here
173 vector<fastjet::PseudoJet> corrected_jets(jets.size());
174
175 // get median pt per unit area -- you must decide the rapidity range
176 // in which you measure the activity. (If doing a ghost-based active
177 // area, make sure that your ghosts go up at least to \sim range+R).
178 double range = 5.0;
179 double median_pt_per_area = clust_seq.median_pt_per_unit_area(range);
180 double median_pt_per_area4vector = clust_seq.median_pt_per_unit_area_4vector(range);
181
182 printf(" ijet rap phi Pt area Pt corr (rap corr phi corr Pt corr)ext\n");
183 for (size_t j = 0; j < jets.size(); j++) {
184
185 // get area of each jet
186 double area = jets[j].area();
187
188 // "standard" correction. Subtract only the Pt
189 double pt_corr = jets[j].perp() - area*median_pt_per_area;
190
191 // "extended" correction
192 fastjet::PseudoJet sub_4vect =
193 median_pt_per_area4vector*jets[j].area_4vector();
194 if (sub_4vect.perp2() >= jets[j].perp2() ||
195 sub_4vect.E() >= jets[j].E()) {
196 // if the correction is too large, set the jet to zero
197 corrected_jets[j] = 0.0 * jets[j];
198 } else {
199 // otherwise do an E-scheme subtraction
200 corrected_jets[j] = jets[j] - sub_4vect;
201 }
202 // NB We could also leave out the above "if":
203 // corrected_jets[j] = jets[j] - sub_4vect;
204 // but the result would be different, since we would not avoid
205 // jets with negative Pt or energy
206
207 printf("%5u %7.3f %7.3f %9.3f %7.3f %9.3f %7.3f %7.3f %9.3f\n",
208 j,jets[j].rap(),jets[j].phi(),jets[j].perp(), area, pt_corr,
209 corrected_jets[j].rap(),corrected_jets[j].phi(), corrected_jets[j].perp());
210 }
211
212 cout << endl;
213 cout << "median pt_over_area = " << median_pt_per_area << endl;
214 cout << "median pt_over_area4vector = " << median_pt_per_area4vector << endl << endl;
215
216
217}
218
219
220
221
222
223
224
int main()
an example program showing how to use fastjet
Definition: 01-basic.cc:52
class that holds a generic area definition
base class that sets interface for extensions of ClusterSequence that provide information about the a...
double median_pt_per_unit_area_4vector(const Selector &selector) const
the median of (pt/area_4vector) for jets contained within the selector range, making use also of the ...
double median_pt_per_unit_area(const Selector &selector) const
the median of (pt/area) for jets contained within the selector range, making use also of the info on ...
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 perp2() const
returns the squared transverse momentum
Definition: PseudoJet.hh:156
virtual PseudoJet area_4vector() const
return the jet 4-vector area.
Definition: PseudoJet.cc:833
void print_jets(const vector< fastjet::PseudoJet > &)
a function that pretty prints a list of jets
Strategy
the various options for the algorithmic strategy to adopt in clustering events with kt and cambridge ...
@ Best
automatic selection of the best (based on N), including the LazyTiled strategies that are new to FJ3....
@ 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