46#include <fastjet/PseudoJet.hh>
47#include <fastjet/ClusterSequence.hh>
48#include <fastjet/Selector.hh>
50#include "fastjet/tools/Filter.hh"
51#include "fastjet/tools/Pruner.hh"
62 vector<PseudoJet> input_particles;
64 double px, py , pz, E;
65 while (cin >> px >> py >> pz >> E) {
68 input_particles.push_back(
PseudoJet(px,py,pz,E));
75 vector<PseudoJet> inclusive_jets =
79 printf(
"%5s %15s %15s %15s %15s\n",
"jet #",
"rapidity",
"phi",
"pt",
"mass");
82 for (
unsigned int i = 0; i < inclusive_jets.size(); i++) {
83 printf(
"%5u %15.8f %15.8f %15.8f %15.8f\n",
84 i, inclusive_jets[i].rap(), inclusive_jets[i].phi(),
85 inclusive_jets[i].perp(),inclusive_jets[i].m());
90 if (inclusive_jets.size()<2){
91 cout <<
"Please provide an event with at least 2 jets above 5 GeV" << endl;
101 vector<Transformer *> groomers;
107 unsigned int nfilt = 3;
114 double ptfrac = 0.03;
121 double rcut_factor = 0.5;
129 cout <<
"\nOriginal jets that will be grooomed: " << endl;
130 for (vector<PseudoJet>::iterator jit=candidates.begin(); jit!=candidates.end(); jit++){
132 cout <<
" rap = " << c.
rap() <<
", phi = " << c.
phi() <<
", pt = " << c.
perp()
133 <<
", mass = " << c.
m()
138 for (
unsigned int i=0; i < groomers.size(); i++){
140 cout <<
"\nUsing groomer: " << f.
description() << endl;
143 for (vector<PseudoJet>::iterator jit=candidates.begin(); jit!=candidates.end(); jit++){
158 n_rejected = fj_struct.
rejected().size();
162 n_rejected = fj_struct.
rejected().size();
166 cout <<
" rap = " << j.
rap() <<
", phi = " << j.
phi() <<
", pt = " << j.
perp()
167 <<
" mass = " << j.
m()
168 <<
" [kept: " << j.
pieces().size()
169 <<
", rejected: " << n_rejected;
178 for (
unsigned int i=0; i < groomers.size(); i++)
delete groomers[i];
int main()
an example program showing how to use Filter and Pruner in FastJet
std::vector< PseudoJet > inclusive_jets(const double ptmin=0.0) const
return a vector of all jets (in the sense of the inclusive algorithm) with pt >= ptmin.
Class to contain structure information for a filtered jet.
const std::vector< PseudoJet > & rejected() const
returns the subjets that were not kept during the filtering procedure (subtracted if the filter reque...
Class that helps perform filtering (Butterworth, Davison, Rubin and Salam, arXiv:0802....
class that is intended to hold a full definition of the jet clusterer
The structure associated with a PseudoJet thas has gone through a Pruner transformer.
std::vector< PseudoJet > rejected() const
return the constituents that have been rejected
Transformer that prunes a jet.
Class to contain pseudojets, including minimal information of use to jet-clustering routines.
double rap() const
returns the rapidity or some large value when the rapidity is infinite
const TransformerType::StructureType & structure_of() const
this is a helper to access any structure created by a Transformer (that is, of type Transformer::Stru...
double phi() const
returns phi (in the range 0..2pi)
double perp() const
returns the scalar transverse momentum
bool has_structure_of() const
check if the PseudoJet has the structure resulting from a Transformer (that is, its structure is comp...
double m() const
returns the invariant mass (If m2() is negative then -sqrt(-m2()) is returned, as in CLHEP)
std::string description() const
return a string describing what kind of PseudoJet we are dealing with
virtual std::vector< PseudoJet > pieces() const
retrieve the pieces that make up the jet.
Selector SelectorNHardest(unsigned int n)
select the n hardest objects
Selector SelectorPtFractionMin(double fraction)
select objects that carry at least a fraction "fraction" of the reference jet.
@ cambridge_algorithm
the longitudinally invariant variant of the cambridge algorithm (aka Aachen algoithm).
@ kt_algorithm
the longitudinally invariant kt algorithm
vector< PseudoJet > sorted_by_pt(const vector< PseudoJet > &jets)
return a vector of jets sorted into decreasing kt2