1#ifndef __FASTJET_TOOLS_PRUNER_HH__
2#define __FASTJET_TOOLS_PRUNER_HH__
34#include "fastjet/ClusterSequence.hh"
35#include "fastjet/WrappedStructure.hh"
36#include "fastjet/tools/Transformer.hh"
40FASTJET_BEGIN_NAMESPACE
45class PruningRecombiner;
51#define FASTJET_PRUNER_STRUCTURE_STORES_RCUT
120 _zcut(zcut), _Rcut_factor(Rcut_factor),
121 _zcut_dyn(0), _Rcut_dyn(0), _get_recombiner_from_jet(true) {}
133 _zcut(zcut), _Rcut_factor(Rcut_factor),
134 _zcut_dyn(0), _Rcut_dyn(0), _get_recombiner_from_jet(false) {}
151 virtual std::string description()
const;
159 bool _check_explicit_ghosts(
const PseudoJet &jet)
const;
167 bool _check_common_recombiner(
const PseudoJet &jet,
169 bool assigned=
false)
const;
176 bool _get_recombiner_from_jet;
197 virtual std::string
description()
const{
return "Pruned PseudoJet";}
201 return validated_cs()->childless_pseudojets();
207 std::vector<PseudoJet> extra_jets()
const;
210 double Rcut()
const {
return _Rcut;}
213 double zcut()
const {
return _zcut;}
250 : _zcut2(zcut*zcut), _Rcut2(Rcut*Rcut),
251 _recombiner(recombiner){}
255 virtual void recombine(
const PseudoJet &pa,
260 virtual std::string description()
const;
263 const std::vector<unsigned int> &
rejected()
const{
return _rejected;}
276 mutable std::vector<unsigned int> _rejected;
303 : _jet_def(jet_def), _zcut(zcut), _Rcut(Rcut){}
309 virtual std::string description()
const;
312 virtual double R()
const {
return _jet_def.R();}
a class that allows a user to introduce their own "plugin" jet finder
An abstract base class that will provide the recombination scheme facilities and/or allow a user to e...
class that is intended to hold a full definition of the jet clusterer
The structure associated with a PseudoJet thas has gone through a Pruner transformer.
double zcut() const
return the value of Rcut that was used for this specific pruning.
PrunerStructure(const PseudoJet &result_jet)
default ctor
double Rcut() const
return the value of Rcut that was used for this specific pruning.
virtual std::string description() const
description
std::vector< PseudoJet > rejected() const
return the constituents that have been rejected
Transformer that prunes a jet.
Pruner(const JetAlgorithm jet_alg, double zcut, double Rcut_factor)
minimal constructor, which takes a jet algorithm, sets the radius to JetDefinition::max_allowable_R (...
Pruner(const JetDefinition &jet_def, double zcut, double Rcut_factor)
alternative ctor in which the full reclustering jet definition can be specified.
virtual double R() const
returns the radius
PruningPlugin(const JetDefinition &jet_def, double zcut, double Rcut)
ctor
const std::vector< unsigned int > & rejected() const
return the history indices that have been pruned away
PruningRecombiner(double zcut, double Rcut, const JetDefinition::Recombiner *recombiner)
ctor
void clear_rejected()
clears the list of rejected indices
Class to contain pseudojets, including minimal information of use to jet-clustering routines.
This wraps a (shared) pointer to an underlying structure.
JetAlgorithm
the various families of jet-clustering algorithm