FastJet 3.4.2
13-boosted_top.cc
Go to the documentation of this file.
1//----------------------------------------------------------------------
2/// \file
3/// \page Example13 13 - boosted top tagging
4///
5/// fastjet example program, illustration of carrying out boosted
6/// top subjet ID analysis using the Johns Hopkins top tagger as
7/// introduced in arXiv:0806.0848 (Kaplan, Rehermann, Schwartz
8/// and Tweedie)
9///
10/// run it with : ./13-boosted_top < data/boosted_top_event.dat
11///
12/// Source code: 13-boosted_top.cc
13//----------------------------------------------------------------------
14
15
16//FJSTARTHEADER
17// $Id$
18//
19// Copyright (c) 2005-2023, Matteo Cacciari, Gavin P. Salam and Gregory Soyez
20//
21//----------------------------------------------------------------------
22// This file is part of FastJet.
23//
24// FastJet is free software; you can redistribute it and/or modify
25// it under the terms of the GNU General Public License as published by
26// the Free Software Foundation; either version 2 of the License, or
27// (at your option) any later version.
28//
29// The algorithms that underlie FastJet have required considerable
30// development. They are described in the original FastJet paper,
31// hep-ph/0512210 and in the manual, arXiv:1111.6097. If you use
32// FastJet as part of work towards a scientific publication, please
33// quote the version you use and include a citation to the manual and
34// optionally also to hep-ph/0512210.
35//
36// FastJet is distributed in the hope that it will be useful,
37// but WITHOUT ANY WARRANTY; without even the implied warranty of
38// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
39// GNU General Public License for more details.
40//
41// You should have received a copy of the GNU General Public License
42// along with FastJet. If not, see <http://www.gnu.org/licenses/>.
43//----------------------------------------------------------------------
44//FJENDHEADER
45
46#include <iostream> // needed for io
47#include <sstream> // needed for internal io
48#include <iomanip>
49#include <cmath>
50
51#include <fastjet/ClusterSequence.hh>
52#include <fastjet/Selector.hh>
53#include <fastjet/tools/JHTopTagger.hh>
54
55using namespace std;
56using namespace fastjet;
57
58
59//----------------------------------------------------------------------
60// forward declaration for printing out info about a jet
61//----------------------------------------------------------------------
62ostream & operator<<(ostream &, const PseudoJet &);
63
64//----------------------------------------------------------------------
65// core of the program
66//----------------------------------------------------------------------
67int main(){
68
69 vector<PseudoJet> particles;
70
71 // read in data in format px py pz E b-tag [last of these is optional]
72 // lines starting with "#" are considered as comments and discarded
73 //----------------------------------------------------------
74 string line;
75 while (getline(cin,line)) {
76 if (line.substr(0,1) == "#") {continue;}
77 istringstream linestream(line);
78 double px,py,pz,E;
79 linestream >> px >> py >> pz >> E;
80
81 // construct the particle
82 particles.push_back(PseudoJet(px,py,pz,E));
83 }
84
85 // compute the parameters to be used through the analysis
86 // ----------------------------------------------------------
87 double Et=0;
88 for (unsigned int i=0; i<particles.size(); i++)
89 Et += particles[i].perp();
90
91 double R, delta_p, delta_r;
92 if (Et>2600){ R=0.4; delta_p=0.05; delta_r=0.19;}
93 else if (Et>1600){ R=0.6; delta_p=0.05; delta_r=0.19;}
94 else if (Et>1000){ R=0.8; delta_p=0.10; delta_r=0.19;}
95 else{ cerr << "Et has to be at least 1 TeV"<< endl; return 1;}
96
97 double ptmin = min(500.0, 0.7*Et/2);
98
99 // find the jets
100 // ----------------------------------------------------------
101 JetDefinition jet_def(cambridge_algorithm, R);
102 ClusterSequence cs(particles, jet_def);
103 vector<PseudoJet> jets = sorted_by_pt(cs.inclusive_jets());
104
105 cout << "Ran: " << jet_def.description() << endl << endl;
106 cout << "2 Hardest jets: " << jets[0] << endl
107 << " " << jets[1] << endl << endl;
108
109 if (jets[0].perp()<ptmin){
110 cout << "No jet above the ptmin threshold" << endl;
111 return 2;
112 }
113
114 // now do jet tagging using the Johns Hopkins top tagger
115 // For simplicity, we just apply it to the hardest jet.
116 //
117 // In addition to delta_p and delta_r, note that there are two
118 // further parameters to the JH top tagger that here are implicitly
119 // set to their defaults:
120 //
121 // - cos_theta_W_max (defaults to 0.7)
122 // - mW (defaults to 80.4).
123 //
124 // The value for mW implicitly assumes that momenta are passed in
125 // GeV.
126 // ----------------------------------------------------------
127 JHTopTagger top_tagger(delta_p, delta_r);
128 top_tagger.set_top_selector(SelectorMassRange(150,200));
129 top_tagger.set_W_selector (SelectorMassRange( 65, 95));
130
131 PseudoJet tagged = top_tagger(jets[0]);
132
133 cout << "Ran the following top tagger: " << top_tagger.description() << endl;
134
135 if (tagged == 0){
136 cout << "No top substructure found" << endl;
137 return 0;
138 }
139
140 cout << "Found top substructure from the hardest jet:" << endl;
141 cout << " top candidate: " << tagged << endl;
142 cout << " |_ W candidate: " << tagged.structure_of<JHTopTagger>().W() << endl;
143 cout << " | |_ W subjet 1: " << tagged.structure_of<JHTopTagger>().W1() << endl;
144 cout << " | |_ W subjet 2: " << tagged.structure_of<JHTopTagger>().W2() << endl;
145 cout << " | cos(theta_W) = " << tagged.structure_of<JHTopTagger>().cos_theta_W() << endl;
146 cout << " |_ non-W subjet: " << tagged.structure_of<JHTopTagger>().non_W() << endl;
147}
148
149
150//----------------------------------------------------------------------
151// does the actual work for printing out a jet
152//----------------------------------------------------------------------
153ostream & operator<<(ostream & ostr, const PseudoJet & jet) {
154 ostr << "pt, y, phi =" << setprecision(6)
155 << " " << setw(9) << jet.perp()
156 << " " << setw(9) << jet.rap()
157 << " " << setw(9) << jet.phi()
158 << ", mass = " << setw(9) << jet.m();
159 return ostr;
160}
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 boosted top tagging using the "Johns Hopkins" method from arXiv:0806....
Definition: JHTopTagger.hh:120
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
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 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
Selector SelectorMassRange(double mmin, double mmax)
select objects with Mmin <= Mass <= Mmax
Definition: Selector.cc:766
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