FastJet  3.3.2
Namespaces | Classes | Typedefs | Enumerations | Functions | Variables
fastjet Namespace Reference

the FastJet namespace More...

Namespaces

 gas
 namespace to hold default parameters for the active area spec
 

Classes

class  _NoInfo
 internal dummy class, used as a default template argument More...
 
class  AreaDefinition
 class that holds a generic area definition More...
 
class  ATLASConePlugin
 Implementation of the ATLAS Cone (plugin for fastjet v2.4 upwards) More...
 
class  BackgroundEstimatorBase
 Abstract base class that provides the basic interface for classes that estimate levels of background radiation in hadron and heavy-ion collider events. More...
 
class  BackgroundJetPtDensity
 Class that implements pt/area_4vector.perp() for background estimation (this is a preliminary class). More...
 
class  BackgroundJetPtMDensity
 Class that implements $ \frac{1}{A} \sum_{i \in jet} (\sqrt{p_{ti}^2+m^2} - p_{ti}) $ for background estimation (this is a preliminary class). More...
 
class  BackgroundJetScalarPtDensity
 Class that implements (scalar pt sum of jet)/(scalar area of jet) for background estimation (this is a preliminary class). More...
 
class  BackgroundRescalingYPolynomial
 A background rescaling that is a simple polynomial in y. More...
 
class  Boost
 Class to boost a PseudoJet. More...
 
class  CASubJetTagger
 clean (almost parameter-free) tagger searching for the element in the clustering history that maximises a chosen distance More...
 
class  CASubJetTaggerStructure
 the structure returned by a CASubJetTagger More...
 
class  CDFJetCluPlugin
 Implementation of the JetClu algorithm from CDF (plugin for fastjet-v2.1 upwards) More...
 
class  CDFMidPointPlugin
 Implementation of the MidPoint algorithm from CDF (plugin for fastjet-v2.1 upwards) More...
 
class  ClusterSequence
 deals with clustering More...
 
class  ClusterSequence1GhostPassiveArea
 Like ClusterSequence with computation of the passive jet area by adding a single ghost. More...
 
class  ClusterSequenceActiveArea
 Like ClusterSequence with computation of the active jet area. More...
 
class  ClusterSequenceActiveAreaExplicitGhosts
 Like ClusterSequence with computation of the active jet area with the addition of explicit ghosts. More...
 
class  ClusterSequenceArea
 General class for user to obtain ClusterSequence with additional area information. More...
 
class  ClusterSequenceAreaBase
 base class that sets interface for extensions of ClusterSequence that provide information about the area of each jet More...
 
class  ClusterSequencePassiveArea
 Like ClusterSequence with computation of the passive jet area. More...
 
class  ClusterSequenceStructure
 Contains any information related to the clustering that should be directly accessible to PseudoJet. More...
 
class  ClusterSequenceVoronoiArea
 Like ClusterSequence with computation of the Voronoi jet area. More...
 
class  CMSIterativeConePlugin
 Implementation of the CMS Iterative Cone (plugin for fastjet v2.4 upwards) More...
 
class  CompositeJetStructure
 The structure for a jet made of pieces. More...
 
class  D0RunIBaseConePlugin
 D0RunIBaseConePlugin is base class for a plugin for FastJet (v3.0 or later) that provides an interface to the D0 version of Run-I cone algorithm. More...
 
class  D0RunIConePlugin
 A plugin for FastJet (v3.0 or later) that provides an interface to the D0 version of Run-I cone algorithm. More...
 
class  D0RunIIConePlugin
 Implementation of the D0 Run II Cone (plugin for fastjet v2.1 upwards) More...
 
class  D0RunIpre96ConePlugin
 A plugin for FastJet (v3.0 or later) that provides an interface to the pre 1996 D0 version of Run-I cone algorithm. More...
 
class  EECambridgePlugin
 Implementation of the e+e- Cambridge algorithm (plugin for fastjet v2.4 upwards) More...
 
class  Error
 base class corresponding to errors that can be thrown by FastJet More...
 
class  EtaPhi
 Shortcut for dealing with eta-phi coordinates. More...
 
class  Filter
 Class that helps perform filtering (Butterworth, Davison, Rubin and Salam, arXiv:0802.2470) and trimming (Krohn, Thaler and Wang, arXiv:0912.1342) on jets, optionally in conjunction with subtraction (Cacciari and Salam, arXiv:0707.1378). More...
 
class  FilterStructure
 Class to contain structure information for a filtered jet. More...
 
class  FunctionOfPseudoJet
 base class providing interface for a generic function of a PseudoJet More...
 
class  GhostedAreaSpec
 Parameters to configure the computation of jet areas using ghosts. More...
 
class  GridJetPlugin
 plugin for fastjet (v3.0 upwards) that clusters particles such that all particles in a given cell of a rectangular rapidity-phi grid end up in a common "jet". More...
 
class  GridMedianBackgroundEstimator
 Background Estimator based on the median pt/area of a set of grid cells. More...
 
class  InternalError
 class corresponding to critical internal errors More...
 
class  JadePlugin
 Implementation of the e+e- Jade algorithm (plugin for fastjet v2.4 upwards) More...
 
class  JetDefinition
 class that is intended to hold a full definition of the jet clusterer More...
 
class  JetMedianBackgroundEstimator
 Class to estimate the pt density of the background per unit area, using the median of the distribution of pt/area from jets that pass some selection criterion. More...
 
class  JHTopTagger
 Class that helps perform boosted top tagging using the "Johns Hopkins" method from arXiv:0806.0848 (Kaplan, Rehermann, Schwartz and Tweedie) More...
 
class  JHTopTaggerStructure
 the structure returned by the JHTopTagger transformer. More...
 
class  LimitedWarning
 class to provide facilities for giving warnings up to some maximum number of times and to provide global summaries of warnings that have been issued. More...
 
class  MassDropTagger
 Class that helps perform 2-pronged boosted tagging using the "mass-drop" technique (with asymmetry cut) introduced by Jonathan Butterworth, Adam Davison, Mathieu Rubin and Gavin Salam in arXiv:0802.2470 in the context of a boosted Higgs search. More...
 
class  MassDropTaggerStructure
 the structure returned by the MassDropTagger transformer. More...
 
class  NestedDefsPlugin
 Plugin to run multiple jet definitions successively (plugin for fastjet v2.4 upwards) More...
 
class  NNBase
 Helps solve closest pair problems with generic interparticle and particle-beam distances. More...
 
class  NNFJN2Plain
 Helps solve closest pair problems with factorised interparticle and beam distances (ie satisfying the FastJet lemma) More...
 
class  NNFJN2Tiled
 Helps solve closest pair problems with factorised interparticle and beam distances (ie satisfying the FastJet lemma) that are on a cylindrical geometry and allow tiling. More...
 
class  NNH
 Help solve closest pair problems with generic interparticle and beam distance (generic case) More...
 
class  NNInfo
 internal helper template class to facilitate initialisation of a BJ with a PseudoJet and extra information. More...
 
class  NNInfo< _NoInfo >
 for cases where there is no extra info More...
 
class  Pruner
 Transformer that prunes a jet. More...
 
class  PrunerStructure
 The structure associated with a PseudoJet thas has gone through a Pruner transformer. More...
 
class  PseudoJet
 Class to contain pseudojets, including minimal information of use to jet-clustering routines. More...
 
class  PseudoJetStructureBase
 Contains any information related to the clustering that should be directly accessible to PseudoJet. More...
 
class  PxConePlugin
 Implementation of the PxCone algorithm (plugin for fastjet v2.1 upwards) More...
 
class  RangeDefinition
 class for holding a range definition specification, given by limits on rapidity and azimuth. More...
 
class  Recluster
 Recluster a jet's constituents with a new jet definition. More...
 
class  RectangularGrid
 Class that holds a generic rectangular tiling. More...
 
class  RestFrameNSubjettinessTagger
 Class that helps perform 2-pronged boosted tagging using a reclustering in the jet's rest frame, supplemented with a cut on N-subjettiness (and a decay angle), as discussed by Ji-Hun Kim in arXiv:1011.1493. More...
 
class  RestFrameNSubjettinessTaggerStructure
 the structure returned by the RestFrameNSubjettinessTagger transformer. More...
 
class  Selector
 Class that encodes information about cuts and other selection criteria that can be applied to PseudoJet(s). More...
 
class  SelectorWorker
 default selector worker is an abstract virtual base class More...
 
class  SharedPtr
 an implementation of C++0x shared pointers (or boost's) More...
 
class  SISConeBaseExtras
 Class that provides extra information about a SISCone clustering. More...
 
class  SISConeExtras
 Class that provides extra information about a SISCone clustering. More...
 
class  SISConePlugin
 Implementation of the SISCone algorithm (plugin for fastjet v2.1 upwards) More...
 
class  SISConeSphericalExtras
 Class that provides extra information about a SISCone clustering. More...
 
class  SISConeSphericalPlugin
 Implementation of the spherical version of the SISCone algorithm (plugin for fastjet v2.1 upwards) More...
 
class  Subtractor
 Class that helps perform jet background subtraction. More...
 
class  TiledJet
 structure analogous to BriefJet, but with the extra information needed for dealing with tiles More...
 
class  TilingBase
 Class to indicate generic structure of tilings. More...
 
class  TilingExtent
 class to perform a fast analysis of the appropriate rapidity range in which to perform tiling More...
 
class  TopTaggerBase
 A base class that provides a common interface for top taggers that are able to return a W (in addition to the top itself). More...
 
class  TopTaggerBaseStructure
 class that specifies the structure common to all top taggers More...
 
class  TrackJetPlugin
 Implementation of the TrackJet algorithm (plugin for fastjet v2.4 upwards) More...
 
class  Transformer
 Base (abstract) class for a jet transformer. More...
 
class  Unboost
 Class to un-boost a PseudoJet. More...
 
class  VoronoiAreaSpec
 Specification for the computation of the Voronoi jet area. More...
 
class  WrappedStructure
 This wraps a (shared) pointer to an underlying structure. More...
 

Typedefs

typedef ClusterSequenceActiveAreaExplicitGhosts ClustSeqActAreaEG
 
typedef ClusterSequenceVoronoiArea::VoronoiAreaCalc VAC
 
typedef ClusterSequenceAreaBase ClusterSequenceWithArea
 
typedef GhostedAreaSpec ActiveAreaSpec
 just provide a typedef for backwards compatibility with programs based on versions 2.0 and 2.1 of fastjet. More...
 
typedef JetAlgorithm JetFinder
 make standard Les Houches nomenclature JetAlgorithm (algorithm is general recipe without the parameters) backward-compatible with old JetFinder
 
typedef integral_type< bool, true > true_type
 the bool 'true' value promoted to a type
 
typedef integral_type< bool, false > false_type
 the bool 'false' value promoted to a type
 
typedef char(& __yes_type)[1]
 
typedef char(& __no_type)[2]
 
typedef Tile2Base< 25 > Tile25
 
typedef Tile2Base< 9 > Tile2
 
typedef CGAL::Triangulation_vertex_base_with_info_2< InitialisedInt, K > Vbb
 
typedef CGAL::Triangulation_hierarchy_vertex_base_2< Vbb > Vb
 
typedef CGAL::Triangulation_face_base_2< K > Fb
 
typedef CGAL::Triangulation_data_structure_2< Vb, Fb > Tds
 
typedef CGAL::Delaunay_triangulation_2< K, Tds > Dt
 
typedef CGAL::Triangulation_hierarchy_2< Dt > Triangulation
 
typedef Triangulation::Vertex_handle Vertex_handle
 
typedef Triangulation::Point Point
 
typedef Triangulation::Vertex_circulator Vertex_circulator
 CGAL Point structure.
 
typedef Triangulation::Face_circulator Face_circulator
 
typedef Triangulation::Face_handle Face_handle
 

Enumerations

enum  AreaType {
  invalid_area = -1, active_area = 0, active_area_explicit_ghosts = 1, one_ghost_passive_area = 10,
  passive_area = 11, voronoi_area =20
}
 the different types of area that are supported
 
enum  Strategy {
  N2MHTLazy9AntiKtSeparateGhosts = -10, N2MHTLazy9 = -7, N2MHTLazy25 = -6, N2MHTLazy9Alt = -5,
  N2MinHeapTiled = -4, N2Tiled = -3, N2PoorTiled = -2, N2Plain = -1,
  N3Dumb = 0, Best = 1, NlnN = 2, NlnN3pi = 3,
  NlnN4pi = 4, NlnNCam4pi = 14, NlnNCam2pi2R = 13, NlnNCam = 12,
  BestFJ30 = 21, plugin_strategy = 999
}
 the various options for the algorithmic strategy to adopt in clustering events with kt and cambridge style algorithms. More...
 
enum  JetAlgorithm {
  kt_algorithm =0, cambridge_algorithm =1, antikt_algorithm =2, genkt_algorithm =3,
  cambridge_for_passive_algorithm =11, genkt_for_passive_algorithm =13, ee_kt_algorithm =50, ee_genkt_algorithm =53,
  plugin_algorithm = 99, undefined_jet_algorithm = 999
}
 the various families of jet-clustering algorithm More...
 
enum  RecombinationScheme {
  E_scheme =0, pt_scheme =1, pt2_scheme =2, Et_scheme =3,
  Et2_scheme =4, BIpt_scheme =5, BIpt2_scheme =6, WTA_pt_scheme =7,
  WTA_modp_scheme =8, external_scheme = 99
}
 The various recombination schemes. More...
 

Functions

int __default_random_generator (int *__iseed)
 
string fastjet_version_string ()
 return a string containing information about the release
 
PseudoJet join (const vector< PseudoJet > &pieces, const JetDefinition::Recombiner &recombiner)
 
PseudoJet join (const PseudoJet &j1, const JetDefinition::Recombiner &recombiner)
 build a "CompositeJet" from a single PseudoJet with an extended structure of type T derived from CompositeJetStructure More...
 
PseudoJet join (const PseudoJet &j1, const PseudoJet &j2, const JetDefinition::Recombiner &recombiner)
 build a "CompositeJet" from two PseudoJet with an extended structure of type T derived from CompositeJetStructure More...
 
PseudoJet join (const PseudoJet &j1, const PseudoJet &j2, const PseudoJet &j3, const JetDefinition::Recombiner &recombiner)
 build a "CompositeJet" from 3 PseudoJet with an extended structure of type T derived from CompositeJetStructure More...
 
PseudoJet join (const PseudoJet &j1, const PseudoJet &j2, const PseudoJet &j3, const PseudoJet &j4, const JetDefinition::Recombiner &recombiner)
 build a "CompositeJet" from 4 PseudoJet with an extended structure of type T derived from CompositeJetStructure More...
 
ostream & operator<< (ostream &ostr, const TiledJet &jet)
 
PseudoJet operator+ (const PseudoJet &jet1, const PseudoJet &jet2)
 
PseudoJet operator- (const PseudoJet &jet1, const PseudoJet &jet2)
 
PseudoJet operator* (double coeff, const PseudoJet &jet)
 
PseudoJet operator* (const PseudoJet &jet, double coeff)
 
PseudoJet operator/ (const PseudoJet &jet, double coeff)
 
bool operator== (const PseudoJet &, const PseudoJet &)
 returns true if the 4 momentum components of the two PseudoJets are identical and all the internal indices (user, cluster_history) More...
 
bool operator== (const PseudoJet &jet, const double val)
 Can only be used with val=0 and tests whether all four momentum components are equal to val (=0.0)
 
bool have_same_momentum (const PseudoJet &jeta, const PseudoJet &jetb)
 returns true if the momenta of the two input jets are identical
 
PseudoJet PtYPhiM (double pt, double y, double phi, double m)
 return a pseudojet with the given pt, y, phi and mass More...
 
void sort_indices (vector< int > &indices, const vector< double > &values)
 
vector< PseudoJetsorted_by_pt (const vector< PseudoJet > &jets)
 return a vector of jets sorted into decreasing kt2
 
vector< PseudoJetsorted_by_rapidity (const vector< PseudoJet > &jets)
 return a vector of jets sorted into increasing rapidity
 
vector< PseudoJetsorted_by_E (const vector< PseudoJet > &jets)
 return a vector of jets sorted into decreasing energy
 
vector< PseudoJetsorted_by_pz (const vector< PseudoJet > &jets)
 return a vector of jets sorted into increasing pz
 
PseudoJet join (const vector< PseudoJet > &pieces)
 
PseudoJet join (const PseudoJet &j1)
 build a "CompositeJet" from a single PseudoJet with an extended structure of type T derived from CompositeJetStructure More...
 
PseudoJet join (const PseudoJet &j1, const PseudoJet &j2)
 build a "CompositeJet" from two PseudoJet with an extended structure of type T derived from CompositeJetStructure More...
 
PseudoJet join (const PseudoJet &j1, const PseudoJet &j2, const PseudoJet &j3)
 build a "CompositeJet" from 3 PseudoJet with an extended structure of type T derived from CompositeJetStructure More...
 
PseudoJet join (const PseudoJet &j1, const PseudoJet &j2, const PseudoJet &j3, const PseudoJet &j4)
 build a "CompositeJet" from 4 PseudoJet with an extended structure of type T derived from CompositeJetStructure More...
 
Selector SelectorIdentity ()
 
Selector operator! (const Selector &s)
 logical not applied on a selector More...
 
Selector operator && (const Selector &s1, const Selector &s2)
 logical and between two selectors More...
 
Selector operator|| (const Selector &s1, const Selector &s2)
 logical or between two selectors More...
 
Selector operator* (const Selector &s1, const Selector &s2)
 successive application of 2 selectors More...
 
Selector SelectorPtMin (double ptmin)
 select objects with pt >= ptmin
 
Selector SelectorPtMax (double ptmax)
 select objects with pt <= ptmax
 
Selector SelectorPtRange (double ptmin, double ptmax)
 select objects with ptmin <= pt <= ptmax
 
Selector SelectorEtMin (double Etmin)
 select objects with Et >= Etmin
 
Selector SelectorEtMax (double Etmax)
 select objects with Et <= Etmax
 
Selector SelectorEtRange (double Etmin, double Etmax)
 select objects with Etmin <= Et <= Etmax
 
Selector SelectorEMin (double Emin)
 select objects with E >= Emin
 
Selector SelectorEMax (double Emax)
 select objects with E <= Emax
 
Selector SelectorERange (double Emin, double Emax)
 select objects with Emin <= E <= Emax
 
Selector SelectorMassMin (double Mmin)
 select objects with Mass >= Mmin
 
Selector SelectorMassMax (double Mmax)
 select objects with Mass <= Mmax
 
Selector SelectorMassRange (double Mmin, double Mmax)
 select objects with Mmin <= Mass <= Mmax
 
Selector SelectorRapMin (double rapmin)
 select objects with rap >= rapmin
 
Selector SelectorRapMax (double rapmax)
 select objects with rap <= rapmax
 
Selector SelectorRapRange (double rapmin, double rapmax)
 select objects with rapmin <= rap <= rapmax
 
Selector SelectorAbsRapMin (double absrapmin)
 select objects with |rap| >= absrapmin
 
Selector SelectorAbsRapMax (double absrapmax)
 select objects with |rap| <= absrapmax
 
Selector SelectorAbsRapRange (double absrapmin, double absrapmax)
 select objects with absrapmin <= |rap| <= absrapmax
 
Selector SelectorEtaMin (double etamin)
 select objects with eta >= etamin
 
Selector SelectorEtaMax (double etamax)
 select objects with eta <= etamax
 
Selector SelectorEtaRange (double etamin, double etamax)
 select objects with etamin <= eta <= etamax
 
Selector SelectorAbsEtaMin (double absetamin)
 select objects with |eta| >= absetamin
 
Selector SelectorAbsEtaMax (double absetamax)
 select objects with |eta| <= absetamax
 
Selector SelectorAbsEtaRange (double absetamin, double absetamax)
 select objects with absetamin <= |eta| <= absetamax
 
Selector SelectorPhiRange (double phimin, double phimax)
 select objects with phimin <= phi <= phimax
 
Selector SelectorRapPhiRange (double rapmin, double rapmax, double phimin, double phimax)
 select objects with rapmin <= rap <= rapmax && phimin <= phi <= phimax More...
 
Selector SelectorNHardest (unsigned int n)
 select the n hardest objects
 
Selector SelectorCircle (const double radius)
 select objets within a distance 'radius' from the location of the reference jet, set by Selector::set_reference(...)
 
Selector SelectorDoughnut (const double radius_in, const double radius_out)
 select objets with distance from the reference jet is between 'radius_in' and 'radius_out'; the reference jet is set by Selector::set_reference(...)
 
Selector SelectorStrip (const double half_width)
 select objets within a rapidity distance 'half_width' from the location of the reference jet, set by Selector::set_reference(...)
 
Selector SelectorRectangle (const double half_rap_width, const double half_phi_width)
 select objets within rapidity distance 'half_rap_width' from the reference jet and azimuthal-angle distance within 'half_phi_width'; the reference jet is set by Selector::set_reference(...)
 
Selector SelectorPtFractionMin (double fraction)
 select objects that carry at least a fraction "fraction" of the reference jet. More...
 
Selector SelectorIsZero ()
 select PseudoJet with 0 momentum
 
Selector SelectorIsPureGhost ()
 select objects that are (or are only made of) ghosts. More...
 
int scomp (const void *p1, const void *p2)
 
template<typename T >
PseudoJet join (const std::vector< PseudoJet > &pieces)
 build a "CompositeJet" from the vector of its pieces with an extended structure of type T derived from CompositeJetStructure More...
 
template<typename T >
PseudoJet join (const std::vector< PseudoJet > &pieces, const JetDefinition::Recombiner &recombiner)
 build a "CompositeJet" from the vector of its pieces with an extended structure of type T derived from CompositeJetStructure More...
 
bool operator!= (const PseudoJet &a, const PseudoJet &b)
 inequality test which is exact opposite of operator==
 
bool operator== (const double val, const PseudoJet &jet)
 
bool operator!= (const PseudoJet &a, const double val)
 Can only be used with val=0 and tests whether at least one of the four momentum components is different from val (=0.0)
 
bool operator!= (const double val, const PseudoJet &a)
 
double dot_product (const PseudoJet &a, const PseudoJet &b)
 returns the 4-vector dot product of a and b
 
double cos_theta (const PseudoJet &a, const PseudoJet &b)
 returns the cosine of the angle between a and b
 
double theta (const PseudoJet &a, const PseudoJet &b)
 returns the angle between a and b
 
std::vector< PseudoJetsorted_by_pt (const std::vector< PseudoJet > &jets)
 return a vector of jets sorted into decreasing transverse momentum
 
std::vector< PseudoJetsorted_by_rapidity (const std::vector< PseudoJet > &jets)
 return a vector of jets sorted into increasing rapidity
 
std::vector< PseudoJetsorted_by_E (const std::vector< PseudoJet > &jets)
 return a vector of jets sorted into decreasing energy
 
std::vector< PseudoJetsorted_by_pz (const std::vector< PseudoJet > &jets)
 return a vector of jets sorted into increasing pz
 
void sort_indices (std::vector< int > &indices, const std::vector< double > &values)
 sort the indices so that values[indices[0->n-1]] is sorted into increasing order
 
template<class T >
std::vector< T > objects_sorted_by_values (const std::vector< T > &objects, const std::vector< double > &values)
 given a vector of values with a one-to-one correspondence with the vector of objects, sort objects into an order such that the associated values would be in increasing order (but don't actually touch the values vector in the process). More...
 
Selector operator && (const Selector &s1, const Selector &s2)
 logical and between two selectors More...
 
template<class T , class U >
bool operator== (SharedPtr< T > const &t, SharedPtr< U > const &u)
 comparison: equality
 
template<class T , class U >
bool operator!= (SharedPtr< T > const &t, SharedPtr< U > const &u)
 comparison: difference
 
template<class T , class U >
bool operator< (SharedPtr< T > const &t, SharedPtr< U > const &u)
 comparison: orgering
 
template<class T >
void swap (SharedPtr< T > &a, SharedPtr< T > &b)
 swapping
 
template<class T >
T * get_pointer (SharedPtr< T > const &t)
 getting the pointer
 
bool floor_ln2_less (unsigned x, unsigned y)
 returns true if floor(ln_base2(x)) < floor(ln_base2(y)), using Chan's neat trick... More...
 
template<class B , class D >
B * cast_if_derived (D *d)
 a little helper that returns a pointer to d of type B* if D is derived from B and NULL otherwise
 
double norm (const VPoint p)
 norm of a vector
 
double vector_product (const VPoint &p1, const VPoint &p2)
 2D vector product
 
double scalar_product (const VPoint &p1, const VPoint &p2)
 scalar product
 
template<class T >
deltaPhi (T phi1, T phi2)
 
template<class T >
deltaR2 (T eta1, T phi1, T eta2, T phi2)
 

Variables

BasicRandom< int > _G_random_int
 
BasicRandom< double > _G_random_double
 
const unsigned int twopow31 = 2147483648U
 
const JetAlgorithm aachen_algorithm = cambridge_algorithm
 provide other possible names for the Cambridge/Aachen algorithm
 
const JetAlgorithm cambridge_aachen_algorithm = cambridge_algorithm
 
const double MaxRap = 1e5
 Used to protect against parton-level events where pt can be zero for some partons, giving rapidity=infinity. More...
 
const double pseudojet_invalid_phi = -100.0
 default value for phi, meaning it (and rapidity) have yet to be calculated)
 
const double pseudojet_invalid_rap = -1e200
 
const double tile_edge_security_margin =1.0e-7
 Rounding errors in the Lazy strategies may cause the following problem: when browsing tiles in the vicinity of the particles being clustered in order to decide which of these tiles may contain particles that need to be updated (because theit NN is one of the particles that are currently clustered), we discard tiles that are deemed "too far from the cell" by the "max_NN_dist" criterion. More...
 
const int n_tile_neighbours = 9
 
const double pi = 3.141592653589793238462643383279502884197
 
const double twopi = 6.283185307179586476925286766559005768394
 
const double pisq = 9.869604401089358618834490999876151135314
 
const double zeta2 = 1.644934066848226436472415166646025189219
 
const double zeta3 = 1.202056903159594285399738161511449990765
 
const double eulergamma = 0.577215664901532860606512090082402431042
 
const double ln2 = 0.693147180559945309417232121458176568076
 
const int INFINITE_VERTEX =-1
 
const int NEW_VERTEX =-2
 
const double HUGE_DOUBLE =1e300
 

Detailed Description

the FastJet namespace

all the fastjet-related material is put under that namespace

Typedef Documentation

◆ ActiveAreaSpec

just provide a typedef for backwards compatibility with programs based on versions 2.0 and 2.1 of fastjet.

Since there is no easy way of telling people this is deprecated at compile or run time, we should be careful before removing this in the future.

Definition at line 248 of file GhostedAreaSpec.hh.

Enumeration Type Documentation

◆ Strategy

the various options for the algorithmic strategy to adopt in clustering events with kt and cambridge style algorithms.

Enumerator
N2MHTLazy9AntiKtSeparateGhosts 

Like N2MHTLazy9 in a number of respects, but does not calculate ghost-ghost distances and so does not carry out ghost-ghost recombination.

If you want active ghosted areas, then this is only suitable for use with the anti-kt algorithm (or genkt with negative p), and does not produce any pure ghost jets. If used with active areas with Kt or Cam algorithms it will actually produce a passive area.

Particles are deemed to be ghosts if their pt is below a threshold (currently 1e-50, hard coded as ghost_limit in LazyTiling9SeparateGhosts).

Currently for events with a couple of thousand normal particles and O(10k) ghosts, this can be quicker than N2MHTLazy9, which would otherwise be the best strategy.

New in FJ3.1

N2MHTLazy9 

only looks into a neighbouring tile for a particle's nearest neighbour (NN) if that particle's in-tile NN is further than the distance to the edge of the neighbouring tile.

Uses tiles of size R and a 3x3 tile grid around the particle. New in FJ3.1

N2MHTLazy25 

Similar to N2MHTLazy9, but uses tiles of size R/2 and a 5x5 tile grid around the particle.

New in FJ3.1

N2MHTLazy9Alt 

Like to N2MHTLazy9 but uses slightly different optimizations, e.g.

for calculations of distance to nearest tile; as of 2014-07-18 it is slightly slower and not recommended for production use. To considered deprecated. New in FJ3.1

N2MinHeapTiled 

faster that N2Tiled above about 500 particles; differs from it by retainig the di(closest j) distances in a MinHeap (sort of priority queue) rather than a simple vector.

N2Tiled 

fastest from about 50..500

N2PoorTiled 

legacy

N2Plain 

fastest below 50

N3Dumb 

worse even than the usual N^3 algorithms

Best 

automatic selection of the best (based on N), including the LazyTiled strategies that are new to FJ3.1

NlnN 

best of the NlnN variants – best overall for N>10^4.

(Does not work for R>=2pi)

NlnN3pi 

legacy N ln N using 3pi coverage of cylinder.

(Does not work for R>=2pi)

NlnN4pi 

legacy N ln N using 4pi coverage of cylinder

NlnNCam4pi 

Chan's closest pair method (in a variant with 4pi coverage), for use exclusively with the Cambridge algorithm.

(Does not work for R>=2pi)

NlnNCam2pi2R 

Chan's closest pair method (in a variant with 2pi+2R coverage), for use exclusively with the Cambridge algorithm.

(Does not work for R>=2pi)

NlnNCam 

Chan's closest pair method (in a variant with 2pi+minimal extra variant), for use exclusively with the Cambridge algorithm.

(Does not work for R>=2pi)

BestFJ30 

the automatic strategy choice that was being made in FJ 3.0 (restricted to strategies that were present in FJ 3.0)

plugin_strategy 

the plugin has been used...

Definition at line 51 of file JetDefinition.hh.

◆ JetAlgorithm

the various families of jet-clustering algorithm

Enumerator
kt_algorithm 

the longitudinally invariant kt algorithm

cambridge_algorithm 

the longitudinally invariant variant of the cambridge algorithm (aka Aachen algoithm).

antikt_algorithm 

like the k_t but with distance measures dij = min(1/kti^2,1/ktj^2) Delta R_{ij}^2 / R^2 diB = 1/kti^2

genkt_algorithm 

like the k_t but with distance measures dij = min(kti^{2p},ktj^{2p}) Delta R_{ij}^2 / R^2 diB = 1/kti^{2p} where p = extra_param()

cambridge_for_passive_algorithm 

a version of cambridge with a special distance measure for particles whose pt is < extra_param(); this is not usually intended for end users, but is instead automatically selected when requesting a passive Cambridge area.

genkt_for_passive_algorithm 

a version of genkt with a special distance measure for particles whose pt is < extra_param() [relevant for passive areas when p<=0] ***** NB: THERE IS CURRENTLY NO IMPLEMENTATION FOR THIS ALG *******

ee_kt_algorithm 

the e+e- kt algorithm

ee_genkt_algorithm 

the e+e- genkt algorithm (R > 2 and p=1 gives ee_kt)

plugin_algorithm 

any plugin algorithm supplied by the user

undefined_jet_algorithm 

the value for the jet algorithm in a JetDefinition for which no algorithm has yet been defined

Definition at line 137 of file JetDefinition.hh.

◆ RecombinationScheme

The various recombination schemes.

Note that the schemes that recombine with non-linear weighting of the directions (e.g. pt2, winner-takes-all) are collinear safe only for algorithms with a suitable ordering of the recombinations: orderings in which, for particles of comparable energies, small-angle clusterings take place before large-angle clusterings. This property is satisfied by all gen-kt algorithms.

Enumerator
E_scheme 

summing the 4-momenta

pt_scheme 

pt weighted recombination of y,phi (and summing of pt's) with preprocessing to make things massless by rescaling E=| p|

pt2_scheme 

pt^2 weighted recombination of y,phi (and summing of pt's) with preprocessing to make things massless by rescaling E=| p|

Et_scheme 

pt weighted recombination of y,phi (and summing of pt's) with preprocessing to make things massless by rescaling | p|->=E

Et2_scheme 

pt^2 weighted recombination of y,phi (and summing of pt's) with preprocessing to make things massless by rescaling | p|->=E

BIpt_scheme 

pt weighted recombination of y,phi (and summing of pt's), with no preprocessing

BIpt2_scheme 

pt^2 weighted recombination of y,phi (and summing of pt's) no preprocessing

WTA_pt_scheme 

pt-based Winner-Takes-All (WTA) recombination: the result of the recombination has the rapidity, azimuth and mass of the PseudoJet with the larger pt, and a pt equal to the sum of the two pt's

WTA_modp_scheme 

mod-p-based Winner-Takes-All (WTA) recombination: the result of the recombination gets the 3-vector direction and mass of the PseudoJet with the larger |3-momentum| (modp), and a |3-momentum| equal to the scalar sum of the two |3-momenta|.

external_scheme 

for the user's external scheme

Definition at line 193 of file JetDefinition.hh.

Function Documentation

◆ join() [1/10]

PseudoJet fastjet::join ( const PseudoJet j1,
const JetDefinition::Recombiner recombiner 
)

build a "CompositeJet" from a single PseudoJet with an extended structure of type T derived from CompositeJetStructure

build a MergedJet from a single PseudoJet

Definition at line 475 of file JetDefinition.cc.

◆ join() [2/10]

PseudoJet fastjet::join ( const PseudoJet j1,
const PseudoJet j2,
const JetDefinition::Recombiner recombiner 
)

build a "CompositeJet" from two PseudoJet with an extended structure of type T derived from CompositeJetStructure

build a MergedJet from 2 PseudoJet

Definition at line 481 of file JetDefinition.cc.

◆ join() [3/10]

PseudoJet fastjet::join ( const PseudoJet j1,
const PseudoJet j2,
const PseudoJet j3,
const JetDefinition::Recombiner recombiner 
)

build a "CompositeJet" from 3 PseudoJet with an extended structure of type T derived from CompositeJetStructure

build a MergedJet from 3 PseudoJet

Definition at line 490 of file JetDefinition.cc.

◆ join() [4/10]

PseudoJet fastjet::join ( const PseudoJet j1,
const PseudoJet j2,
const PseudoJet j3,
const PseudoJet j4,
const JetDefinition::Recombiner recombiner 
)

build a "CompositeJet" from 4 PseudoJet with an extended structure of type T derived from CompositeJetStructure

build a MergedJet from 4 PseudoJet

Definition at line 500 of file JetDefinition.cc.

◆ operator==()

bool fastjet::operator== ( const PseudoJet ,
const PseudoJet  
)

returns true if the 4 momentum components of the two PseudoJets are identical and all the internal indices (user, cluster_history)

  • structure and user-info shared pointers are too

Definition at line 255 of file PseudoJet.cc.

◆ PtYPhiM()

PseudoJet fastjet::PtYPhiM ( double  pt,
double  y,
double  phi,
double  m 
)

return a pseudojet with the given pt, y, phi and mass

return a pseudojet with the given pt, y, phi and mass (phi should satisfy -2pi<phi<4pi)

Definition at line 365 of file PseudoJet.cc.

◆ join() [5/10]

PseudoJet fastjet::join ( const PseudoJet j1)

build a "CompositeJet" from a single PseudoJet with an extended structure of type T derived from CompositeJetStructure

build a MergedJet from a single PseudoJet

Definition at line 831 of file PseudoJet.cc.

◆ join() [6/10]

PseudoJet fastjet::join ( const PseudoJet j1,
const PseudoJet j2 
)

build a "CompositeJet" from two PseudoJet with an extended structure of type T derived from CompositeJetStructure

build a MergedJet from 2 PseudoJet

Definition at line 836 of file PseudoJet.cc.

◆ join() [7/10]

PseudoJet fastjet::join ( const PseudoJet j1,
const PseudoJet j2,
const PseudoJet j3 
)

build a "CompositeJet" from 3 PseudoJet with an extended structure of type T derived from CompositeJetStructure

build a MergedJet from 3 PseudoJet

Definition at line 845 of file PseudoJet.cc.

◆ join() [8/10]

PseudoJet fastjet::join ( const PseudoJet j1,
const PseudoJet j2,
const PseudoJet j3,
const PseudoJet j4 
)

build a "CompositeJet" from 4 PseudoJet with an extended structure of type T derived from CompositeJetStructure

build a MergedJet from 4 PseudoJet

Definition at line 855 of file PseudoJet.cc.

◆ join() [9/10]

template<typename T >
PseudoJet fastjet::join ( const std::vector< PseudoJet > &  pieces)

build a "CompositeJet" from the vector of its pieces with an extended structure of type T derived from CompositeJetStructure

build a "CompositeJet" from the vector of its pieces

In this case, E-scheme recombination is assumed to compute the total momentum

Definition at line 155 of file CompositeJetStructure.hh.

◆ join() [10/10]

template<typename T >
PseudoJet fastjet::join ( const std::vector< PseudoJet > &  pieces,
const JetDefinition::Recombiner recombiner 
)

build a "CompositeJet" from the vector of its pieces with an extended structure of type T derived from CompositeJetStructure

build a "CompositeJet" from the vector of its pieces

In this case, E-scheme recombination is assumed to compute the total momentum

Definition at line 218 of file CompositeJetStructure.hh.

◆ objects_sorted_by_values()

template<class T >
std::vector<T> fastjet::objects_sorted_by_values ( const std::vector< T > &  objects,
const std::vector< double > &  values 
)

given a vector of values with a one-to-one correspondence with the vector of objects, sort objects into an order such that the associated values would be in increasing order (but don't actually touch the values vector in the process).

Definition at line 907 of file PseudoJet.hh.

◆ floor_ln2_less()

bool fastjet::floor_ln2_less ( unsigned  x,
unsigned  y 
)
inline

returns true if floor(ln_base2(x)) < floor(ln_base2(y)), using Chan's neat trick...

Definition at line 221 of file ClosestPair2D.hh.

Variable Documentation

◆ MaxRap

const double fastjet::MaxRap = 1e5

Used to protect against parton-level events where pt can be zero for some partons, giving rapidity=infinity.

KtJet fails in those cases.

Definition at line 52 of file PseudoJet.hh.

◆ tile_edge_security_margin

const double fastjet::tile_edge_security_margin =1.0e-7

Rounding errors in the Lazy strategies may cause the following problem: when browsing tiles in the vicinity of the particles being clustered in order to decide which of these tiles may contain particles that need to be updated (because theit NN is one of the particles that are currently clustered), we discard tiles that are deemed "too far from the cell" by the "max_NN_dist" criterion.

Because of rounding error, this condition can sometimes miss cases where an update is needed.

An example of this happens if a particle '1' is, say, at the lower edge of the rapidity of a given tile, with a particle '2' in the tile directly on its left at the same rapidity. Assume also that max_NN_dist in 2's tile corresponds to the distance between 2 and teh tile of 1. If 2 is 1's NN then in case 2 gets clustered, 1's NN needs to be updated. However, rounding errors in the calculation of the distance between 1 and 2 may result is something slightly larger than the max_NN_dist in 2's tile.

This situation corresponds to the bug reported by Jochen Olt on February 12 2015 [see issue-tracker/2015-02-infinite-loop], causing an infinite loop.

To prevent this, the simplest solution is, when looking at tiles to browse for updateds, to add a margin of security close to the edges of the cell, i.e. instead of updating only tiles for which distance<=max_NN_dist, we will update tiles for which distance<=max_NN_dist+tile_edge_security_margin.

Note that this does not need to be done when computing nearest neighbours [rounding errors are tolerated there] but it is critical when tracking points that have to be updated.

Definition at line 71 of file LazyTiling9Alt.hh.