FastJet  3.1.3
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Modules Pages
Classes | Public Types | Public Member Functions | Protected Attributes | Friends | List of all members
fastjet::PseudoJet Class Reference

Class to contain pseudojets, including minimal information of use to jet-clustering routines. More...

#include <fastjet/PseudoJet.hh>

Collaboration diagram for fastjet::PseudoJet:
Collaboration graph
[legend]

Classes

class  InexistentUserInfo
 error class to be thrown if accessing user info when it doesn't exist More...
 
class  UserInfoBase
 a base class to hold extra user information in a PseudoJet More...
 

Public Types

enum  {
  X =0, Y =1, Z =2, T =3,
  NUM_COORDINATES =4, SIZE =NUM_COORDINATES
}
 

Public Member Functions

template<>
 PseudoJet (const siscone::Cmomentum &four_vector)
 shortcut for converting siscone Cmomentum into PseudoJet
 
template<>
 PseudoJet (const siscone_spherical::CSphmomentum &four_vector)
 shortcut for converting siscone CSphmomentum into PseudoJet
 
Constructors and destructor
 PseudoJet ()
 default constructor, which as of FJ3.0 provides an object for which all operations are now valid and which has zero momentum
 
 PseudoJet (const double px, const double py, const double pz, const double E)
 construct a pseudojet from explicit components
 
template<class L >
 PseudoJet (const L &some_four_vector)
 constructor from any object that has px,py,pz,E = some_four_vector[0–3],
 
 PseudoJet (bool)
 
virtual ~PseudoJet ()
 default (virtual) destructor
 
Kinematic access functions
double E () const
 
double e () const
 
double px () const
 
double py () const
 
double pz () const
 
double phi () const
 returns phi (in the range 0..2pi)
 
double phi_std () const
 returns phi in the range -pi..pi
 
double phi_02pi () const
 returns phi in the range 0..2pi
 
double rap () const
 returns the rapidity or some large value when the rapidity is infinite
 
double rapidity () const
 the same as rap()
 
double pseudorapidity () const
 returns the pseudo-rapidity or some large value when the rapidity is infinite
 
double eta () const
 
double pt2 () const
 returns the squared transverse momentum
 
double pt () const
 returns the scalar transverse momentum
 
double perp2 () const
 returns the squared transverse momentum
 
double perp () const
 returns the scalar transverse momentum
 
double kt2 () const
 returns the squared transverse momentum
 
double m2 () const
 returns the squared invariant mass // like CLHEP
 
double m () const
 returns the invariant mass (If m2() is negative then -sqrt(-m2()) is returned, as in CLHEP)
 
double mperp2 () const
 returns the squared transverse mass = kt^2+m^2
 
double mperp () const
 returns the transverse mass = sqrt(kt^2+m^2)
 
double mt2 () const
 returns the squared transverse mass = kt^2+m^2
 
double mt () const
 returns the transverse mass = sqrt(kt^2+m^2)
 
double modp2 () const
 return the squared 3-vector modulus = px^2+py^2+pz^2
 
double modp () const
 return the 3-vector modulus = sqrt(px^2+py^2+pz^2)
 
double Et () const
 return the transverse energy
 
double Et2 () const
 return the transverse energy squared
 
double operator() (int i) const
 returns component i, where X==0, Y==1, Z==2, E==3
 
double operator[] (int i) const
 returns component i, where X==0, Y==1, Z==2, E==3
 
double kt_distance (const PseudoJet &other) const
 returns kt distance (R=1) between this jet and another
 
double plain_distance (const PseudoJet &other) const
 returns squared cylinder (rap-phi) distance between this jet and another
 
double squared_distance (const PseudoJet &other) const
 returns squared cylinder (rap-phi) distance between this jet and another
 
double delta_R (const PseudoJet &other) const
 return the cylinder (rap-phi) distance between this jet and another, $\Delta_R = \sqrt{\Delta y^2 + \Delta \phi^2}$. More...
 
double delta_phi_to (const PseudoJet &other) const
 returns other.phi() - this.phi(), constrained to be in range -pi . More...
 
double beam_distance () const
 returns distance between this jet and the beam
 
std::valarray< double > four_mom () const
 return a valarray containing the four-momentum (components 0-2 are 3-mom, component 3 is energy). More...
 
Kinematic modification functions
PseudoJetboost (const PseudoJet &prest)
 transform this jet (given in the rest frame of prest) into a jet in the lab frame [NOT FULLY TESTED] More...
 
PseudoJetunboost (const PseudoJet &prest)
 transform this jet (given in lab) into a jet in the rest frame of prest [NOT FULLY TESTED] More...
 
void operator*= (double)
 multiply the jet's momentum by the coefficient
 
void operator/= (double)
 divide the jet's momentum by the coefficient
 
void operator+= (const PseudoJet &)
 add the other jet's momentum to this jet
 
void operator-= (const PseudoJet &)
 subtract the other jet's momentum from this jet
 
void reset (double px, double py, double pz, double E)
 reset the 4-momentum according to the supplied components and put the user and history indices back to their default values
 
void reset (const PseudoJet &psjet)
 reset the PseudoJet to be equal to psjet (including its indices); NB if the argument is derived from a PseudoJet then the "reset" used will be the templated version More...
 
template<class L >
void reset (const L &some_four_vector)
 reset the 4-momentum according to the supplied generic 4-vector (accessible via indexing, [0]==px,...[3]==E) and put the user and history indices back to their default values. More...
 
void reset_PtYPhiM (double pt_in, double y_in, double phi_in, double m_in=0.0)
 reset the PseudoJet according to the specified pt, rapidity, azimuth and mass (also resetting indices, etc.) (phi should satisfy -2pi<phi<4pi)
 
void reset_momentum (double px, double py, double pz, double E)
 reset the 4-momentum according to the supplied components but leave all other information (indices, user info, etc.) untouched
 
void reset_momentum (const PseudoJet &pj)
 reset the 4-momentum according to the components of the supplied PseudoJet, including cached components; note that the template version (below) will be called for classes derived from PJ. More...
 
void reset_momentum_PtYPhiM (double pt, double y, double phi, double m=0.0)
 reset the 4-momentum according to the specified pt, rapidity, azimuth and mass (phi should satisfy -2pi<phi<4pi)
 
template<class L >
void reset_momentum (const L &some_four_vector)
 reset the 4-momentum according to the supplied generic 4-vector (accessible via indexing, [0]==px,...[3]==E), but leave all other information (indices, user info, etc.) untouched
 
void set_cached_rap_phi (double rap, double phi)
 in some cases when setting a 4-momentum, the user/program knows what rapidity and azimuth are associated with that 4-momentum; by calling this routine the user can provide the information directly to the PseudoJet and avoid expensive rap-phi recalculations. More...
 
User index functions

To allow the user to set and access an integer index which can be exploited by the user to associate extra information with a particle/jet (for example pdg id, or an indication of a particle's origin within the user's analysis)

int user_index () const
 return the user_index,
 
void set_user_index (const int index)
 set the user_index, intended to allow the user to add simple identifying information to a particle/jet
 
User information types and functions

Allows PseudoJet to carry extra user info (as an object derived from UserInfoBase).

void set_user_info (UserInfoBase *user_info_in)
 sets the internal shared pointer to the user information. More...
 
template<class L >
const L & user_info () const
 returns a reference to the dynamic cast conversion of user_info to type L. More...
 
bool has_user_info () const
 returns true if the PseudoJet has user information
 
template<class L >
bool has_user_info () const
 returns true if the PseudoJet has user information than can be cast to the template argument type. More...
 
const UserInfoBaseuser_info_ptr () const
 retrieve a pointer to the (const) user information
 
const SharedPtr< UserInfoBase > & user_info_shared_ptr () const
 retrieve a (const) shared pointer to the user information
 
SharedPtr< UserInfoBase > & user_info_shared_ptr ()
 retrieve a (non-const) shared pointer to the user information; you can use this, for example, to set the shared pointer, eg More...
 
Description

Since a PseudoJet can have a structure that contains a variety of information, we provide a description that allows one to check exactly what kind of PseudoJet we are dealing with

std::string description () const
 return a string describing what kind of PseudoJet we are dealing with
 
Access to the associated ClusterSequence object.

In addition to having kinematic information, jets may contain a reference to an associated ClusterSequence (this is the case, for example, if the jet has been returned by a ClusterSequence member function).

bool has_associated_cluster_sequence () const
 returns true if this PseudoJet has an associated ClusterSequence.
 
bool has_associated_cs () const
 shorthand for has_associated_cluster_sequence()
 
bool has_valid_cluster_sequence () const
 returns true if this PseudoJet has an associated and still valid(ated) ClusterSequence. More...
 
bool has_valid_cs () const
 shorthand for has_valid_cluster_sequence()
 
const ClusterSequenceassociated_cluster_sequence () const
 get a (const) pointer to the parent ClusterSequence (NULL if inexistent)
 
const ClusterSequenceassociated_cs () const
 
const ClusterSequencevalidated_cluster_sequence () const
 if the jet has a valid associated cluster sequence then return a pointer to it; otherwise throw an error
 
const ClusterSequencevalidated_cs () const
 shorthand for validated_cluster_sequence()
 
const ClusterSequenceAreaBasevalidated_cluster_sequence_area_base () const
 if the jet has valid area information then return a pointer to the associated ClusterSequenceAreaBase object; otherwise throw an error
 
const ClusterSequenceAreaBasevalidated_csab () const
 shorthand for validated_cluster_sequence_area_base()
 
Access to the associated PseudoJetStructureBase object.

In addition to having kinematic information, jets may contain a reference to an associated ClusterSequence (this is the case, for example, if the jet has been returned by a ClusterSequence member function).

void set_structure_shared_ptr (const SharedPtr< PseudoJetStructureBase > &structure)
 set the associated structure
 
bool has_structure () const
 return true if there is some structure associated with this PseudoJet
 
const PseudoJetStructureBasestructure_ptr () const
 return a pointer to the structure (of type PseudoJetStructureBase*) associated with this PseudoJet. More...
 
PseudoJetStructureBasestructure_non_const_ptr ()
 return a non-const pointer to the structure (of type PseudoJetStructureBase*) associated with this PseudoJet. More...
 
const PseudoJetStructureBasevalidated_structure_ptr () const
 return a pointer to the structure (of type PseudoJetStructureBase*) associated with this PseudoJet. More...
 
const SharedPtr
< PseudoJetStructureBase > & 
structure_shared_ptr () const
 return a reference to the shared pointer to the PseudoJetStructureBase associated with this PseudoJet
 
template<typename StructureType >
const StructureType & structure () const
 returns a reference to the structure casted to the requested structure type More...
 
template<typename TransformerType >
bool has_structure_of () const
 check if the PseudoJet has the structure resulting from a Transformer (that is, its structure is compatible with a Transformer::StructureType). More...
 
template<typename TransformerType >
const
TransformerType::StructureType & 
structure_of () const
 this is a helper to access any structure created by a Transformer (that is, of type Transformer::StructureType). More...
 
Methods for access to information about jet structure

These allow access to jet constituents, and other jet subtructure information.

They only work if the jet is associated with a ClusterSequence.

virtual bool has_partner (PseudoJet &partner) const
 check if it has been recombined with another PseudoJet in which case, return its partner through the argument. More...
 
virtual bool has_child (PseudoJet &child) const
 check if it has been recombined with another PseudoJet in which case, return its child through the argument. More...
 
virtual bool has_parents (PseudoJet &parent1, PseudoJet &parent2) const
 check if it is the product of a recombination, in which case return the 2 parents through the 'parent1' and 'parent2' arguments. More...
 
virtual bool contains (const PseudoJet &constituent) const
 check if the current PseudoJet contains the one passed as argument. More...
 
virtual bool is_inside (const PseudoJet &jet) const
 check if the current PseudoJet is contained the one passed as argument. More...
 
virtual bool has_constituents () const
 returns true if the PseudoJet has constituents
 
virtual std::vector< PseudoJetconstituents () const
 retrieve the constituents. More...
 
virtual bool has_exclusive_subjets () const
 returns true if the PseudoJet has support for exclusive subjets
 
std::vector< PseudoJetexclusive_subjets (const double dcut) const
 return a vector of all subjets of the current jet (in the sense of the exclusive algorithm) that would be obtained when running the algorithm with the given dcut. More...
 
int n_exclusive_subjets (const double dcut) const
 return the size of exclusive_subjets(...); still n ln n with same coefficient, but marginally more efficient than manually taking exclusive_subjets.size() More...
 
std::vector< PseudoJetexclusive_subjets (int nsub) const
 return the list of subjets obtained by unclustering the supplied jet down to nsub subjets. More...
 
std::vector< PseudoJetexclusive_subjets_up_to (int nsub) const
 return the list of subjets obtained by unclustering the supplied jet down to nsub subjets (or all constituents if there are fewer than nsub). More...
 
double exclusive_subdmerge (int nsub) const
 Returns the dij that was present in the merging nsub+1 -> nsub subjets inside this jet. More...
 
double exclusive_subdmerge_max (int nsub) const
 Returns the maximum dij that occurred in the whole event at the stage that the nsub+1 -> nsub merge of subjets occurred inside this jet. More...
 
virtual bool has_pieces () const
 returns true if a jet has pieces More...
 
virtual std::vector< PseudoJetpieces () const
 retrieve the pieces that make up the jet. More...
 
virtual bool has_area () const
 check if it has a defined area
 
virtual double area () const
 return the jet (scalar) area. More...
 
virtual double area_error () const
 return the error (uncertainty) associated with the determination of the area of this jet. More...
 
virtual PseudoJet area_4vector () const
 return the jet 4-vector area. More...
 
virtual bool is_pure_ghost () const
 true if this jet is made exclusively of ghosts. More...
 
Members mainly intended for internal use
int cluster_hist_index () const
 return the cluster_hist_index, intended to be used by clustering routines. More...
 
void set_cluster_hist_index (const int index)
 set the cluster_hist_index, intended to be used by clustering routines.
 
int cluster_sequence_history_index () const
 alternative name for cluster_hist_index() [perhaps more meaningful]
 
void set_cluster_sequence_history_index (const int index)
 alternative name for set_cluster_hist_index(...) [perhaps more meaningful]
 

Protected Attributes

SharedPtr< PseudoJetStructureBase_structure
 
SharedPtr< UserInfoBase_user_info
 

Friends

PseudoJet operator* (double, const PseudoJet &)
 

Detailed Description

Class to contain pseudojets, including minimal information of use to jet-clustering routines.

Definition at line 67 of file PseudoJet.hh.

Member Function Documentation

double fastjet::PseudoJet::delta_R ( const PseudoJet other) const
inline

return the cylinder (rap-phi) distance between this jet and another, $\Delta_R = \sqrt{\Delta y^2 + \Delta \phi^2}$.

Definition at line 191 of file PseudoJet.hh.

double fastjet::PseudoJet::delta_phi_to ( const PseudoJet other) const

returns other.phi() - this.phi(), constrained to be in range -pi .

returns other.phi() - this.phi(), i.e.

. pi

the phi distance to other, constrained to be in range -pi .. pi

Definition at line 401 of file PseudoJet.cc.

valarray< double > fastjet::PseudoJet::four_mom ( ) const

return a valarray containing the four-momentum (components 0-2 are 3-mom, component 3 is energy).

Definition at line 116 of file PseudoJet.cc.

PseudoJet & fastjet::PseudoJet::boost ( const PseudoJet prest)

transform this jet (given in the rest frame of prest) into a jet in the lab frame [NOT FULLY TESTED]

transform this jet (given in the rest frame of prest) into a jet in the lab frame

Definition at line 282 of file PseudoJet.cc.

PseudoJet & fastjet::PseudoJet::unboost ( const PseudoJet prest)

transform this jet (given in lab) into a jet in the rest frame of prest [NOT FULLY TESTED]

transform this jet (given in lab) into a jet in the rest frame of prest

Definition at line 309 of file PseudoJet.cc.

void fastjet::PseudoJet::reset ( const PseudoJet psjet)
inline

reset the PseudoJet to be equal to psjet (including its indices); NB if the argument is derived from a PseudoJet then the "reset" used will be the templated version

Note: this is included on top of the templated version because PseudoJet is not "derived" from PseudoJet, so the templated reset would not handle this case properly.

Definition at line 244 of file PseudoJet.hh.

template<class L >
void fastjet::PseudoJet::reset ( const L &  some_four_vector)
inline

reset the 4-momentum according to the supplied generic 4-vector (accessible via indexing, [0]==px,...[3]==E) and put the user and history indices back to their default values.

Definition at line 251 of file PseudoJet.hh.

void fastjet::PseudoJet::reset_momentum ( const PseudoJet pj)
inline

reset the 4-momentum according to the components of the supplied PseudoJet, including cached components; note that the template version (below) will be called for classes derived from PJ.

Definition at line 943 of file PseudoJet.hh.

void fastjet::PseudoJet::set_cached_rap_phi ( double  rap,
double  phi 
)

in some cases when setting a 4-momentum, the user/program knows what rapidity and azimuth are associated with that 4-momentum; by calling this routine the user can provide the information directly to the PseudoJet and avoid expensive rap-phi recalculations.

  • Parameters
    raprapidity
  • Parameters
    phi(in range -twopi...4*pi)
    USE WITH CAUTION: there are no checks that the rapidity and azimuth supplied are sensible, nor does this reset the 4-momentum components if things don't match.

Definition at line 340 of file PseudoJet.cc.

void fastjet::PseudoJet::set_user_info ( UserInfoBase user_info_in)
inline

sets the internal shared pointer to the user information.

Note that the PseudoJet will now own the pointer, and delete the corresponding object when it (the jet, and any copies of the jet) goes out of scope.

Definition at line 387 of file PseudoJet.hh.

template<class L >
const L& fastjet::PseudoJet::user_info ( ) const
inline

returns a reference to the dynamic cast conversion of user_info to type L.

Usage: suppose you have previously set the user info with a pointer to an object of type MyInfo,

class MyInfo: public PseudoJet::UserInfoBase { MyInfo(int id) : _pdg_id(id); int pdg_id() const {return _pdg_id;} int _pdg_id; };

PseudoJet particle(...); particle.set_user_info(new MyInfo(its_pdg_id));

Then you would access that pdg_id() as

particle.user_info<MyInfo>().pdg_id();

It's overkill for just a single integer, but scales easily to more extensive information.

Note that user_info() throws an InexistentUserInfo() error if there is no user info; throws a std::bad_cast if the conversion doesn't work

If this behaviour does not fit your needs, use instead the the user_info_ptr() or user_info_shared_ptr() member functions.

Definition at line 420 of file PseudoJet.hh.

template<class L >
bool fastjet::PseudoJet::has_user_info ( ) const
inline

returns true if the PseudoJet has user information than can be cast to the template argument type.

Definition at line 433 of file PseudoJet.hh.

SharedPtr<UserInfoBase>& fastjet::PseudoJet::user_info_shared_ptr ( )
inline

retrieve a (non-const) shared pointer to the user information; you can use this, for example, to set the shared pointer, eg

p2.user_info_shared_ptr() = p1.user_info_shared_ptr();

or

SharedPtr<PseudoJet::UserInfoBase> info_shared(new MyInfo(...));
p2.user_info_shared_ptr() = info_shared;

Definition at line 462 of file PseudoJet.hh.

bool fastjet::PseudoJet::has_valid_cluster_sequence ( ) const

returns true if this PseudoJet has an associated and still valid(ated) ClusterSequence.

Definition at line 447 of file PseudoJet.cc.

const PseudoJetStructureBase * fastjet::PseudoJet::structure_ptr ( ) const

return a pointer to the structure (of type PseudoJetStructureBase*) associated with this PseudoJet.

return NULL if there is no associated structure

Definition at line 479 of file PseudoJet.cc.

PseudoJetStructureBase * fastjet::PseudoJet::structure_non_const_ptr ( )

return a non-const pointer to the structure (of type PseudoJetStructureBase*) associated with this PseudoJet.

return NULL if there is no associated structure

Only use this if you know what you are doing. In any case, prefer the 'structure_ptr()' (the const version) to this method, unless you really need a write access to the PseudoJet's underlying structure.

Definition at line 494 of file PseudoJet.cc.

const PseudoJetStructureBase * fastjet::PseudoJet::validated_structure_ptr ( ) const

return a pointer to the structure (of type PseudoJetStructureBase*) associated with this PseudoJet.

throw an error if there is no associated structure

Definition at line 504 of file PseudoJet.cc.

template<typename StructureType >
const StructureType & fastjet::PseudoJet::structure ( ) const

returns a reference to the structure casted to the requested structure type

If there is no structure associated, an Error is thrown. If the type is not met, a std::bad_cast error is thrown.

Definition at line 963 of file PseudoJet.hh.

template<typename TransformerType >
bool fastjet::PseudoJet::has_structure_of ( ) const

check if the PseudoJet has the structure resulting from a Transformer (that is, its structure is compatible with a Transformer::StructureType).

If there is no structure, false is returned.

Definition at line 971 of file PseudoJet.hh.

template<typename TransformerType >
const TransformerType::StructureType & fastjet::PseudoJet::structure_of ( ) const

this is a helper to access any structure created by a Transformer (that is, of type Transformer::StructureType).

If there is no structure, or if the structure is not compatible with TransformerType, an error is thrown.

Definition at line 981 of file PseudoJet.hh.

bool fastjet::PseudoJet::has_partner ( PseudoJet partner) const
virtual

check if it has been recombined with another PseudoJet in which case, return its partner through the argument.

Otherwise, 'partner' is set to 0.

an Error is thrown if this PseudoJet has no currently valid associated ClusterSequence

Definition at line 525 of file PseudoJet.cc.

bool fastjet::PseudoJet::has_child ( PseudoJet child) const
virtual

check if it has been recombined with another PseudoJet in which case, return its child through the argument.

Otherwise, 'child' is set to 0.

an Error is thrown if this PseudoJet has no currently valid associated ClusterSequence

Definition at line 536 of file PseudoJet.cc.

bool fastjet::PseudoJet::has_parents ( PseudoJet parent1,
PseudoJet parent2 
) const
virtual

check if it is the product of a recombination, in which case return the 2 parents through the 'parent1' and 'parent2' arguments.

Otherwise, set these to 0.

an Error is thrown if this PseudoJet has no currently valid associated ClusterSequence

Definition at line 547 of file PseudoJet.cc.

bool fastjet::PseudoJet::contains ( const PseudoJet constituent) const
virtual

check if the current PseudoJet contains the one passed as argument.

an Error is thrown if this PseudoJet has no currently valid associated ClusterSequence

Definition at line 557 of file PseudoJet.cc.

bool fastjet::PseudoJet::is_inside ( const PseudoJet jet) const
virtual

check if the current PseudoJet is contained the one passed as argument.

an Error is thrown if this PseudoJet has no currently valid associated ClusterSequence

Definition at line 567 of file PseudoJet.cc.

vector< PseudoJet > fastjet::PseudoJet::constituents ( ) const
virtual

retrieve the constituents.

an Error is thrown if this PseudoJet has no currently valid associated ClusterSequence or other substructure information

Definition at line 580 of file PseudoJet.cc.

std::vector< PseudoJet > fastjet::PseudoJet::exclusive_subjets ( const double  dcut) const

return a vector of all subjets of the current jet (in the sense of the exclusive algorithm) that would be obtained when running the algorithm with the given dcut.

Time taken is O(m ln m), where m is the number of subjets that are found. If m gets to be of order of the total number of constituents in the jet, this could be substantially slower than just getting that list of constituents.

an Error is thrown if this PseudoJet has no currently valid associated ClusterSequence

Definition at line 603 of file PseudoJet.cc.

int fastjet::PseudoJet::n_exclusive_subjets ( const double  dcut) const

return the size of exclusive_subjets(...); still n ln n with same coefficient, but marginally more efficient than manually taking exclusive_subjets.size()

an Error is thrown if this PseudoJet has no currently valid associated ClusterSequence

Definition at line 614 of file PseudoJet.cc.

std::vector< PseudoJet > fastjet::PseudoJet::exclusive_subjets ( int  nsub) const

return the list of subjets obtained by unclustering the supplied jet down to nsub subjets.

Throws an error if there are fewer than nsub particles in the jet.

For ClusterSequence type jets, requires nsub ln nsub time

An Error is thrown if this PseudoJet has no currently valid associated ClusterSequence

Definition at line 634 of file PseudoJet.cc.

std::vector< PseudoJet > fastjet::PseudoJet::exclusive_subjets_up_to ( int  nsub) const

return the list of subjets obtained by unclustering the supplied jet down to nsub subjets (or all constituents if there are fewer than nsub).

For ClusterSequence type jets, requires nsub ln nsub time

An Error is thrown if this PseudoJet has no currently valid associated ClusterSequence

Definition at line 627 of file PseudoJet.cc.

double fastjet::PseudoJet::exclusive_subdmerge ( int  nsub) const

Returns the dij that was present in the merging nsub+1 -> nsub subjets inside this jet.

Returns 0 if there were nsub or fewer constituents in the jet.

an Error is thrown if this PseudoJet has no currently valid associated ClusterSequence

Definition at line 651 of file PseudoJet.cc.

double fastjet::PseudoJet::exclusive_subdmerge_max ( int  nsub) const

Returns the maximum dij that occurred in the whole event at the stage that the nsub+1 -> nsub merge of subjets occurred inside this jet.

Returns 0 if there were nsub or fewer constituents in the jet.

an Error is thrown if this PseudoJet has no currently valid associated ClusterSequence

Definition at line 662 of file PseudoJet.cc.

bool fastjet::PseudoJet::has_pieces ( ) const
virtual

returns true if a jet has pieces

By default a single particle or a jet coming from a ClusterSequence have no pieces and this methos will return false.

In practice, this is equivalent to have an structure of type CompositeJetStructure.

Definition at line 671 of file PseudoJet.cc.

std::vector< PseudoJet > fastjet::PseudoJet::pieces ( ) const
virtual

retrieve the pieces that make up the jet.

If the jet does not support pieces, an error is throw

Definition at line 680 of file PseudoJet.cc.

double fastjet::PseudoJet::area ( ) const
virtual

return the jet (scalar) area.

throws an Error if there is no support for area in the parent CS

Definition at line 717 of file PseudoJet.cc.

double fastjet::PseudoJet::area_error ( ) const
virtual

return the error (uncertainty) associated with the determination of the area of this jet.

throws an Error if there is no support for area in the parent CS

Definition at line 725 of file PseudoJet.cc.

PseudoJet fastjet::PseudoJet::area_4vector ( ) const
virtual

return the jet 4-vector area.

throws an Error if there is no support for area in the parent CS

Definition at line 732 of file PseudoJet.cc.

bool fastjet::PseudoJet::is_pure_ghost ( ) const
virtual

true if this jet is made exclusively of ghosts.

throws an Error if there is no support for area in the parent CS

Definition at line 739 of file PseudoJet.cc.

int fastjet::PseudoJet::cluster_hist_index ( ) const
inline

return the cluster_hist_index, intended to be used by clustering routines.

Definition at line 770 of file PseudoJet.hh.


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