FastJet  3.3.0
Public Types | Public Member Functions | List of all members
fastjet::Filter Class Reference

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). More...

#include <fastjet/tools/Filter.hh>

Inheritance diagram for fastjet::Filter:
Inheritance graph
[legend]
Collaboration diagram for fastjet::Filter:
Collaboration graph
[legend]

Public Types

typedef FilterStructure StructureType
 
- Public Types inherited from fastjet::Transformer
typedef PseudoJetStructureBase StructureType
 A typedef that is needed to ensure that the PseudoJet::structure_of() template function works.
 

Public Member Functions

 Filter ()
 trivial ctor Note: this is just for derived classes a Filter initialised through this constructor will not work!
 
 Filter (JetDefinition subjet_def, Selector selector, double rho=0.0)
 define a filter that decomposes a jet into subjets using a generic JetDefinition and then keeps only a subset of these subjets according to a Selector. More...
 
 Filter (double Rfilt, Selector selector, double rho=0.0)
 Same as the full constructor (see above) but just specifying the radius By default, Cambridge-Aachen is used If the jet (or all its pieces) is obtained with a non-default recombiner, that one will be used. More...
 
 Filter (FunctionOfPseudoJet< double > *Rfilt_func, Selector selector, double rho=0.0)
 Same as the full constructor (see above) but just specifying a filtering radius that will depend on the jet being filtered As for the previous case, Cambridge-Aachen is used If the jet (or all its pieces) is obtained with a non-default recombiner, that one will be used. More...
 
virtual ~Filter ()
 default dtor
 
void set_subtractor (const FunctionOfPseudoJet< PseudoJet > *subtractor_in)
 Set a subtractor that is applied to all individual subjets before deciding which ones to keep. More...
 
const FunctionOfPseudoJet< PseudoJet > * subtractor () const
 Set a subtractor that is applied to all individual subjets before deciding which ones to keep. More...
 
virtual PseudoJet result (const PseudoJet &jet) const
 runs the filtering and sets kept and rejected to be the jets of interest (with non-zero rho, they will have been subtracted). More...
 
virtual std::string description () const
 class description
 
- Public Member Functions inherited from fastjet::Transformer
 Transformer ()
 default ctor
 
virtual ~Transformer ()
 default dtor
 
- Public Member Functions inherited from fastjet::FunctionOfPseudoJet< PseudoJet >
 FunctionOfPseudoJet ()
 default ctor
 
virtual ~FunctionOfPseudoJet ()
 default dtor (virtual to allow safe polymorphism)
 
PseudoJet operator() (const PseudoJet &pj) const
 apply the function using the "traditional" () operator. More...
 
std::vector< PseudoJetoperator() (const std::vector< PseudoJet > &pjs) const
 apply the function on a vector of PseudoJet, returning a vector of the results. More...
 

Detailed Description

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).

For example, to apply filtering that reclusters a jet's constituents with the Cambridge/Aachen jet algorithm with R=0.3 and then selects the 3 hardest subjets, one can use the following code:

Filter filter(JetDefinition(cambridge_algorithm, 0.3), SelectorNHardest(3));
PseudoJet filtered_jet = filter(original_jet);

To obtain trimming, involving for example the selection of all subjets carrying at least 3% of the original jet's pt, the selector would be replaced by SelectorPtFractionMin(0.03).

To additionally perform subtraction on the subjets prior to selection, either include a 3rd argument specifying the background density rho, or call the set_subtractor(...) member function. If subtraction is requested, the original jet must be the result of a clustering with active area with explicit ghosts support or a merging of such pieces.

The information on the subjets that were kept and rejected can be obtained using:

vector<PseudoJet> kept_subjets = filtered_jet.pieces();
vector<PseudoJet> rejected_subjets = filtered_jet.structure_of<Filter>().rejected();

Implementation Note

If the original jet was defined with the Cambridge/Aachen algorithm (or is made of pieces each of which comes from the C/A alg) and the filtering definition is C/A, then the filter does not rerun the C/A algorithm on the constituents, but instead makes use of the existent C/A cluster sequence in the original jet. This increases the speed of the filter.

See also 11 - use of filtering for a further usage example.

Support for areas, reuse of C/A cluster sequences, etc., considerably complicates the implementation of Filter. For an explanation of how a simpler filter might be coded, see the "User-defined transformers" appendix of the manual.

Definition at line 97 of file Filter.hh.

Constructor & Destructor Documentation

◆ Filter() [1/3]

fastjet::Filter::Filter ( JetDefinition  subjet_def,
Selector  selector,
double  rho = 0.0 
)
inline

define a filter that decomposes a jet into subjets using a generic JetDefinition and then keeps only a subset of these subjets according to a Selector.

Optionally, each subjet may be internally bakground-subtracted prior to selection.

Parameters
subjet_defthe jet definition applied to obtain the subjets
selectorthe Selector applied to compute the kept subjets
rhoif non-zero, backgruond-subtract each subjet befor selection

Note: internal subtraction only applies on jets that are obtained with a cluster sequence with area support and explicit ghosts

Definition at line 116 of file Filter.hh.

◆ Filter() [2/3]

fastjet::Filter::Filter ( double  Rfilt,
Selector  selector,
double  rho = 0.0 
)
inline

Same as the full constructor (see above) but just specifying the radius By default, Cambridge-Aachen is used If the jet (or all its pieces) is obtained with a non-default recombiner, that one will be used.

Parameters
Rfiltthe filtering radius

Definition at line 124 of file Filter.hh.

◆ Filter() [3/3]

fastjet::Filter::Filter ( FunctionOfPseudoJet< double > *  Rfilt_func,
Selector  selector,
double  rho = 0.0 
)
inline

Same as the full constructor (see above) but just specifying a filtering radius that will depend on the jet being filtered As for the previous case, Cambridge-Aachen is used If the jet (or all its pieces) is obtained with a non-default recombiner, that one will be used.

Parameters
Rfilt_functhe filtering radius function of a PseudoJet

Definition at line 136 of file Filter.hh.

Member Function Documentation

◆ set_subtractor()

void fastjet::Filter::set_subtractor ( const FunctionOfPseudoJet< PseudoJet > *  subtractor_in)
inline

Set a subtractor that is applied to all individual subjets before deciding which ones to keep.

It takes precedence over a non-zero rho.

Definition at line 144 of file Filter.hh.

◆ subtractor()

const FunctionOfPseudoJet<PseudoJet>* fastjet::Filter::subtractor ( ) const
inline

Set a subtractor that is applied to all individual subjets before deciding which ones to keep.

It takes precedence over a non-zero rho.

Definition at line 148 of file Filter.hh.

◆ result()

PseudoJet fastjet::Filter::result ( const PseudoJet jet) const
virtual

runs the filtering and sets kept and rejected to be the jets of interest (with non-zero rho, they will have been subtracted).

Parameters
jetthe jet that gets filtered
Returns
the filtered jet

Implements fastjet::Transformer.

Definition at line 82 of file Filter.cc.


The documentation for this class was generated from the following files: