FastJet
3.3.3
|
clean (almost parameter-free) tagger searching for the element in the clustering history that maximises a chosen distance More...
#include <fastjet/tools/CASubJetTagger.hh>
Classes | |
class | JetAux |
class that contains the result internally More... | |
Public Types | |
enum | ScaleChoice { kt2_distance, jade_distance, jade2_distance, plain_distance, mass_drop_distance, dot_product_distance } |
the available choices of auxiliary scale with respect to which to order the splittings | |
typedef CASubJetTaggerStructure | StructureType |
the type of Structure returned | |
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 | |
CASubJetTagger (ScaleChoice scale_choice=jade_distance, double z_threshold=0.1) | |
just constructs | |
void | set_dr_min (double drmin) |
sets a minimum delta R below which spliting will be ignored (only relevant if set prior to calling run()) | |
virtual std::string | description () const |
returns a textual description of the tagger | |
void | set_absolute_z_cut (bool abs_z_cut=true) |
If (abs_z_cut) is set to false (the default) then for a splitting to be considered, each subjet must satisfy. More... | |
virtual PseudoJet | result (const PseudoJet &jet) const |
runs the tagger on the given jet and returns the tagged PseudoJet if successful, or a PseudoJet==0 otherwise (standard access is through operator()). More... | |
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< PseudoJet > | operator() (const std::vector< PseudoJet > &pjs) const |
apply the function on a vector of PseudoJet, returning a vector of the results. More... | |
Protected Member Functions | |
void | _recurse_through_jet (const PseudoJet ¤t_jet, JetAux &aux_max, const PseudoJet &original_jet) const |
work through the jet, establishing a distance at each branching More... | |
Protected Attributes | |
ScaleChoice | _scale_choice |
double | _z_threshold |
double | _dr2_min |
bool | _absolute_z_cut |
Static Protected Attributes | |
static LimitedWarning | _non_ca_warnings |
clean (almost parameter-free) tagger searching for the element in the clustering history that maximises a chosen distance
class to help us get a clean (almost parameter-free) handle on substructure inside a C/A jet. It follows the logic described in arXiv:0906.0728 (and is inspired by the original Cambridge algorithm paper in its use of separate angular and dimensionful distances), but provides some extra flexibility.
It searches for all splittings that pass a symmetry cut (zcut) and then selects the one with the largest auxiliary scale choice (e.g. jade distance of the splitting, kt distance of the splitting, etc.)
By default, the zcut is calculated from the fraction of the child pt carried by the parent jet. If one calls set_absolute_z_cut the fraction of transverse momentum will be computed wrt the original jet.
original code copyright (C) 2009 by Gavin Salam, released under the GPL.
Definition at line 106 of file CASubJetTagger.hh.
|
inline |
If (abs_z_cut) is set to false (the default) then for a splitting to be considered, each subjet must satisfy.
p_{t,sub} > z_threshold * p_{t,parent}
whereas if it is set to true, then each subject must satisfy
p_{t,sub} > z_threshold * p_{t,original-jet}
where parent is the immediate parent of the splitting, and original jet is the one supplied to the run() function.
Only relevant is called prior to run().
Definition at line 146 of file CASubJetTagger.hh.
runs the tagger on the given jet and returns the tagged PseudoJet if successful, or a PseudoJet==0 otherwise (standard access is through operator()).
Implements fastjet::Transformer.
Definition at line 69 of file CASubJetTagger.cc.
|
inlineprotected |
work through the jet, establishing a distance at each branching
make sure the objects are not too close together
Definition at line 105 of file CASubJetTagger.cc.