FastJet 3.4.2
12-boosted_higgs.cc
Go to the documentation of this file.
1//----------------------------------------------------------------------
2/// \file
3/// \page Example12 12 - boosted Higgs tagging
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/// - using a boosted higgs tagger
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/// run it with : ./12-boosted_higgs < data/HZ-event-Hmass115.dat
20///
21/// Source code: 12-boosted_higgs.cc
22//----------------------------------------------------------------------
23
24
25//FJSTARTHEADER
26// $Id$
27//
28// Copyright (c) 2005-2023, Matteo Cacciari, Gavin P. Salam and Gregory Soyez
29//
30//----------------------------------------------------------------------
31// This file is part of FastJet.
32//
33// FastJet is free software; you can redistribute it and/or modify
34// it under the terms of the GNU General Public License as published by
35// the Free Software Foundation; either version 2 of the License, or
36// (at your option) any later version.
37//
38// The algorithms that underlie FastJet have required considerable
39// development. They are described in the original FastJet paper,
40// hep-ph/0512210 and in the manual, arXiv:1111.6097. If you use
41// FastJet as part of work towards a scientific publication, please
42// quote the version you use and include a citation to the manual and
43// optionally also to hep-ph/0512210.
44//
45// FastJet is distributed in the hope that it will be useful,
46// but WITHOUT ANY WARRANTY; without even the implied warranty of
47// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
48// GNU General Public License for more details.
49//
50// You should have received a copy of the GNU General Public License
51// along with FastJet. If not, see <http://www.gnu.org/licenses/>.
52//----------------------------------------------------------------------
53//FJENDHEADER
54
55#include <iostream> // needed for io
56#include <sstream> // needed for internal io
57#include <iomanip>
58#include <cmath>
59
60#include <fastjet/ClusterSequence.hh>
61#include <fastjet/tools/MassDropTagger.hh>
62#include <fastjet/tools/Filter.hh>
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 &, const PseudoJet &);
110
111
112//----------------------------------------------------------------------
113// core of the program
114//----------------------------------------------------------------------
115int main(){
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 jet tagging using a mass drop tagger
159 //
160 // Note: if you prefer, you may as well use a CASubJetTagger
161 // CASubJetTagger ca_tagger;
162 // PseudoJet tagged = ca_tagger(jets[0]);
163 // This requires including fastjet/tools/CASubJetTagger.hh
164 // You also need to adapt the 2 lines below accessing
165 // the extra structural information provided by the tagger
166 //----------------------------------------------------------
167 MassDropTagger md_tagger(0.667, 0.09);
168 PseudoJet tagged = md_tagger(jets[0]);
169
170 if (tagged == 0){
171 cout << "No substructure found" << endl;
172 return 0;
173 }
174
175 PseudoJet parent1 = tagged.pieces()[0];
176 PseudoJet parent2 = tagged.pieces()[1];
177 cout << "Found suitable pair of subjets: " << endl;
178 cout << " " << parent1 << endl;
179 cout << " " << parent2 << endl;
180 cout << "Total = " << endl;
181 cout << " " << tagged << endl;
182 cout << "(mass drop = " << tagged.structure_of<MassDropTagger>().mu()
183 << ", y = " << tagged.structure_of<MassDropTagger>().y() << ")"
184 << endl << endl;
185
186 // next we "filter" it, to remove UE & pileup contamination
187 //----------------------------------------------------------
188 //
189 // [there are two ways of doing this; here we directly use the
190 // existing cluster sequence and find the exclusive subjets of
191 // this_jet (i.e. work backwards within the cs starting from
192 // this_jet); alternatively one can recluster just the
193 // constituents of the jet]
194 //
195 // first get separation between the subjets (called Rbb -- assuming
196 // it's a Higgs!)
197 //
198 // See example 11-filter.cc for another way of implementing the dynamic
199 // Rfilt used below
200 double Rbb = parent1.delta_R(parent2);
201 double Rfilt = min(Rbb/2, 0.3); // somewhat arbitrary choice
202 unsigned nfilt = 3; // number of pieces we'll take
203 cout << "Subjet separation (Rbb) = " << Rbb << ", Rfilt = " << Rfilt << endl;
204
205 Filter filter(JetDefinition(cambridge_algorithm, Rfilt, &flav_recombiner),
206 SelectorNHardest(nfilt));
207 PseudoJet filtered = filter(tagged);
208
209 // now print out the filtered jets and reconstruct total
210 // at the same time
211 const vector<PseudoJet> & filtered_pieces = filtered.pieces();
212 cout << "Filtered pieces are " << endl;
213 for (unsigned i = 0; i < nfilt && i < filtered_pieces.size(); i++) {
214 cout << " " << filtered_pieces[i] << endl;
215 }
216 cout << "Filtered total is " << endl;
217 cout << " " << filtered << endl;
218
219}
220
221
222//----------------------------------------------------------------------
223// does the actual work for printing out a jet
224//----------------------------------------------------------------------
225ostream & operator<<(ostream & ostr, const PseudoJet & jet) {
226 ostr << "pt, y, phi ="
227 << " " << setw(10) << jet.perp()
228 << " " << setw(6) << jet.rap()
229 << " " << setw(6) << jet.phi()
230 << ", mass = " << setw(10) << jet.m()
231 << ", btag = " << jet.user_index();
232 return ostr;
233}
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
Class that helps perform filtering (Butterworth, Davison, Rubin and Salam, arXiv:0802....
Definition: Filter.hh:97
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 that helps perform 2-pronged boosted tagging using the "mass-drop" technique (with asymmetry cu...
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
const TransformerType::StructureType & structure_of() const
this is a helper to access any structure created by a Transformer (that is, of type Transformer::Stru...
Definition: PseudoJet.hh:1146
double delta_R(const PseudoJet &other) const
return the cylinder (rap-phi) distance between this jet and another, .
Definition: PseudoJet.hh:214
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 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
virtual std::vector< PseudoJet > pieces() const
retrieve the pieces that make up the jet.
Definition: PseudoJet.cc:781
Selector SelectorNHardest(unsigned int n)
select the n hardest objects
Definition: Selector.cc:1074
the FastJet namespace
RecombinationScheme
The various recombination schemes.
vector< PseudoJet > sorted_by_pt(const vector< PseudoJet > &jets)
return a vector of jets sorted into decreasing kt2
Definition: PseudoJet.cc:871