FastJet 3.0.4
Classes | Public Types | Public Member Functions | Protected Member Functions | Protected Attributes | Static Protected Attributes
fastjet::CASubJetTagger Class Reference

clean (almost parameter-free) tagger searching for the element in the clustering history that maximises a chosen distance More...

#include <fastjet/tools/CASubJetTagger.hh>

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

List of all members.

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 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.
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()).

Protected Member Functions

void _recurse_through_jet (const PseudoJet &current_jet, JetAux &aux_max, const PseudoJet &original_jet) const

Protected Attributes

ScaleChoice _scale_choice
double _z_threshold
double _dr2_min
bool _absolute_z_cut

Static Protected Attributes

static LimitedWarning _non_ca_warnings

Detailed Description

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.

Options

Input conditions

Output/structure

Definition at line 104 of file CASubJetTagger.hh.


Member Function Documentation

void fastjet::CASubJetTagger::set_absolute_z_cut ( bool  abs_z_cut = true) [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 144 of file CASubJetTagger.hh.

PseudoJet fastjet::CASubJetTagger::result ( const PseudoJet jet) const [virtual]

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 67 of file CASubJetTagger.cc.

void fastjet::CASubJetTagger::_recurse_through_jet ( const PseudoJet jet,
JetAux aux,
const PseudoJet original_jet 
) const [inline, protected]

---------------------------------------------------------------------- work through the jet, establishing a distance at each branching

make sure the objects are not _too_ close together

Definition at line 103 of file CASubJetTagger.cc.


The documentation for this class was generated from the following files:
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends