FastJet 3.4.1
Filter.cc
1//FJSTARTHEADER
2// $Id$
3//
4// Copyright (c) 2005-2023, Matteo Cacciari, Gavin P. Salam and Gregory Soyez
5//
6//----------------------------------------------------------------------
7// This file is part of FastJet.
8//
9// FastJet is free software; you can redistribute it and/or modify
10// it under the terms of the GNU General Public License as published by
11// the Free Software Foundation; either version 2 of the License, or
12// (at your option) any later version.
13//
14// The algorithms that underlie FastJet have required considerable
15// development. They are described in the original FastJet paper,
16// hep-ph/0512210 and in the manual, arXiv:1111.6097. If you use
17// FastJet as part of work towards a scientific publication, please
18// quote the version you use and include a citation to the manual and
19// optionally also to hep-ph/0512210.
20//
21// FastJet is distributed in the hope that it will be useful,
22// but WITHOUT ANY WARRANTY; without even the implied warranty of
23// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
24// GNU General Public License for more details.
25//
26// You should have received a copy of the GNU General Public License
27// along with FastJet. If not, see <http://www.gnu.org/licenses/>.
28//----------------------------------------------------------------------
29//FJENDHEADER
30
31#include "fastjet/tools/Filter.hh"
32#include "fastjet/tools/Recluster.hh"
33#include "fastjet/tools/Subtractor.hh"
34#include "fastjet/ClusterSequenceActiveAreaExplicitGhosts.hh"
35#include <cassert>
36#include <algorithm>
37#include <sstream>
38#include <typeinfo>
39
40using namespace std;
41
42
43FASTJET_BEGIN_NAMESPACE // defined in fastjet/internal/base.hh
44
45//----------------------------------------------------------------------
46// Filter class implementation
47//----------------------------------------------------------------------
48
49// class description
50string Filter::description() const {
51 if (!_initialised){
52 return "uninitialised Filter";
53 }
54
55 ostringstream ostr;
56 ostr << "Filter with subjet_def = ";
57 if (_Rfiltfunc) {
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 = "
62 << _Rfilt
63 << " (recomb. scheme deduced from jet, or E-scheme if not unique)";
64 } else {
65 ostr << _subjet_def.description();
66 }
67 ostr<< ", selection " << _selector.description();
68 const FunctionOfPseudoJet<PseudoJet> *local_subtractor = _subtractor;
69 if (local_subtractor) {
70 ostr << ", subtractor: " << local_subtractor->description();
71 } else if (_rho != 0) {
72 ostr << ", subtracting with rho = " << _rho;
73 }
74 return ostr.str();
75}
76
77
78// core functions
79//----------------------------------------------------------------------
80
81// return a vector of subjets, which are the ones that would be kept
82// by the filtering
83PseudoJet Filter::result(const PseudoJet &jet) const {
84 if (!_initialised){
85 //Q: do we throw or do we return an empty PJ?
86 throw Error("uninitialised Filter");
87 }
88
89 // start by getting the list of subjets (including a list of sanity
90 // checks)
91 // NB: subjets is empty to begin with (see the comment for
92 // _set_filtered_elements_cafilt)
93 vector<PseudoJet> subjets;
94 //JetDefinition subjet_def;
95 bool ca_optimised = _set_filtered_elements(jet, subjets);
96
97 // apply subtraction if needed:
98 const FunctionOfPseudoJet<PseudoJet> *local_subtractor = _subtractor;
99 if (local_subtractor){
100 subjets = (*local_subtractor)(subjets);
101 } else if (_rho!=0){
102 if (subjets.size()>0){
103 //const ClusterSequenceAreaBase *csab = subjets[0].validated_csab();
104 for (unsigned int i=0;i<subjets.size();i++){
105 //subjets[i]=csab->subtracted_jet(subjets[i], _rho);
106 subjets[i]=Subtractor(_rho)(subjets[i]);
107 }
108 }
109 }
110
111 // now build the vector of kept and rejected subjets
112 vector<PseudoJet> kept, rejected;
113 // Note that in the following line we make a copy of the _selector
114 // to avoid issues with needing a mutable _selector
115 Selector selector_copy = _selector;
116 if (selector_copy.takes_reference()) selector_copy.set_reference(jet);
117 selector_copy.sift(subjets, kept, rejected);
118
119 // gather the info under the form of a PseudoJet
120 return _finalise(jet, kept, rejected, ca_optimised);
121}
122
123
124// sets filtered_elements to be all the subjets on which filtering will work
125//
126// return true when the subjets have been optained using teh optimised
127// method for C/A
128bool Filter::_set_filtered_elements(const PseudoJet & jet,
129 vector<PseudoJet> & filtered_elements) const {
130 // create the recluster instance
131 Recluster recluster;
132 if ((_Rfilt>=0) || (_Rfiltfunc))
133 recluster = Recluster(cambridge_algorithm, (_Rfiltfunc) ? (*_Rfiltfunc)(jet) : _Rfilt, Recluster::keep_all);
134 else
135 recluster = Recluster(_subjet_def, false, Recluster::keep_all);
136
137 // get the subjets
138 //JetDefinition subjet_def;
139 return recluster.get_new_jets_and_def(jet, filtered_elements);
140}
141
142// gather the information about what is kept and rejected under the
143// form of a PseudoJet with a special ClusterSequenceInfo
144PseudoJet Filter::_finalise(const PseudoJet & /*jet*/,
145 vector<PseudoJet> & kept,
146 vector<PseudoJet> & rejected,
147 bool ca_optimisation_used) const {
148 PseudoJet filtered_jet;
149
150 if (kept.size()+rejected.size()>0){
151 // figure out which recombiner to use
152 const JetDefinition::Recombiner &rec = (kept.size()>0)
153 ? *(kept[0].associated_cs()->jet_def().recombiner())
154 : *(rejected[0].associated_cs()->jet_def().recombiner());
155
156 // create an appropriate structure and transfer the info to it
157 filtered_jet = join<StructureType>(kept, rec);
158 } else {
159 filtered_jet = join<StructureType>(kept);
160 }
161 StructureType *fs = (StructureType*) filtered_jet.structure_non_const_ptr();
162 fs->_rejected = rejected;
163
164 // if we've used C/A optimisation, we need to get rid of the area
165 // information if it comes from a non-explicit-ghost clustering.
166 // (because in that case it can be erroneous due the lack of
167 // information about empty areas)
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)
173 fs->discard_area();
174 }
175
176 return filtered_jet;
177}
178
179
180FASTJET_END_NAMESPACE // defined in fastjet/internal/base.hh
base class corresponding to errors that can be thrown by FastJet
Definition: Error.hh:52
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.
Definition: PseudoJet.hh:68
Recluster a jet's constituents with a new jet definition.
Definition: Recluster.hh:73
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.
Definition: Recluster.cc:108
Class that encodes information about cuts and other selection criteria that can be applied to PseudoJ...
Definition: Selector.hh:149
bool takes_reference() const
returns true if this can be applied jet by jet
Definition: Selector.hh:287
const Selector & set_reference(const PseudoJet &reference)
set the reference jet for this Selector
Definition: Selector.hh:292
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
Definition: Selector.cc:153
Class that helps perform jet background subtraction.
Definition: Subtractor.hh:62