31 #include "fastjet/tools/Filter.hh"
32 #include "fastjet/tools/Recluster.hh"
33 #include "fastjet/tools/Subtractor.hh"
34 #include "fastjet/ClusterSequenceActiveAreaExplicitGhosts.hh"
43 FASTJET_BEGIN_NAMESPACE
50 string Filter::description()
const {
52 return "uninitialised Filter";
56 ostr <<
"Filter with subjet_def = ";
58 ostr <<
"Cambridge/Aachen algorithm with dynamic Rfilt"
59 <<
" (recomb. scheme deduced from jet, or E-scheme if not unique)";
60 }
else if (_Rfilt > 0) {
61 ostr <<
"Cambridge/Aachen algorithm with Rfilt = "
63 <<
" (recomb. scheme deduced from jet, or E-scheme if not unique)";
65 ostr << _subjet_def.description();
67 ostr<<
", selection " << _selector.description();
69 if (local_subtractor) {
70 ostr <<
", subtractor: " << local_subtractor->
description();
71 }
else if (_rho != 0) {
72 ostr <<
", subtracting with rho = " << _rho;
86 throw Error(
"uninitialised Filter");
93 vector<PseudoJet> subjets;
95 bool ca_optimised = _set_filtered_elements(jet, subjets);
99 if (local_subtractor){
100 subjets = (*local_subtractor)(subjets);
102 if (subjets.size()>0){
104 for (
unsigned int i=0;i<subjets.size();i++){
112 vector<PseudoJet> kept, rejected;
117 selector_copy.
sift(subjets, kept, rejected);
120 return _finalise(jet, kept, rejected, ca_optimised);
128 bool Filter::_set_filtered_elements(
const PseudoJet & jet,
129 vector<PseudoJet> & filtered_elements)
const {
132 if ((_Rfilt>=0) || (_Rfiltfunc))
133 recluster =
Recluster(cambridge_algorithm, (_Rfiltfunc) ? (*_Rfiltfunc)(jet) : _Rfilt, Recluster::keep_all);
135 recluster =
Recluster(_subjet_def,
false, Recluster::keep_all);
144 PseudoJet Filter::_finalise(
const PseudoJet & ,
145 vector<PseudoJet> & kept,
146 vector<PseudoJet> & rejected,
147 bool ca_optimisation_used)
const {
148 PseudoJet filtered_jet;
150 if (kept.size()+rejected.size()>0){
152 const JetDefinition::Recombiner &rec = (kept.size()>0)
153 ? *(kept[0].associated_cs()->jet_def().recombiner())
154 : *(rejected[0].associated_cs()->jet_def().recombiner());
157 filtered_jet = join<StructureType>(kept, rec);
159 filtered_jet = join<StructureType>(kept);
161 StructureType *fs = (StructureType*) filtered_jet.structure_non_const_ptr();
162 fs->_rejected = rejected;
168 if ((ca_optimisation_used) && (kept.size()+rejected.size()>0)){
169 bool has_non_explicit_ghost_area = (kept.size()>0)
170 ? (kept[0].has_area() && (!(kept[0].validated_csab()->has_explicit_ghosts())))
171 : (rejected[0].has_area() && (!(rejected[0].validated_csab()->has_explicit_ghosts())));
172 if (has_non_explicit_ghost_area)
180 FASTJET_END_NAMESPACE
base class corresponding to errors that can be thrown by FastJet
virtual std::string description() const
returns a description of the function (an empty string by default)
Class to contain pseudojets, including minimal information of use to jet-clustering routines.
Recluster a jet's constituents with a new jet definition.
bool get_new_jets_and_def(const PseudoJet &input_jet, std::vector< PseudoJet > &output_jets) const
A lower-level method that does the actual work of reclustering the input jet.
Class that encodes information about cuts and other selection criteria that can be applied to PseudoJ...
const Selector & set_reference(const PseudoJet &reference)
set the reference jet for this Selector
bool takes_reference() const
returns true if this can be applied jet by jet
void sift(const std::vector< PseudoJet > &jets, std::vector< PseudoJet > &jets_that_pass, std::vector< PseudoJet > &jets_that_fail) const
sift the input jets into two vectors – those that pass the selector and those that do not
Class that helps perform jet background subtraction.