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

base class for user-defined ordering of stable cones (used for prorgessive removal) More...

#include <fastjet/SISConeBasePlugin.hh>

Inheritance diagram for fastjet::SISConeBasePlugin::UserScaleBase:
Inheritance graph
[legend]
Collaboration diagram for fastjet::SISConeBasePlugin::UserScaleBase:
Collaboration graph
[legend]

Classes

class  StructureType
 the structure that allows to store the information contained into a siscone::Cjet (built internally in SISCone from a stable cone) into a PseudoJet More...
 

Public Member Functions

virtual ~UserScaleBase ()
 empty virtual dtor
 
virtual double result (const PseudoJet &jet) const =0
 returns the scale associated with a given jet More...
 
virtual bool is_larger (const PseudoJet &a, const PseudoJet &b) const
 returns true when the scale associated with jet a is larger than the scale associated with jet b More...
 
- Public Member Functions inherited from fastjet::FunctionOfPseudoJet< double >
 FunctionOfPseudoJet ()
 default ctor
 
virtual ~FunctionOfPseudoJet ()
 default dtor (virtual to allow safe polymorphism)
 
virtual std::string description () const
 returns a description of the function (an empty string by default)
 
double operator() (const PseudoJet &pj) const
 apply the function using the "traditional" () operator. More...
 
std::vector< double > operator() (const std::vector< PseudoJet > &pjs) const
 apply the function on a vector of PseudoJet, returning a vector of the results. More...
 

Detailed Description

base class for user-defined ordering of stable cones (used for prorgessive removal)

derived classes have to implement the () operator that returns the scale associated with a given jet.

It is also highly recommended to implement the is_larger() method whenever possible, in order to avoid rounding issues known to lead to possible infrared unsafeties.

The jets that are passed to this class will carry the structure of type SISConePlugin::StructureType which allows to retreive easily the following information:

vector<PseudoJet> constituents = jet.constituents(); unsigned int n_constituents = jet.structure_of<SISConePlugin::UserScaleBase>().size(); int index = jet.structure_of<SISConePlugin::UserScaleBase>().constituent_index(index i); const PseudoJet & p = jet.structure_of<SISConePlugin::UserScaleBase>().constituent(index i); double scalar_pt = jet.structure_of<SISConePlugin::UserScaleBase>().pt_tilde();

see SISConePlugin::StructureType below for further details

Definition at line 154 of file SISConeBasePlugin.hh.

Member Function Documentation

virtual double fastjet::SISConeBasePlugin::UserScaleBase::result ( const PseudoJet jet) const
pure virtual

returns the scale associated with a given jet

"progressive removal" iteratively removes the stable cone with the largest scale

Implements fastjet::FunctionOfPseudoJet< double >.

bool fastjet::SISConeBasePlugin::UserScaleBase::is_larger ( const PseudoJet a,
const PseudoJet b 
) const
virtual

returns true when the scale associated with jet a is larger than the scale associated with jet b

By default this does a simple direct comparison but it can be overloaded for higher precision [recommended if possible]

Definition at line 11 of file SISConeBasePlugin.cc.


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