FastJet 3.4.2
10-subjets.cc
Go to the documentation of this file.
1//----------------------------------------------------------------------
2/// \file
3/// \page Example10 10 - extracting subjets
4///
5/// fastjet example program to show how to access subjets;
6///
7/// See also 12-boosted_higgs.cc to see the use of subjets for
8/// identifying boosted higgs (and other objects)
9///
10/// run it with : ./10-subjets < data/single-event.dat
11///
12/// Source code: 10-subjets.cc
13//----------------------------------------------------------------------
14
15//FJSTARTHEADER
16// $Id$
17//
18// Copyright (c) 2005-2023, Matteo Cacciari, Gavin P. Salam and Gregory Soyez
19//
20//----------------------------------------------------------------------
21// This file is part of FastJet.
22//
23// FastJet is free software; you can redistribute it and/or modify
24// it under the terms of the GNU General Public License as published by
25// the Free Software Foundation; either version 2 of the License, or
26// (at your option) any later version.
27//
28// The algorithms that underlie FastJet have required considerable
29// development. They are described in the original FastJet paper,
30// hep-ph/0512210 and in the manual, arXiv:1111.6097. If you use
31// FastJet as part of work towards a scientific publication, please
32// quote the version you use and include a citation to the manual and
33// optionally also to hep-ph/0512210.
34//
35// FastJet is distributed in the hope that it will be useful,
36// but WITHOUT ANY WARRANTY; without even the implied warranty of
37// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
38// GNU General Public License for more details.
39//
40// You should have received a copy of the GNU General Public License
41// along with FastJet. If not, see <http://www.gnu.org/licenses/>.
42//----------------------------------------------------------------------
43//FJENDHEADER
44
45#include "fastjet/ClusterSequence.hh"
46#include <iostream> // needed for io
47#include <cstdio> // needed for io
48
49using namespace std;
50using namespace fastjet;
51
52int main(){
53
54 // read in input particles
55 //----------------------------------------------------------
56 vector<PseudoJet> input_particles;
57
58 double px, py , pz, E;
59 while (cin >> px >> py >> pz >> E) {
60 // create a PseudoJet with these components and put it onto
61 // back of the input_particles vector
62 input_particles.push_back(PseudoJet(px,py,pz,E));
63 }
64
65
66 // create a jet definition:
67 // for subjet studies, Cambridge/Aachen is the natural algorithm
68 //----------------------------------------------------------
69 double R = 1.0;
70 JetDefinition jet_def(cambridge_algorithm, R);
71
72
73 // run the jet clustering with the above jet definition
74 // and get the jets above 5 GeV
75 //----------------------------------------------------------
76 ClusterSequence clust_seq(input_particles, jet_def);
77 double ptmin = 6.0;
78 vector<PseudoJet> inclusive_jets = sorted_by_pt(clust_seq.inclusive_jets(ptmin));
79
80 // extract the subjets at a smaller angular scale (Rsub=0.5)
81 //
82 // This is done by using ClusterSequence::exclusive_subjets(dcut):
83 // for the Cambridge/Aachen algorithm, running with R and then
84 // asking for exclusive subjets with dcut should give the same
85 // subjets as rerunning the algorithm with R'=R*sqrt(dcut) on the
86 // jet's constituents.
87 //
88 // At the same time we output a summary of what has been done and the
89 // resulting subjets
90 //----------------------------------------------------------
91 double Rsub = 0.5;
92 double dcut = pow(Rsub/R,2);
93
94 // a "header" for the output
95 cout << "Ran " << jet_def.description() << endl;
96 cout << "Showing the jets above " << ptmin << " GeV" << endl;
97 cout << "And their subjets for Rsub = " << Rsub << endl;
98 printf("%10s %15s %15s %15s %15s\n","jet #", "rapidity", "phi", "pt", "n constituents");
99
100 // show the jets and their subjets
101 for (unsigned int i = 0; i < inclusive_jets.size(); i++) {
102 // get the subjets
103 vector<PseudoJet> subjets = sorted_by_pt(inclusive_jets[i].exclusive_subjets(dcut));
104
105 cout << endl;
106 // print the jet and its subjets
107 printf("%5u %15.8f %15.8f %15.8f %8d\n", i,
108 inclusive_jets[i].rap(), inclusive_jets[i].phi(),
109 inclusive_jets[i].perp(), int(inclusive_jets[i].constituents().size()));
110
111 for (unsigned int j=0; j<subjets.size(); j++)
112 printf(" sub%4u %15.8f %15.8f %15.8f %8u\n", j,
113 subjets[j].rap(), subjets[j].phi(),
114 subjets[j].perp(),
115 (unsigned int) subjets[j].constituents().size());
116 }
117
118 return 0;
119}
int main()
an example program showing how to use fastjet
Definition: 01-basic.cc:52
deals with clustering
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
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