50 #include <fastjet/PseudoJet.hh>    51 #include <fastjet/ClusterSequence.hh>    52 #include <fastjet/Selector.hh>    54 #include "fastjet/tools/Filter.hh"    56 #include "fastjet/tools/GridMedianBackgroundEstimator.hh"    57 #include "fastjet/ClusterSequenceArea.hh"    58 #include "fastjet/tools/Subtractor.hh"    72   DynamicRfilt(
double Rmax, 
double deltaR_factor) : _Rmax(Rmax), _deltaR_factor(deltaR_factor){}
    78     vector<PseudoJet> pieces = j.
pieces();
    79     if (pieces.size() != 2) 
return _Rmax;
    81     double deltaR = pieces[0].delta_R(pieces[1]);
    82     return min(_Rmax, _deltaR_factor * deltaR);
    86   double _Rmax, _deltaR_factor;
    93   vector<PseudoJet> input_particles;
    95   double px, py , pz, E;
    96   while (cin >> px >> py >> pz >> E) {
    99     input_particles.push_back(
PseudoJet(px,py,pz,E)); 
   113   printf(
"%5s %15s %15s %15s\n",
"jet #", 
"rapidity", 
"phi", 
"pt");
   116   for (
unsigned int i = 0; i < inclusive_jets.size(); i++) {
   117     printf(
"%5u %15.8f %15.8f %15.8f\n",
   118            i, inclusive_jets[i].rap(), inclusive_jets[i].phi(),
   119            inclusive_jets[i].perp());
   124   if (inclusive_jets.size()<3){
   125     cout << 
"Please provide an event with at least 3 jets above 5 GeV" << endl;
   134   vector<PseudoJet> candidates;
   135   candidates.push_back(inclusive_jets[0]);
   136   candidates.push_back(join(inclusive_jets[1],inclusive_jets[2]));
   141   vector<Filter> filters;
   161   double rho = bkgd.
rho();
   173   filters.push_back(filt);
   181   cout << 
"\nOriginal jets that will be filtered: " << endl;
   182   for (vector<PseudoJet>::iterator jit=candidates.begin(); jit!=candidates.end(); jit++){
   184     cout << 
"  rap = " << c.
rap() << 
", phi = " << c.
phi() << 
", pt = " << c.
perp() 
   189   for (vector<Filter>::iterator it=filters.begin(); it!=filters.end(); it++){
   191     cout << 
"\nUsing filter: " << f.
description() << endl;
   194     for (vector<PseudoJet>::iterator jit=candidates.begin(); jit!=candidates.end(); jit++){
   210       cout << 
"  rap = " << j.
rap() << 
", phi = " << j.
phi() << 
", pt = " << j.
perp() 
   211            << 
"  [kept: " << j.
pieces().size() << 
", rejected: "   212            << fj_struct.
rejected().size() << 
"]" << endl;
 Class that helps perform jet background subtraction. 
 
Selector SelectorPtFractionMin(double fraction)
select objects that carry at least a fraction "fraction" of the reference jet. 
 
vector< PseudoJet > sorted_by_pt(const vector< PseudoJet > &jets)
return a vector of jets sorted into decreasing kt2 
 
General class for user to obtain ClusterSequence with additional area information. 
 
T * get() const
get the stored pointer 
 
Selector SelectorNHardest(unsigned int n)
select the n hardest objects 
 
void set_subtractor(const FunctionOfPseudoJet< PseudoJet > *subtractor_in)
Set a subtractor that is applied to all individual subjets before deciding which ones to keep...
 
Class that helps perform filtering (Butterworth, Davison, Rubin and Salam, arXiv:0802.2470) and trimming (Krohn, Thaler and Wang, arXiv:0912.1342) on jets, optionally in conjunction with subtraction (Cacciari and Salam, arXiv:0707.1378). 
 
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. 
 
the longitudinally invariant kt algorithm 
 
class that holds a generic area definition 
 
virtual bool has_pieces() const
returns true if a jet has pieces 
 
double perp() const
returns the scalar transverse momentum 
 
virtual std::vector< PseudoJet > pieces() const
retrieve the pieces that make up the jet. 
 
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...
 
an implementation of C++0x shared pointers (or boost's) 
 
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 rap() const
returns the rapidity or some large value when the rapidity is infinite 
 
base class providing interface for a generic function of a PseudoJet 
 
int main()
an example program showing how to use Filter in FastJet 
 
the longitudinally invariant variant of the cambridge algorithm (aka Aachen algoithm). 
 
bool has_structure_of() const
check if the PseudoJet has the structure resulting from a Transformer (that is, its structure is comp...
 
std::string description() const
return a string describing what kind of PseudoJet we are dealing with 
 
Class to contain pseudojets, including minimal information of use to jet-clustering routines...
 
virtual std::string description() const
class description 
 
class that is intended to hold a full definition of the jet clusterer