FastJet 3.4.2
12-boosted_higgs-old.cc
Go to the documentation of this file.
1//----------------------------------------------------------------------
2/// \file
3/// \page Example12old 12 - boosted Higgs tagging (old version)
4///
5/// fastjet example program, illustration of carrying out boosted
6/// Higgs subjet ID analysis
7///
8/// It illustrates two kinds of functionality:
9///
10/// - following the decomposition of a jet into pieces
11/// - following information on a b-tag through the jet
12///
13/// This kind of functionality was used in arXiv:0802.2470
14/// (Butterworth, Davison, Rubin & Salam) for boosted Higgs searches,
15/// and related functionality was used in arXiv:0806.0848 (Kaplan,
16/// Rehermann, Schwartz & Tweedie) in searching for boosted tops
17/// (without b-tag assumptions).
18///
19/// Note that this example is deprecated --- see 12-boosted_higgs.cc
20/// for the newest version --- so it is not built by default
21///
22/// run it with : ./12-boosted_higgs-old < data/HZ-event-Hmass115.dat
23///
24/// Source code: 12-boosted_higgs-old.cc
25//----------------------------------------------------------------------
26
27
28//FJSTARTHEADER
29// $Id$
30//
31// Copyright (c) 2005-2023, Matteo Cacciari, Gavin P. Salam and Gregory Soyez
32//
33//----------------------------------------------------------------------
34// This file is part of FastJet.
35//
36// FastJet is free software; you can redistribute it and/or modify
37// it under the terms of the GNU General Public License as published by
38// the Free Software Foundation; either version 2 of the License, or
39// (at your option) any later version.
40//
41// The algorithms that underlie FastJet have required considerable
42// development. They are described in the original FastJet paper,
43// hep-ph/0512210 and in the manual, arXiv:1111.6097. If you use
44// FastJet as part of work towards a scientific publication, please
45// quote the version you use and include a citation to the manual and
46// optionally also to hep-ph/0512210.
47//
48// FastJet is distributed in the hope that it will be useful,
49// but WITHOUT ANY WARRANTY; without even the implied warranty of
50// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
51// GNU General Public License for more details.
52//
53// You should have received a copy of the GNU General Public License
54// along with FastJet. If not, see <http://www.gnu.org/licenses/>.
55//----------------------------------------------------------------------
56//FJENDHEADER
57
58#include "fastjet/ClusterSequence.hh"
59#include <iostream> // needed for io
60#include <sstream> // needed for internal io
61#include <iomanip>
62#include <cmath>
63
64using namespace std;
65using namespace fastjet;
66
67
68//----------------------------------------------------------------------
69// set up a class to give standard (by default E-scheme)
70// recombination, with additional tracking of flavour information in
71// the user_index.
72//
73// b-tagged particles are assumed to have their user_index set to 1,
74// and other particles should have user_index to 0.
75//
76// Watch out however that, by default, the user_index of a particle is
77// set to -1 and you may not have control over that (e.g. if you
78// compute the jet area using explicit ghosts, the ghosts will have a
79// default user_index of -1). For that reason, if one of the particle
80// being combined has a user index of -1, we assume it is not b-tagged
81// (i.e. we count it as 0 in the recombination)
82//
83// This will work for native algorithms, but not for all plugins
84//----------------------------------------------------------------------
86
87class FlavourRecombiner : public DefRecomb {
88public:
89 FlavourRecombiner(RecombinationScheme recomb_scheme = E_scheme) :
90 DefRecomb(recomb_scheme) {};
91
92 virtual std::string description() const {
93 return DefRecomb::description()+" (with user index addition)";}
94
95 /// recombine pa and pb and put result into pab
96 virtual void recombine(const PseudoJet & pa, const PseudoJet & pb,
97 PseudoJet & pab) const {
98 DefRecomb::recombine(pa,pb,pab);
99 // Note: see the above discussion for the fact that we consider
100 // negative user indices as "0"
101 pab.set_user_index(max(pa.user_index(),0) + max(pb.user_index(),0));
102 }
103};
104
105
106//----------------------------------------------------------------------
107// forward declaration for printing out info about a jet
108//----------------------------------------------------------------------
109ostream & operator<<(ostream &, PseudoJet &);
110
111
112//----------------------------------------------------------------------
113// core of the program
114//----------------------------------------------------------------------
115int main (int argc, char ** argv) {
116
117 vector<PseudoJet> particles;
118
119 // read in data in format px py pz E b-tag [last of these is optional]
120 // lines starting with "#" are considered as comments and discarded
121 //----------------------------------------------------------
122
123 string line;
124 while (getline(cin,line)) {
125 if (line.substr(0,1) == "#") {continue;}
126 istringstream linestream(line);
127 double px,py,pz,E;
128 linestream >> px >> py >> pz >> E;
129
130 // optionally read in btag information
131 int btag;
132 if (! (linestream >> btag)) btag = 0;
133
134 // construct the particle
135 PseudoJet particle(px,py,pz,E);
136 particle.set_user_index(btag); // btag info goes in user index, for flavour tracking
137 particles.push_back(particle);
138 }
139
140
141 // set up the jet finding
142 //
143 // This also shows how to use the "FlavourRecombiner" user-defined
144 // recombiner
145 // ----------------------------------------------------------
146 double R = 1.2;
147 FlavourRecombiner flav_recombiner; // for tracking flavour
148 JetDefinition jet_def(cambridge_algorithm, R, &flav_recombiner);
149
150
151 // run the jet finding; find the hardest jet
152 ClusterSequence cs(particles, jet_def);
153 vector<PseudoJet> jets = sorted_by_pt(cs.inclusive_jets());
154
155 cout << "Ran: " << jet_def.description() << endl << endl;
156 cout << "Hardest jet: " << jets[0] << endl << endl;
157
158 // now do the subjet decomposition
159 //----------------------------------------------------------
160 //
161 // when unpeeling a C/A jet, often only a very soft piece may break off;
162 // the mass_drop_threshold indicates how much "lighter" the heavier of the two
163 // resulting pieces must be in order for us to consider that we've really
164 // seen some form of substructure
165 double mass_drop_threshold = 0.667;
166 // QCD backgrounds that give larger jet masses have a component
167 // where a quite soft gluon is emitted; to eliminate part of this
168 // one can place a cut on the asymmetry of the branching;
169 //
170 // Here the cut is expressed in terms of y, the kt-distance scaled
171 // to the squared jet mass; an easier way to see it is in terms of
172 // a requirement on the momentum fraction in the splitting: z/(1-z)
173 // and (1-z)/z > rtycut^2 [the correspondence holds only at LO]
174 double rtycut = 0.3;
175
176 PseudoJet this_jet = jets[0];
177 PseudoJet parent1, parent2;
178 bool had_parents;
179
180 while ((had_parents = this_jet.has_parents(parent1,parent2))) {
181 // make parent1 the more massive jet
182 if (parent1.m() < parent2.m()) swap(parent1,parent2);
183
184 // if we pass the conditions on the mass drop and its degree of
185 // asymmetry (z/(1-z) \sim kt_dist/m^2 > rtycut), then we've found
186 // something interesting, so exit the loop
187 if (parent1.m() < mass_drop_threshold * this_jet.m() &&
188 parent1.kt_distance(parent2) > pow(rtycut,2) * this_jet.m2()) {
189 break;
190 } else {
191 // otherwise try a futher decomposition on the more massive jet
192 this_jet = parent1;
193 }
194 }
195
196 // look to see what we found
197 if (!had_parents) {
198 cout << "Did not find suitable hard substructure in this event." << endl;
199 return 0;
200 }
201
202 cout << "Found suitable pair of subjets: " << endl;
203 cout << " " << parent1 << endl;
204 cout << " " << parent2 << endl;
205 cout << "Total = " << endl;
206 cout << " " << this_jet << endl << endl;
207
208 // next we "filter" it, to remove UE & pileup contamination
209 //----------------------------------------------------------
210 //
211 // [there are two ways of doing this; here we directly use the
212 // exsiting cluster sequence and find the exclusive subjets of
213 // this_jet (i.e. work backwards within the cs starting from
214 // this_jet); alternatively one can recluster just the
215 // constituents of the jet]
216 //
217 // first get separation between the subjets (called Rbb -- assuming it's a Higgs!)
218 double Rbb = sqrt(parent1.squared_distance(parent2));
219 double Rfilt = min(Rbb/2, 0.3); // somewhat arbitrary choice
220 unsigned nfilt = 3; // number of pieces we'll take
221 cout << "Subjet separation (Rbb) = " << Rbb << ", Rfilt = " << Rfilt << endl;
222
223 double dcut = pow(Rfilt/R,2); // for C/A get a view at Rfilt by
224 // using a dcut=(Rfilt/R)^2
225 vector<PseudoJet> filt_subjets = sorted_by_pt(this_jet.exclusive_subjets(dcut));
226
227 // now print out the filtered jets and reconstruct total
228 // at the same time
229 cout << "Filtered pieces are " << endl;
230 cout << " " << filt_subjets[0] << endl;
231 PseudoJet filtered_total = filt_subjets[0];
232 for (unsigned i = 1; i < nfilt && i < filt_subjets.size(); i++) {
233 cout << " " << filt_subjets[i] << endl;
234 flav_recombiner.plus_equal(filtered_total, filt_subjets[i]);
235 }
236 cout << "Filtered total is " << endl;
237 cout << " " << filtered_total << endl;
238
239}
240
241
242//----------------------------------------------------------------------
243// does the actual work for printing out a jet
244//----------------------------------------------------------------------
245ostream & operator<<(ostream & ostr, PseudoJet & jet) {
246 ostr << "pt, y, phi ="
247 << " " << setw(10) << jet.perp()
248 << " " << setw(6) << jet.rap()
249 << " " << setw(6) << jet.phi()
250 << ", mass = " << setw(10) << jet.m()
251 << ", btag = " << jet.user_index();
252 return ostr;
253}
int main()
an example program showing how to use fastjet
Definition: 01-basic.cc:52
ostream & operator<<(ostream &, PseudoJet &)
does the actual work for printing out a jet
deals with clustering
A class that will provide the recombination scheme facilities and/or allow a user to extend these fac...
virtual void recombine(const PseudoJet &pa, const PseudoJet &pb, PseudoJet &pab) const override
recombine pa and pb and put result into pab
virtual std::string description() const override
return a textual description of the recombination scheme implemented here
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
void set_user_index(const int index)
set the user_index, intended to allow the user to add simple identifying information to a particle/je...
Definition: PseudoJet.hh:392
double squared_distance(const PseudoJet &other) const
returns squared cylinder (rap-phi) distance between this jet and another
Definition: PseudoJet.hh:209
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 m2() const
returns the squared invariant mass // like CLHEP
Definition: PseudoJet.hh:163
double kt_distance(const PseudoJet &other) const
returns kt distance (R=1) between this jet and another
Definition: PseudoJet.cc:475
virtual bool has_parents(PseudoJet &parent1, PseudoJet &parent2) const
check if it is the product of a recombination, in which case return the 2 parents through the 'parent...
Definition: PseudoJet.cc:648
std::vector< PseudoJet > exclusive_subjets(const double dcut) const
return a vector of all subjets of the current jet (in the sense of the exclusive algorithm) that woul...
Definition: PseudoJet.cc:704
double m() const
returns the invariant mass (If m2() is negative then -sqrt(-m2()) is returned, as in CLHEP)
Definition: PseudoJet.hh:1064
int user_index() const
return the user_index,
Definition: PseudoJet.hh:389
the FastJet namespace
RecombinationScheme
The various recombination schemes.
void swap(SharedPtr< T > &a, SharedPtr< T > &b)
swapping
Definition: SharedPtr.hh:619
vector< PseudoJet > sorted_by_pt(const vector< PseudoJet > &jets)
return a vector of jets sorted into decreasing kt2
Definition: PseudoJet.cc:871