Class that helps perform filtering/trimming on jets, and optionally subtraction (if rho > 0). More...
#include <Filter.hh>
Public Types | |
typedef FilteredJetStructure | 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. | |
virtual | ~Filter () |
default dtor | |
virtual PseudoJet | operator() (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::vector< PseudoJet > | operator() (const std::vector< PseudoJet > &originals) const |
action of the transformer on each jet from the vector this has to be repeated because it shares the same name as the operator()(PseudoJet) | |
std::string | description () const |
class description | |
Protected Member Functions | |
void | _set_filtered_elements (const PseudoJet &jet, std::vector< PseudoJet > &filtered_elements, ClusterSequence *&internal_cs) const |
sets filtered_elements to be all the subjets on which filtering will work [NB: this routine is work in progress as part of a transition to a Filter that also works on jet collections] | |
PseudoJet | _finalise (const PseudoJet &jet, std::vector< PseudoJet > &kept, std::vector< PseudoJet > &rejected, ClusterSequence *&internal_cs) const |
gather the information about what is kept and rejected under the form of a PseudoJet with a special ClusterSequenceInfo | |
bool | _recursively_check_ca (const PseudoJet &jet) const |
check if the jet is obtained from C/A or a superposition of C/A pieces | |
void | _set_filtered_elements_cafilt (const PseudoJet &jet, std::vector< PseudoJet > &filtered_elements, double Rfilt) const |
set the filtered elements in the simple case of C/A+C/A | |
ClusterSequence * | _set_filtered_elements_generic_unsubtracted (const PseudoJet &jet, std::vector< PseudoJet > &filtered_elements) const |
set the filtered elements in the generic re-clustering case (wo subtraction) | |
ClusterSequence * | _set_filtered_elements_generic_subtracted (const PseudoJet &jet, std::vector< PseudoJet > &filtered_elements) const |
set the filtered elements in the generic re-clustering case (with subtraction) | |
Protected Attributes | |
JetDefinition | _subjet_def |
the jet definition to use to extract the subjets | |
Selector | _selector |
the subjet selection criterium | |
double | _rho |
the background density (used for subtraction when possible) |
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, 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()
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 12 - use of filtering for a usage example.
Definition at line 131 of file Filter.hh.
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
Definition at line 150 of file Filter.hh.
: _subjet_def(subjet_def), _selector(selector), _rho(rho) {}
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 |
Reimplemented from fastjet::Transformer.
Definition at line 58 of file Filter.cc.
References _finalise(), _selector, _set_filtered_elements(), fastjet::Selector::nullify_non_selected(), fastjet::Selector::set_reference(), and fastjet::Selector::takes_reference().
{ ClusterSequence *internal_cs; // start by getting the list of subjets (including a sanity check // that all the jets in the argument vector share the same // underlying ClusterSequence) vector<PseudoJet> subjets; // NB: empty to begin with (see the comment for _set_filtered_elements_cafilt) _set_filtered_elements(jet, subjets, internal_cs); // decide what to keep and what to reject // we first make a copy of the pointers then apply the selector vector<const PseudoJet*> subjet_pointers; for (unsigned int i=0;i<subjets.size(); i++) subjet_pointers.push_back(&(subjets[i])); // Note that the following line is the one requiring that _seelector be declared as mutable if (_selector.takes_reference()) _selector.set_reference(jet); _selector.nullify_non_selected(subjet_pointers); // now build the vector of kept and rejected subjets vector<PseudoJet> kept, rejected; for (unsigned int i=0;i<subjets.size(); i++){ if (subjet_pointers[i]==0) rejected.push_back(subjets[i]); else kept.push_back(subjets[i]); } // gather the info under the form of a PseudoJet return _finalise(jet, kept, rejected, internal_cs); }