Main Page | Namespace List | Class Hierarchy | Class List | Directories | File List | Namespace Members | Class Members | File Members

fastjet Namespace Reference


Classes

class  ActiveAreaSpec
 Class that defines the parameters that go into the measurement of active jet areas. More...
class  ClusterSequence
 deals with clustering More...
class  ClusterSequenceActiveArea
 Class that behaves essentially like ClusterSequence except that it also provides access to the area of a jet (which will be a random quantity. More...
class  ClusterSequenceActiveAreaExplicitGhosts
 Class that behaves essentially like ClusterSequence except that it also provides access to the area of a jet (which will be a random quantity. More...
class  ClusterSequenceWithArea
 base class that sets interface for extensions of ClusterSequence that provide information about the area of each jet; More...
class  Error
 class corresponding to errors that will be thrown by fastjet More...
class  JetDefinition
 class that is intended to hold a full definition of the jet clusterer More...
class  PseudoJet
 Class to contain pseudojets, including minimal information of use to to jet-clustering routines. More...
class  IndexedSortHelper
 a class that helps us carry out indexed sorting. More...
class  BasicRandom
class  BasicRandom< int >
class  BasicRandom< double >
 template specialization (double) for the BasicRandom template class. More...
class  ClosestPair2D
 concrete implementation for finding closest pairs in 2D -- will use Chan's (hopefully efficient) shuffle based structures More...
class  Coord2D
 class for representing 2d coordinates and carrying out some basic operations on them More...
class  ClosestPair2DBase
 abstract base class for finding closest pairs in 2D More...
class  Dnn2piCylinder
 class derived from DynamicNearestNeighbours that provides an implementation for the surface of cylinder (using one DnnPlane object spanning 0--2pi). More...
class  Dnn3piCylinder
 class derived from DynamicNearestNeighbours that provides an implementation for the surface of cylinder (using one DnnPlane object spanning 0--3pi). More...
class  Dnn4piCylinder
 class derived from DynamicNearestNeighbours that provides an implementation for the surface of cylinder (using two copies of DnnPlane, one running from 0--2pi, the other from pi--3pi). More...
class  DnnPlane
 class derived from DynamicNearestNeighbours that provides an implementation for the Euclidean plane More...
class  EtaPhi
 use a class instead of a pair so that phi can be sanitized and put into proper range on initialization. More...
class  DnnError
 class corresponding to errors that will be thrown by Dynamic Nearest Neighbours code More...
class  DynamicNearestNeighbours
 Abstract base class for quick location of nearest neighbours in a set of points, with facilities for adding and removing points from the set after initialisation. More...
class  MinHeap
 A class which provides a "heap"-like structure that allows access to a the minimal value of a dynamically changing set of numbers. More...
class  SearchTree
 This is the class for a search tree designed to be especially efficient when looking for successors and predecessors (to be used in Chan's CP algorithm). More...
struct  K
 the basic geometrical kernel that lies at the base of all CGAL operations More...
class  InitialisedInt
 A class to provide an "int" with an initial value. More...
class  SISConePlugin
 SISConePlugin is a plugin for fastjet (v2.1 upwards) that provides an interface to the seedless infrared safe cone jet finder by Gregory Soyez and Gavin Salam. More...
class  SISConeExtras
 Class that provides extra information about a SISCone clustering. More...
class  PxConePlugin
 PxConePlugin is a plugin for fastjet (v2.1 upwards) that provides an interface to the fortran pxcone iterative cone algorithm with midpoint seeds. More...
class  CDFJetCluPlugin
 a plugin for fastjet-v2.1 that provides an interface to the CDF jetclu algorithm More...
class  CDFMidPointPlugin
 CDFMidPointPlugin is a plugin for fastjet (v2.1 upwards) that provides an interface to the CDF version of Run-II iterative cone algorithm with midpoint seeds (also known as the Iterative Legacy Cone Algorithm, ILCA). More...

Namespaces

namespace  Private

Typedefs

typedef ClusterSequenceActiveAreaExplicitGhosts ClustSeqActAreaEG
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  Strategy {
  N2MinHeapTiled = -4, N2Tiled = -3, N2PoorTiled = -2, N2Plain = -1,
  N3Dumb = 0, Best = 1, NlnN = 2, NlnN3pi = 3,
  NlnN4pi = 4, NlnNCam4pi = 14, NlnNCam2pi2R = 13, NlnNCam = 12,
  plugin_strategy = 999
}
 the various options for the algorithmic strategy to adopt in clustering events with kt and cambridge style algorithms. More...
enum  JetFinder { kt_algorithm = 0, cambridge_algorithm = 1, plugin_algorithm = 99 }
 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, external_scheme = 99
}
 the various recombination schemes More...

Functions

int __default_random_generator (int *__iseed)
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 have_same_momentum (const PseudoJet &jeta, const PseudoJet &jetb)
 returns true if the momenta of the two input jets are identical
void sort_indices (vector< int > &indices, const vector< double > &values)
template<class T>
vector< T > objects_sorted_by_values (const vector< T > &objects, const 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
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
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
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).
bool floor_ln2_less (unsigned x, unsigned y)
 returns true if floor(ln_base2(x)) < floor(ln_base2(y)), using Chan's neat trick.

Variables

BasicRandom< int > _G_random_int
BasicRandom< double > _G_random_double
const unsigned int huge_unsigned = 4294967295U
const unsigned int twopow31 = 2147483648U
const double MaxRap = 1e5
 Used to protect against parton-level events where pt can be zero for some partons, giving rapidity=infinity.
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


Typedef Documentation

typedef ClusterSequenceActiveAreaExplicitGhosts fastjet::ClustSeqActAreaEG
 

Definition at line 38 of file ClusterSequenceActiveAreaExplicitGhosts.cc.

typedef CGAL::Delaunay_triangulation_2<K,Tds> fastjet::Dt
 

Definition at line 85 of file Triangulation.hh.

typedef Triangulation::Face_circulator fastjet::Face_circulator
 

Definition at line 92 of file Triangulation.hh.

typedef Triangulation::Face_handle fastjet::Face_handle
 

Definition at line 93 of file Triangulation.hh.

typedef CGAL::Triangulation_face_base_2<K> fastjet::Fb
 

Definition at line 83 of file Triangulation.hh.

typedef Triangulation::Point fastjet::Point
 

Definition at line 90 of file Triangulation.hh.

typedef CGAL::Triangulation_data_structure_2<Vb,Fb> fastjet::Tds
 

Definition at line 84 of file Triangulation.hh.

typedef CGAL::Triangulation_hierarchy_2<Dt> fastjet::Triangulation
 

Definition at line 86 of file Triangulation.hh.

typedef CGAL::Triangulation_hierarchy_vertex_base_2<Vbb> fastjet::Vb
 

Definition at line 82 of file Triangulation.hh.

typedef CGAL::Triangulation_vertex_base_with_info_2<InitialisedInt,K> fastjet::Vbb
 

Definition at line 81 of file Triangulation.hh.

typedef Triangulation::Vertex_circulator fastjet::Vertex_circulator
 

CGAL Point structure.

Definition at line 91 of file Triangulation.hh.

typedef Triangulation::Vertex_handle fastjet::Vertex_handle
 

Definition at line 89 of file Triangulation.hh.


Enumeration Type Documentation

enum fastjet::JetFinder
 

the various families of jet-clustering algorithm

Enumeration values:
kt_algorithm  the longitudinally invariant kt algorithm
cambridge_algorithm  the longitudinally invariant variant of the cambridge algorithm (aka Aachen algoithm).
plugin_algorithm  any plugin algorithm supplied by the user

Definition at line 76 of file JetDefinition.hh.

00076                {
00078   kt_algorithm=0,
00081   cambridge_algorithm=1,
00083   plugin_algorithm = 99
00084 };

enum fastjet::RecombinationScheme
 

the various recombination schemes

Enumeration values:
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
external_scheme  for the user's external scheme

Definition at line 89 of file JetDefinition.hh.

00089                          {
00091   E_scheme=0,
00094   pt_scheme=1,
00097   pt2_scheme=2,
00100   Et_scheme=3,
00103   Et2_scheme=4,
00106   BIpt_scheme=5,
00109   BIpt2_scheme=6,
00111   external_scheme = 99
00112 };

enum fastjet::Strategy
 

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

Enumeration values:
N2MinHeapTiled  experimental ...
N2Tiled  fastest from about 50..10^4
N2PoorTiled  legacy
N2Plain  fastest below 50
N3Dumb  worse even than the usual N^3 algorithms
Best  automatic selection of the best (based on N)
NlnN  best of the NlnN variants -- best overall for N>10^4
NlnN3pi  legacy N ln N using 3pi coverage of cylinder
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.
NlnNCam2pi2R 
NlnNCam 
plugin_strategy  the plugin has been used...

Definition at line 45 of file JetDefinition.hh.

00045               {
00047   N2MinHeapTiled   = -4, 
00049   N2Tiled     = -3, 
00051   N2PoorTiled = -2, 
00053   N2Plain     = -1, 
00055   N3Dumb      =  0, 
00057   Best        =  1, 
00059   NlnN        =  2, 
00061   NlnN3pi     =  3, 
00063   NlnN4pi     =  4,
00066   NlnNCam4pi   = 14,
00067   NlnNCam2pi2R = 13,
00068   NlnNCam      = 12, // 2piMultD
00070   plugin_strategy = 999
00071 };


Function Documentation

int fastjet::__default_random_generator int *  __iseed  ) 
 

Definition at line 31 of file BasicRandom.cc.

Referenced by fastjet::BasicRandom< double >::operator()(), and fastjet::BasicRandom< int >::operator()().

00032 {
00033   int __k = __iseed[0]/53668;
00034   __iseed[0] = (__iseed[0] - __k*53668)*40014 - __k*12211;
00035   if(__iseed[0] < 0) __iseed[0] += 2147483563;
00036   
00037   __k = __iseed[1]/52774;
00038   __iseed[1] = (__iseed[1] - __k*52774)*40692 - __k*3791;
00039   if(__iseed[1] < 0) __iseed[1] += 2147483399;
00040   
00041   int __iz = __iseed[0] - __iseed[1];
00042   if(__iz < 1) __iz += 2147483562;
00043   
00044   return __iz;
00045 }

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 212 of file ClosestPair2D.hh.

Referenced by fastjet::ClosestPair2D::Shuffle::operator<().

00212                                                    {
00213   if (x>y) return false;
00214   return (x < (x^y)); // beware of operator precedence...
00215 }

bool fastjet::have_same_momentum const PseudoJet &  jeta,
const PseudoJet &  jetb
 

returns true if the momenta of the two input jets are identical

Definition at line 272 of file PseudoJet.cc.

References fastjet::PseudoJet::E(), fastjet::PseudoJet::px(), fastjet::PseudoJet::py(), and fastjet::PseudoJet::pz().

Referenced by fastjet::SISConePlugin::run_clustering().

00272                                                                         {
00273   return jeta.px() == jetb.px()
00274     &&   jeta.py() == jetb.py()
00275     &&   jeta.pz() == jetb.pz()
00276     &&   jeta.E()  == jetb.E();
00277 }

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

template<class T>
vector<T> fastjet::objects_sorted_by_values const vector< T > &  objects,
const 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

Definition at line 315 of file PseudoJet.cc.

References sort_indices().

Referenced by sorted_by_E(), sorted_by_pt(), and sorted_by_rapidity().

00317                                                       {
00318 
00319   assert(objects.size() == values.size());
00320 
00321   // get a vector of indices
00322   vector<int> indices(values.size());
00323   for (size_t i = 0; i < indices.size(); i++) {indices[i] = i;}
00324   
00325   // sort the indices
00326   sort_indices(indices, values);
00327   
00328   // copy the objects 
00329   vector<T> objects_sorted(objects.size());
00330   
00331   // place the objects in the correct order
00332   for (size_t i = 0; i < indices.size(); i++) {
00333     objects_sorted[i] = objects[indices[i]];
00334   }
00335 
00336   return objects_sorted;
00337 }

PseudoJet fastjet::operator * const PseudoJet &  jet,
double  coeff
 

Definition at line 167 of file PseudoJet.cc.

00167                                                           {
00168   return coeff*jet;
00169 } 

PseudoJet fastjet::operator * double  coeff,
const PseudoJet &  jet
 

Definition at line 157 of file PseudoJet.cc.

00157                                                           {
00158   //return PseudoJet(coeff*jet.four_mom());
00159   // the following code is hopefully more efficient
00160   PseudoJet coeff_times_jet(jet);
00161   coeff_times_jet *= coeff;
00162   return coeff_times_jet;
00163 } 

PseudoJet fastjet::operator+ const PseudoJet &  jet1,
const PseudoJet &  jet2
 

Definition at line 137 of file PseudoJet.cc.

References fastjet::PseudoJet::E(), fastjet::PseudoJet::px(), fastjet::PseudoJet::py(), and fastjet::PseudoJet::pz().

00137                                                                      {
00138   //return PseudoJet(jet1.four_mom()+jet2.four_mom());
00139   return PseudoJet(jet1.px()+jet2.px(),
00140                    jet1.py()+jet2.py(),
00141                    jet1.pz()+jet2.pz(),
00142                    jet1.E() +jet2.E()  );
00143 } 

PseudoJet fastjet::operator- const PseudoJet &  jet1,
const PseudoJet &  jet2
 

Definition at line 147 of file PseudoJet.cc.

References fastjet::PseudoJet::E(), fastjet::PseudoJet::px(), fastjet::PseudoJet::py(), and fastjet::PseudoJet::pz().

00147                                                                      {
00148   //return PseudoJet(jet1.four_mom()-jet2.four_mom());
00149   return PseudoJet(jet1.px()-jet2.px(),
00150                    jet1.py()-jet2.py(),
00151                    jet1.pz()-jet2.pz(),
00152                    jet1.E() -jet2.E()  );
00153 } 

PseudoJet fastjet::operator/ const PseudoJet &  jet,
double  coeff
 

Definition at line 173 of file PseudoJet.cc.

00173                                                           {
00174   return (1.0/coeff)*jet;
00175 } 

void fastjet::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

void fastjet::sort_indices vector< int > &  indices,
const vector< double > &  values
 

Definition at line 305 of file PseudoJet.cc.

Referenced by objects_sorted_by_values().

00306                                                         {
00307   IndexedSortHelper index_sort_helper(&values);
00308   sort(indices.begin(), indices.end(), index_sort_helper);
00309 }

std::vector<PseudoJet> fastjet::sorted_by_E const std::vector< PseudoJet > &  jets  ) 
 

return a vector of jets sorted into decreasing energy

vector<PseudoJet> fastjet::sorted_by_E const vector< PseudoJet > &  jets  ) 
 

return a vector of jets sorted into decreasing energy

Definition at line 357 of file PseudoJet.cc.

References objects_sorted_by_values().

Referenced by main().

00357                                                               {
00358   vector<double> energies(jets.size());
00359   for (size_t i = 0; i < jets.size(); i++) {energies[i] = -jets[i].E();}
00360   return objects_sorted_by_values(jets, energies);
00361 }

std::vector<PseudoJet> fastjet::sorted_by_pt const std::vector< PseudoJet > &  jets  ) 
 

return a vector of jets sorted into decreasing transverse momentum

vector<PseudoJet> fastjet::sorted_by_pt const vector< PseudoJet > &  jets  ) 
 

return a vector of jets sorted into decreasing kt2

Definition at line 341 of file PseudoJet.cc.

References objects_sorted_by_values().

Referenced by main(), and print_jets().

00341                                                                {
00342   vector<double> minus_kt2(jets.size());
00343   for (size_t i = 0; i < jets.size(); i++) {minus_kt2[i] = -jets[i].kt2();}
00344   return objects_sorted_by_values(jets, minus_kt2);
00345 }

std::vector<PseudoJet> fastjet::sorted_by_rapidity const std::vector< PseudoJet > &  jets  ) 
 

return a vector of jets sorted into increasing rapidity

vector<PseudoJet> fastjet::sorted_by_rapidity const vector< PseudoJet > &  jets  ) 
 

return a vector of jets sorted into increasing rapidity

Definition at line 349 of file PseudoJet.cc.

References objects_sorted_by_values().

00349                                                                      {
00350   vector<double> rapidities(jets.size());
00351   for (size_t i = 0; i < jets.size(); i++) {rapidities[i] = jets[i].rap();}
00352   return objects_sorted_by_values(jets, rapidities);
00353 }


Variable Documentation

BasicRandom< double > fastjet::_G_random_double
 

Definition at line 49 of file BasicRandom.cc.

BasicRandom< int > fastjet::_G_random_int
 

Definition at line 48 of file BasicRandom.cc.

const double fastjet::eulergamma = 0.577215664901532860606512090082402431042
 

Definition at line 47 of file numconsts.hh.

const double fastjet::HUGE_DOUBLE = 1e300
 

Definition at line 56 of file Triangulation.hh.

Referenced by fastjet::DnnPlane::_SetAndUpdateNearest(), and fastjet::DnnPlane::_SetNearest().

const unsigned int fastjet::huge_unsigned = 4294967295U
 

Definition at line 39 of file ClosestPair2D.cc.

const int fastjet::INFINITE_VERTEX = -1
 

Definition at line 54 of file Triangulation.hh.

Referenced by fastjet::DnnPlane::_SetAndUpdateNearest(), fastjet::DnnPlane::_SetNearest(), fastjet::DnnPlane::DnnPlane(), and fastjet::DnnPlane::RemoveAndAddPoints().

const double fastjet::ln2 = 0.693147180559945309417232121458176568076
 

Definition at line 48 of file numconsts.hh.

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 50 of file PseudoJet.hh.

Referenced by fastjet::PseudoJet::_finish_init().

const int fastjet::NEW_VERTEX = -2
 

Definition at line 55 of file Triangulation.hh.

Referenced by fastjet::DnnPlane::_CrashIfVertexPresent(), and fastjet::InitialisedInt::InitialisedInt().

const double fastjet::pi = 3.141592653589793238462643383279502884197
 

Definition at line 42 of file numconsts.hh.

Referenced by fastjet::ClusterSequence::_bj_dist(), fastjet::ClusterSequence::_CP2DChan_limited_cluster(), fastjet::Dnn3piCylinder::_RegisterCylinderPoint(), fastjet::Dnn2piCylinder::_RegisterCylinderPoint(), fastjet::Dnn4piCylinder::_remap_phi(), fastjet::Dnn3piCylinder::_remap_phi(), fastjet::Dnn2piCylinder::_remap_phi(), fastjet::Dnn4piCylinder::Dnn4piCylinder(), fastjet::JetDefinition::JetDefinition(), fastjet::PseudoJet::kt_distance(), fastjet::PseudoJet::phi_std(), fastjet::PseudoJet::plain_distance(), and fastjet::JetDefinition::DefaultRecombiner::recombine().

const double fastjet::pisq = 9.869604401089358618834490999876151135314
 

Definition at line 44 of file numconsts.hh.

const double fastjet::twopi = 6.283185307179586476925286766559005768394
 

Definition at line 43 of file numconsts.hh.

Referenced by fastjet::ClusterSequence::_bj_dist(), fastjet::ClusterSequence::_CP2DChan_cluster(), fastjet::Dnn2piCylinder::_CreateNecessaryMirrorPoints(), fastjet::ClusterSequence::_delaunay_cluster(), fastjet::PseudoJet::_finish_init(), fastjet::ClusterSequence::_initialise_tiles(), fastjet::ActiveAreaSpec::_initialize(), fastjet::ClusterSequence::_tile_index(), fastjet::PseudoJet::kt_distance(), main(), fastjet::PseudoJet::plain_distance(), and fastjet::JetDefinition::DefaultRecombiner::recombine().

const unsigned int fastjet::twopow31 = 2147483648U
 

Definition at line 40 of file ClosestPair2D.cc.

Referenced by fastjet::ClosestPair2D::_initialize(), and fastjet::ClosestPair2D::_point2shuffle().

const double fastjet::zeta2 = 1.644934066848226436472415166646025189219
 

Definition at line 45 of file numconsts.hh.

const double fastjet::zeta3 = 1.202056903159594285399738161511449990765
 

Definition at line 46 of file numconsts.hh.


Generated on Mon Apr 2 20:58:16 2007 for fastjet by  doxygen 1.4.2