52#include <fastjet/PseudoJet.hh>
53#include <fastjet/ClusterSequence.hh>
54#include <fastjet/Selector.hh>
56#include "fastjet/tools/Filter.hh"
58#include "fastjet/tools/GridMedianBackgroundEstimator.hh"
59#include "fastjet/ClusterSequenceArea.hh"
60#include "fastjet/tools/Subtractor.hh"
74 DynamicRfilt(
double Rmax,
double deltaR_factor) : _Rmax(Rmax), _deltaR_factor(deltaR_factor){}
80 vector<PseudoJet> pieces = j.
pieces();
81 if (pieces.size() != 2)
return _Rmax;
83 double deltaR = pieces[0].delta_R(pieces[1]);
84 return min(_Rmax, _deltaR_factor * deltaR);
88 double _Rmax, _deltaR_factor;
95 vector<PseudoJet> input_particles;
97 double px, py , pz, E;
98 while (cin >> px >> py >> pz >> E) {
101 input_particles.push_back(
PseudoJet(px,py,pz,E));
115 printf(
"%5s %15s %15s %15s\n",
"jet #",
"rapidity",
"phi",
"pt");
118 for (
unsigned int i = 0; i < inclusive_jets.size(); i++) {
119 printf(
"%5u %15.8f %15.8f %15.8f\n",
120 i, inclusive_jets[i].rap(), inclusive_jets[i].phi(),
121 inclusive_jets[i].perp());
126 if (inclusive_jets.size()<3){
127 cout <<
"Please provide an event with at least 3 jets above 5 GeV" << endl;
136 vector<PseudoJet> candidates;
137 candidates.push_back(inclusive_jets[0]);
138 candidates.push_back(join(inclusive_jets[1],inclusive_jets[2]));
143 vector<Filter> filters;
163 double rho = bkgd.
rho();
175 filters.push_back(filt);
183 cout <<
"\nOriginal jets that will be filtered: " << endl;
184 for (vector<PseudoJet>::iterator jit=candidates.begin(); jit!=candidates.end(); jit++){
186 cout <<
" rap = " << c.
rap() <<
", phi = " << c.
phi() <<
", pt = " << c.
perp()
191 for (vector<Filter>::iterator it=filters.begin(); it!=filters.end(); it++){
193 cout <<
"\nUsing filter: " << f.
description() << endl;
196 for (vector<PseudoJet>::iterator jit=candidates.begin(); jit!=candidates.end(); jit++){
212 cout <<
" rap = " << j.
rap() <<
", phi = " << j.
phi() <<
", pt = " << j.
perp()
213 <<
" [kept: " << j.
pieces().size() <<
", rejected: "
214 << fj_struct.
rejected().size() <<
"]" << endl;
int main()
an example program showing how to use Filter in FastJet
class that holds a generic area definition
General class for user to obtain ClusterSequence with additional area information.
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....
void set_subtractor(const FunctionOfPseudoJet< PseudoJet > *subtractor_in)
Set a subtractor that is applied to all individual subjets before deciding which ones to keep.
virtual std::string description() const
class description
base class providing interface for a generic function of a PseudoJet
virtual TOut result(const PseudoJet &pj) const =0
the action of the function this has to be overloaded in derived classes
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.
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...
virtual bool has_pieces() const
returns true if a jet has pieces
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.
An implementation of shared pointers that is broadly similar to C++11 shared_ptr (https://en....
T * get() const
get the stored pointer
Class that helps perform jet background subtraction.
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