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