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 ostr <<
", subtractor: " << _subtractor->description();
70 }
else if (_rho != 0) {
71 ostr <<
", subtracting with rho = " << _rho;
85 throw Error(
"uninitialised Filter");
92 vector<PseudoJet> subjets;
94 bool ca_optimised = _set_filtered_elements(jet, subjets);
98 subjets = (*_subtractor)(subjets);
100 if (subjets.size()>0){
102 for (
unsigned int i=0;i<subjets.size();i++){
110 vector<PseudoJet> kept, rejected;
115 selector_copy.
sift(subjets, kept, rejected);
118 return _finalise(jet, kept, rejected, ca_optimised);
126 bool Filter::_set_filtered_elements(
const PseudoJet & jet,
127 vector<PseudoJet> & filtered_elements)
const {
130 if ((_Rfilt>=0) || (_Rfiltfunc))
131 recluster =
Recluster(cambridge_algorithm, (_Rfiltfunc) ? (*_Rfiltfunc)(jet) : _Rfilt, Recluster::keep_all);
133 recluster =
Recluster(_subjet_def,
false, Recluster::keep_all);
142 PseudoJet Filter::_finalise(
const PseudoJet & ,
143 vector<PseudoJet> & kept,
144 vector<PseudoJet> & rejected,
145 bool ca_optimisation_used)
const {
146 PseudoJet filtered_jet;
148 if (kept.size()+rejected.size()>0){
150 const JetDefinition::Recombiner &rec = (kept.size()>0)
151 ? *(kept[0].associated_cs()->jet_def().recombiner())
152 : *(rejected[0].associated_cs()->jet_def().recombiner());
155 filtered_jet = join<StructureType>(kept, rec);
157 filtered_jet = join<StructureType>(kept);
159 StructureType *fs = (StructureType*) filtered_jet.structure_non_const_ptr();
160 fs->_rejected = rejected;
166 if ((ca_optimisation_used) && (kept.size()+rejected.size()>0)){
167 bool has_non_explicit_ghost_area = (kept.size()>0)
168 ? (kept[0].has_area() && (!(kept[0].validated_csab()->has_explicit_ghosts())))
169 : (rejected[0].has_area() && (!(rejected[0].validated_csab()->has_explicit_ghosts())));
170 if (has_non_explicit_ghost_area)
178 FASTJET_END_NAMESPACE
Class that helps perform jet background subtraction.
bool takes_reference() const
returns true if this can be applied jet by jet
base class corresponding to errors that can be thrown by FastJet
const Selector & set_reference(const PseudoJet &reference)
set the reference jet for this Selector
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 encodes information about cuts and other selection criteria that can be applied to PseudoJ...
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 to contain pseudojets, including minimal information of use to jet-clustering routines...