FastJet 3.0beta1
|
Class that helps perform filtering/trimming on jets, and optionally subtraction (if rho > 0). More...
#include <fastjet/tools/Filter.hh>
Public Types | |
typedef FilterStructure | StructureType |
information about the associated structure type | |
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. | |
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. | |
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. | |
virtual | ~Filter () |
default dtor | |
void | set_subtractor (const Transformer *subtractor) |
Set a subtractor that is applied to all individual subjets before deciding which ones to keep. | |
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). | |
virtual std::string | description () const |
class description | |
Protected Attributes | |
std::vector< PseudoJet > | all_pieces |
Class that helps perform filtering/trimming on jets, and optionally subtraction (if rho > 0).
Though the original version was applied on Cambridge/Aachen jets, this one takes any jet (that has constituents) and reclusters it with a given algorithm. A user-provided Selector is applied to decide which of the subjets are kept to produce the filtered jet (others are discarded).
The constructor has the following arguments:
Filtering as proposed in arXiv:0802.2470 for boosted object reconstruction (and used also in arXiv:0810.1304 for dijet reconstructions) involves two parameters, the filtering radius, Rfilt, and the number of subjets you wish to keep, nfilt. To get a filter of this kind define
Filter filter(JetDefinition(cambridge_algorithm,Rfilt), SelectorNHardest(nfilt));
You apply it as follows
PseudoJet filtered_jet = filter(jet);
To get trimming defined with respect to a jet's pt, arXiv:0912.1342, you need an Rtrim to define subjets and a pt_fraction_min to decide which subjets to keep:
Filter trimmer(JetDefinition(cambridge_algorithm,Rfilt), SelectorPtFractionMin(pt_fraction_min));
You then apply it as before
PseudoJet trimmed_jet = trimmer(jet);
You can then find out which pieces were filtered or trimmed jet is made of by calling
trimmed_jet.pieces()
Trimming defined with respect to an event's effective mass can be carried out with a SelectorPtMin(...) selector.
More sophisticated filters/trimmers can easily be obtained by combining Selectors.
[MORE INFO, E.G. ON PIECES REJECTED, SHOULD FOLLOW]
If the 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.
See also 11 - use of filtering for a usage example.
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.
subjet_def | the jet definition applied to obtain the subjets |
selector | the Selector applied to compute the kept subjets |
rho | if 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
fastjet::Filter::Filter | ( | double | Rfilt, |
Selector | selector, | ||
double | rho = 0.0 |
||
) | [inline] |
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.
Rfilt_func | the filtering radius function of a PseudoJet |
void fastjet::Filter::set_subtractor | ( | const Transformer * | subtractor | ) | [inline] |
runs the filtering and sets kept and rejected to be the jets of interest (with non-zero rho, they will have been subtracted).
jet | the jet that gets filtered |
Implements fastjet::Transformer.