Class to contain pseudojets, including minimal information of use to to jet-clustering routines. More...
#include <PseudoJet.hh>
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 leaves PseudoJet unusable | |
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], | |
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 | kt2 () const |
returns the squared transverse momentum | |
double | perp2 () const |
returns the squared transverse momentum | |
double | perp () const |
returns the scalar transverse momentum | |
double | m2 () const |
returns the squared invariant mass // like 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 | m () const |
returns the invariant mass (If m2() is negative then -sqrt(-m2()) is returned, as in CLHEP) | |
double | modp2 () const |
return 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_phi_to (const PseudoJet &other) const |
returns other.phi() - this.phi(), constrained to be in range -pi . | |
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). | |
Kinematic modification functions | |
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] | |
PseudoJet & | unboost (const PseudoJet &prest) |
transform this jet (given in lab) into a jet in the rest frame of prest [NOT FULLY TESTED] | |
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 (which does not know about indices...) | |
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. | |
User index functions | |
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. | |
template<class L > | |
const L & | user_info () const |
returns a reference to the dynamic cast conversion of user_info to type L. | |
const UserInfoBase * | user_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 | |
Description | |
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 (and still valid) ClusterSequence. | |
const ClusterSequence * | associated_cluster_sequence () const |
get a (const) pointer to the parent ClusterSequence (NULL if inexistent) | |
const ClusterSequence * | validated_cs () const |
if the jet has a valid associated cluster sequence then return a pointer to it; otherwise throw an error | |
const ClusterSequenceAreaBase * | validated_csab () const |
if the jet has valid area information then return a pointer to the associated ClusterSequenceAreaBase object; otherwise throw an error | |
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 strusture associated with this PseudoJet | |
const PseudoJetStructureBase * | structure_ptr () const |
return a pointer to the structure (of type PseudoJetStructureBase*) associated wioth this PseudoJet. | |
const PseudoJetStructureBase * | validated_structure_ptr () const |
return a pointer to the structure (of type PseudoJetStructureBase*) associated wioth this PseudoJet. | |
const SharedPtr < PseudoJetStructureBase > & | structure_shared_ptr () const |
return a reference to the shared pointer to the PseudoJetStructureBase associated wioth this PseudoJet | |
template<typename StructureType > | |
const StructureType & | structure () const |
returns a reference to the structure casted to the requested structure type | |
template<typename TransformerType > | |
bool | has_properties_of () const |
check if the PseudoJet has the properties of the result of a Transformer (that is, its structure is compatible with a Transformer::StructureType) if there is no structure, false is returned | |
template<typename TransformerType > | |
const TransformerType::StructureType * | extra_properties () const |
this is a helper to access an structuree created by a Transformer (that is, of type Transformer::StructureType) NULL is returned if the corresponding type is not met if there is no structure, an error is thrown | |
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. | |
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. | |
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. | |
virtual bool | contains (const PseudoJet &constituent) const |
check if the current PseudoJet contains the one passed as argument. | |
virtual bool | is_inside (const PseudoJet &jet) const |
check if the current PseudoJet is contained the one passed as argument. | |
virtual bool | has_constituents () const |
returns true if the PseudoJet has constituents | |
virtual std::vector< PseudoJet > | constituents () const |
retrieve the constituents. | |
virtual bool | has_exclusive_subjets () const |
returns true if the PseudoJet has support for exclusive subjets | |
std::vector< 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. | |
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() | |
std::vector< PseudoJet > | exclusive_subjets (int nsub) const |
return the list of subjets obtained by unclustering the supplied jet down to n subjets (or all constituents if there are fewer than n). | |
double | exclusive_subdmerge (int nsub) const |
return the dij that was present in the merging nsub+1 -> nsub subjets inside this jet. | |
double | exclusive_subdmerge_max (int nsub) const |
return the maximum dij that occurred in the whole event at the stage that the nsub+1 -> nsub merge of subjets occurred inside this jet. | |
virtual bool | has_pieces () const |
returns true if a jet has pieces | |
virtual std::vector< PseudoJet > | pieces () const |
retrieve the pieces that make up the jet. | |
virtual bool | has_area () const |
check if it has a defined area | |
virtual double | area () const |
return the jet (scalar) area. | |
virtual double | area_error () const |
return the error (uncertainty) associated with the determination of the area of this jet. | |
virtual PseudoJet | area_4vector () const |
return the jet 4-vector area. | |
virtual bool | is_pure_ghost () const |
true if this jet is made exclusively of ghosts. | |
Members mainly intended for internal use | |
int | cluster_hist_index () const |
return the cluster_hist_index, intended to be used by clustering routines. | |
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 |
Class to contain pseudojets, including minimal information of use to to jet-clustering routines.
Definition at line 65 of file PseudoJet.hh.
double fastjet::PseudoJet::m | ( | ) | const [inline] |
returns the invariant mass (If m2() is negative then -sqrt(-m2()) is returned, as in CLHEP)
specialization of the "reset" template for case where something is reset to a pseudojet -- it then takes the user and history indices from the psjet
Definition at line 753 of file PseudoJet.hh.
References m2().
Referenced by boost(), operator<<(), and unboost().
{ double mm = m2(); return mm < 0.0 ? -std::sqrt(-mm) : std::sqrt(mm); }
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 326 of file PseudoJet.cc.
References phi().
{ double dphi = other.phi() - phi(); if (dphi > pi) dphi -= twopi; if (dphi < -pi) dphi += twopi; return dphi; }
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 106 of file PseudoJet.cc.
{
valarray<double> mom(4);
mom[0] = _px;
mom[1] = _py;
mom[2] = _pz;
mom[3] = _E ;
return mom;
}
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 lab) into a jet in the rest frame of prest
Definition at line 235 of file PseudoJet.cc.
References m().
{ if (prest.px() == 0.0 && prest.py() == 0.0 && prest.pz() == 0.0) return *this; double m = prest.m(); assert(m != 0); double pf4 = ( px()*prest.px() + py()*prest.py() + pz()*prest.pz() + E()*prest.E() )/m; double fn = (pf4 + E()) / (prest.E() + m); _px += fn*prest.px(); _py += fn*prest.py(); _pz += fn*prest.pz(); _E = pf4; _finish_init(); // we need to recalculate phi,rap,kt2 return *this; }
transform this jet (given in lab) into a jet in the rest frame of prest [NOT FULLY TESTED]
transform this jet (given in the rest frame of prest) into a jet in the lab frame;
Definition at line 262 of file PseudoJet.cc.
References m().
{ if (prest.px() == 0.0 && prest.py() == 0.0 && prest.pz() == 0.0) return *this; double m = prest.m(); assert(m != 0); double pf4 = ( -px()*prest.px() - py()*prest.py() - pz()*prest.pz() + E()*prest.E() )/m; double fn = (pf4 + E()) / (prest.E() + m); _px -= fn*prest.px(); _py -= fn*prest.py(); _pz -= fn*prest.pz(); _E = pf4; _finish_init(); // we need to recalculate phi,rap,kt2 return *this; }
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 213 of file PseudoJet.hh.
{ reset(some_four_vector[0], some_four_vector[1], some_four_vector[2], some_four_vector[3]); }
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 287 of file PseudoJet.hh.
{ _user_info.reset(user_info_in); }
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 320 of file PseudoJet.hh.
{ if (_user_info.get() == 0) throw InexistentUserInfo(); return dynamic_cast<const L &>(* _user_info.get()); }
bool fastjet::PseudoJet::has_associated_cluster_sequence | ( | ) | const |
returns true if this PseudoJet has an associated (and still valid) ClusterSequence.
Definition at line 355 of file PseudoJet.cc.
Referenced by fastjet::Filter::_recursively_check_ca(), associated_cluster_sequence(), has_area(), and fastjet::ClusterSequenceStructure::object_in_jet().
{
return (_structure()) && (_structure->has_associated_cluster_sequence());
}
const PseudoJetStructureBase * fastjet::PseudoJet::structure_ptr | ( | ) | const |
return a pointer to the structure (of type PseudoJetStructureBase*) associated wioth this PseudoJet.
return NULL if there is no associated structure
Definition at line 397 of file PseudoJet.cc.
{ if (!_structure()) return NULL; return _structure(); }
const PseudoJetStructureBase * fastjet::PseudoJet::validated_structure_ptr | ( | ) | const |
return a pointer to the structure (of type PseudoJetStructureBase*) associated wioth this PseudoJet.
throw an error if there is no associated structure
Definition at line 407 of file PseudoJet.cc.
Referenced by area(), area_4vector(), area_error(), constituents(), contains(), exclusive_subdmerge(), exclusive_subdmerge_max(), exclusive_subjets(), has_area(), has_child(), has_parents(), has_partner(), is_inside(), is_pure_ghost(), n_exclusive_subjets(), structure(), and validated_cs().
{ if (!_structure()) throw Error("Trying to access the structure of a PseudoJet which has no associated structure"); return _structure(); }
const StructureType & fastjet::PseudoJet::structure | ( | ) | const |
returns a reference to the structure casted to the requested structure type
If there is no sructure associated, an Error is thrown. If the type is not met, a std::bad_cast error is thrown.
Definition at line 779 of file PseudoJet.hh.
References validated_structure_ptr().
{ return dynamic_cast<const StructureType &>(* validated_structure_ptr()); }
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 428 of file PseudoJet.cc.
References fastjet::PseudoJetStructureBase::has_partner(), and validated_structure_ptr().
{ return validated_structure_ptr()->has_partner(*this, partner); }
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 439 of file PseudoJet.cc.
References fastjet::PseudoJetStructureBase::has_child(), and validated_structure_ptr().
{ return validated_structure_ptr()->has_child(*this, child); }
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 450 of file PseudoJet.cc.
References fastjet::PseudoJetStructureBase::has_parents(), and validated_structure_ptr().
{ return validated_structure_ptr()->has_parents(*this, parent1, parent2); }
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 460 of file PseudoJet.cc.
References fastjet::PseudoJetStructureBase::object_in_jet(), and validated_structure_ptr().
{ return validated_structure_ptr()->object_in_jet(constituent, *this); }
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 470 of file PseudoJet.cc.
References fastjet::PseudoJetStructureBase::object_in_jet(), and validated_structure_ptr().
{ return validated_structure_ptr()->object_in_jet(*this, jet); }
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 483 of file PseudoJet.cc.
References fastjet::PseudoJetStructureBase::constituents(), and validated_structure_ptr().
{ return validated_structure_ptr()->constituents(*this); }
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 506 of file PseudoJet.cc.
References fastjet::PseudoJetStructureBase::exclusive_subjets(), and validated_structure_ptr().
{ return validated_structure_ptr()->exclusive_subjets(*this, dcut); }
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 517 of file PseudoJet.cc.
References fastjet::PseudoJetStructureBase::n_exclusive_subjets(), and validated_structure_ptr().
{ return validated_structure_ptr()->n_exclusive_subjets(*this, dcut); }
std::vector< PseudoJet > fastjet::PseudoJet::exclusive_subjets | ( | int | nsub | ) | const |
return the list of subjets obtained by unclustering the supplied jet down to n subjets (or all constituents if there are fewer than n).
requires n ln n time
an Error is thrown if this PseudoJet has no currently valid associated ClusterSequence
Definition at line 530 of file PseudoJet.cc.
References fastjet::PseudoJetStructureBase::exclusive_subjets(), and validated_structure_ptr().
{ return validated_structure_ptr()->exclusive_subjets(*this, nsub); }
double fastjet::PseudoJet::exclusive_subdmerge | ( | int | nsub | ) | const |
return the dij that was present in the merging nsub+1 -> nsub subjets inside this jet.
an Error is thrown if this PseudoJet has no currently valid associated ClusterSequence
Definition at line 540 of file PseudoJet.cc.
References fastjet::PseudoJetStructureBase::exclusive_subdmerge(), and validated_structure_ptr().
{ return validated_structure_ptr()->exclusive_subdmerge(*this, nsub); }
double fastjet::PseudoJet::exclusive_subdmerge_max | ( | int | nsub | ) | const |
return the maximum dij that occurred in the whole event at the stage that the nsub+1 -> nsub merge of subjets occurred inside this jet.
an Error is thrown if this PseudoJet has no currently valid associated ClusterSequence
Definition at line 551 of file PseudoJet.cc.
References fastjet::PseudoJetStructureBase::exclusive_subdmerge_max(), and validated_structure_ptr().
{ return validated_structure_ptr()->exclusive_subdmerge_max(*this, nsub); }
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 560 of file PseudoJet.cc.
Referenced by fastjet::Filter::_recursively_check_ca(), and pieces().
{
return ((_structure()) && (_structure->has_pieces()));
}
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 569 of file PseudoJet.cc.
References has_pieces().
Referenced by fastjet::Filter::_recursively_check_ca(), and main().
{ if (!has_pieces()) throw Error("Trying to retrieve the pieces of a PseudoJet that has no support for pieces."); return _structure->pieces(*this); }
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 602 of file PseudoJet.cc.
References fastjet::PseudoJetStructureBase::area(), and validated_structure_ptr().
{ return validated_structure_ptr()->area(*this); }
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 610 of file PseudoJet.cc.
References fastjet::PseudoJetStructureBase::area_error(), and validated_structure_ptr().
{ return validated_structure_ptr()->area_error(*this); }
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 617 of file PseudoJet.cc.
References fastjet::PseudoJetStructureBase::area_4vector(), and validated_structure_ptr().
Referenced by fastjet::ClusterSequenceArea::area_4vector().
{ return validated_structure_ptr()->area_4vector(*this); }
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 624 of file PseudoJet.cc.
References fastjet::PseudoJetStructureBase::is_pure_ghost(), and validated_structure_ptr().
{ return validated_structure_ptr()->is_pure_ghost(*this); }
int fastjet::PseudoJet::cluster_hist_index | ( | ) | const [inline] |
return the cluster_hist_index, intended to be used by clustering routines.
Definition at line 600 of file PseudoJet.hh.
Referenced by fastjet::ClusterSequenceActiveAreaExplicitGhosts::area(), fastjet::ClusterSequenceVoronoiArea::area(), fastjet::ClusterSequenceActiveArea::area(), fastjet::ClusterSequenceActiveAreaExplicitGhosts::area_4vector(), fastjet::ClusterSequenceVoronoiArea::area_4vector(), fastjet::ClusterSequenceActiveArea::area_4vector(), fastjet::ClusterSequenceActiveArea::area_error(), fastjet::ClusterSequenceActiveAreaExplicitGhosts::is_pure_ghost(), fastjet::SISConeBaseExtras::pass(), and fastjet::ClusterSequenceAreaBase::subtracted_jet().
{return _cluster_hist_index;}