FastJet  3.1.3
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Modules Pages
Public Types | Public Member Functions | List of all members
fastjet::Pruner Class Reference

Transformer that prunes a jet. More...

#include <fastjet/tools/Pruner.hh>

Inheritance diagram for fastjet::Pruner:
Inheritance graph
[legend]
Collaboration diagram for fastjet::Pruner:
Collaboration graph
[legend]

Public Types

typedef PrunerStructure StructureType
 
- Public Types inherited from fastjet::Transformer
typedef PseudoJetStructureBase StructureType
 A typedef that is needed to ensure that the PseudoJet::structure_of() template function works.
 

Public Member Functions

 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 (practically equivalent to infinity) and also tries to use a recombiner based on the one in the jet definition of the particular jet being pruned. More...
 
 Pruner (const JetDefinition &jet_def, double zcut, double Rcut_factor)
 alternative ctor in which the full reclustering jet definition can be specified. More...
 
 Pruner (const JetDefinition &jet_def, const FunctionOfPseudoJet< double > *zcut_dyn, const FunctionOfPseudoJet< double > *Rcut_dyn)
 alternative ctor in which the pt-fraction cut and angular distance cut are functions of the jet being pruned. More...
 
virtual PseudoJet result (const PseudoJet &jet) const
 action on a single jet
 
virtual std::string description () const
 description
 
- Public Member Functions inherited from fastjet::Transformer
 Transformer ()
 default ctor
 
virtual ~Transformer ()
 default dtor
 
- Public Member Functions inherited from fastjet::FunctionOfPseudoJet< PseudoJet >
 FunctionOfPseudoJet ()
 default ctor
 
virtual ~FunctionOfPseudoJet ()
 default dtor (virtual to allow safe polymorphism)
 
PseudoJet operator() (const PseudoJet &pj) const
 apply the function using the "traditional" () operator. More...
 
std::vector< PseudoJetoperator() (const std::vector< PseudoJet > &pjs) const
 apply the function on a vector of PseudoJet, returning a vector of the results. More...
 

Detailed Description

Transformer that prunes a jet.

This transformer prunes a jet according to the ideas presented in arXiv:0903.5081 (S.D. Ellis, C.K. Vermilion and J.R. Walsh).

The jet's constituents are reclustered with a user-specified jet definition, with the modification that objects i and j are only recombined if at least one of the following two criteria is satisfied:

If both these criteria fail, i and j are not recombined, the harder of i and j is kept, and the softer is rejected.

Usage:

Pruner pruner(jet_def, zcut, Rcut_factor);
PseudoJet pruned_jet = pruner(jet);

The pruned_jet has a valid associated cluster sequence. In addition the subjets of the original jet that have been vetoed by pruning (i.e. have been 'pruned away') can be accessed using

vector<PseudoJet> rejected_subjets = pruned_jet.structure_of<Pruner>().rejected();

If the re-clustering happens to find more than a single inclusive jet (this should normally not happen if the radius of the jet definition used for the reclustering was set large enough), the hardest of these jets is retured as the result of the Pruner. The other jets can be accessed through

vector<PseudoJet> extra_jets = pruned_jet.structure_of<Pruner>().extra_jets();

Instead of using Rcut_factor and zcut, one can alternatively construct a Pruner by passing two (pointers to) functions of PseudoJet that dynamically compute the Rcut and zcut to be used for the jet being pruned.

When the jet being pruned has area support and explicit ghosts, the resulting pruned jet will likewise have area.

Definition at line 107 of file Pruner.hh.

Constructor & Destructor Documentation

fastjet::Pruner::Pruner ( const JetAlgorithm  jet_alg,
double  zcut,
double  Rcut_factor 
)
inline

minimal constructor, which takes a jet algorithm, sets the radius to JetDefinition::max_allowable_R (practically equivalent to infinity) and also tries to use a recombiner based on the one in the jet definition of the particular jet being pruned.

Parameters
jet_algthe jet algorithm for the internal clustering
zcutpt-fraction cut in the pruning
Rcut_factorthe angular distance cut in the pruning will be Rcut_factor * 2m/pt

Definition at line 118 of file Pruner.hh.

fastjet::Pruner::Pruner ( const JetDefinition jet_def,
double  zcut,
double  Rcut_factor 
)
inline

alternative ctor in which the full reclustering jet definition can be specified.

Parameters
jet_defthe jet definition for the internal clustering
zcutpt-fraction cut in the pruning
Rcut_factorthe angular distance cut in the pruning will be Rcut_factor * 2m/pt

Definition at line 131 of file Pruner.hh.

fastjet::Pruner::Pruner ( const JetDefinition jet_def,
const FunctionOfPseudoJet< double > *  zcut_dyn,
const FunctionOfPseudoJet< double > *  Rcut_dyn 
)

alternative ctor in which the pt-fraction cut and angular distance cut are functions of the jet being pruned.

Parameters
jet_defthe jet definition for the internal clustering
zcut_dyndynamic pt-fraction cut in the pruning
Rcut_dyndynamic angular distance cut in the pruning

Definition at line 54 of file Pruner.cc.


The documentation for this class was generated from the following files: