fastjet 2.4.5
|
Class to contain pseudojets, including minimal information of use to to jet-clustering routines. More...
#include <PseudoJet.hh>
Public Types | |
enum | { X = 0, Y = 1, Z = 2, T = 3, NUM_COORDINATES = 4, SIZE = NUM_COORDINATES } |
Public Member Functions | |
PseudoJet () | |
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], | |
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 | |
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] | |
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] | |
int | user_index () const |
return the user_index, intended to allow the user to "add" information | |
void | set_user_index (const int index) |
set the user_index, intended to allow the user to "add" information | |
std::valarray< double > | four_mom () const |
return a valarray containing the four-momentum (components 0-2 are 3-mom, component 3 is energy). | |
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 | |
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. | |
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 | |
Private Member Functions | |
void | _finish_init () |
calculate phi, rap, kt2 based on the 4-momentum components | |
void | _reset_indices () |
set the indices to default values | |
Private Attributes | |
double | _px |
double | _py |
double | _pz |
double | _E |
double | _phi |
double | _rap |
double | _kt2 |
int | _cluster_hist_index |
int | _user_index |
Class to contain pseudojets, including minimal information of use to to jet-clustering routines.
Definition at line 52 of file PseudoJet.hh.
anonymous enum |
Definition at line 123 of file PseudoJet.hh.
{ X=0, Y=1, Z=2, T=3, NUM_COORDINATES=4, SIZE=NUM_COORDINATES };
fastjet::PseudoJet::PseudoJet | ( | ) | [inline] |
Definition at line 55 of file PseudoJet.hh.
{};
fastjet::PseudoJet::PseudoJet | ( | const double | px, |
const double | py, | ||
const double | pz, | ||
const double | E | ||
) |
construct a pseudojet from explicit components
Definition at line 47 of file PseudoJet.cc.
{ _E = E ; _px = px; _py = py; _pz = pz; this->_finish_init(); // some default values for the history and user indices _reset_indices(); }
fastjet::PseudoJet::PseudoJet | ( | const L & | some_four_vector | ) | [inline] |
constructor from any object that has px,py,pz,E = some_four_vector[0--3],
Definition at line 286 of file PseudoJet.hh.
{ _px = some_four_vector[0]; _py = some_four_vector[1]; _pz = some_four_vector[2]; _E = some_four_vector[3]; _finish_init(); // some default values for these two indices _reset_indices(); }
fastjet::PseudoJet::PseudoJet | ( | const siscone::Cmomentum & | four_vector | ) |
shortcut for converting siscone Cmomentum into PseudoJet
Definition at line 19 of file SISConePlugin.cc.
{ (*this) = PseudoJet(four_vector.px,four_vector.py,four_vector.pz, four_vector.E); }
fastjet::PseudoJet::PseudoJet | ( | const siscone_spherical::CSphmomentum & | four_vector | ) |
shortcut for converting siscone CSphmomentum into PseudoJet
Definition at line 19 of file SISConeSphericalPlugin.cc.
{ (*this) = PseudoJet(four_vector.px,four_vector.py,four_vector.pz, four_vector.E); }
void fastjet::PseudoJet::_finish_init | ( | ) | [private] |
calculate phi, rap, kt2 based on the 4-momentum components
do standard end of initialisation
Definition at line 64 of file PseudoJet.cc.
References fastjet::MaxRap, and twopi.
{ _kt2 = this->px()*this->px() + this->py()*this->py(); if (_kt2 == 0.0) { _phi = 0.0; } else { _phi = atan2(this->py(),this->px()); } if (_phi < 0.0) {_phi += twopi;} if (_phi >= twopi) {_phi -= twopi;} // can happen if phi=-|eps<1e-15|? if (this->E() == abs(this->pz()) && _kt2 == 0) { // Point has infinite rapidity -- convert that into a very large // number, but in such a way that different 0-pt momenta will have // different rapidities (so as to lift the degeneracy between // them) [this can be relevant at parton-level] double MaxRapHere = MaxRap + abs(this->pz()); if (this->pz() >= 0.0) {_rap = MaxRapHere;} else {_rap = -MaxRapHere;} } else { // get the rapidity in a way that's modestly insensitive to roundoff // error when things pz,E are large (actually the best we can do without // explicit knowledge of mass) double effective_m2 = max(0.0,m2()); // force non tachyonic mass double E_plus_pz = _E + abs(_pz); // the safer of p+, p- // p+/p- = (p+ p-) / (p-)^2 = (kt^2+m^2)/(p-)^2 _rap = 0.5*log((_kt2 + effective_m2)/(E_plus_pz*E_plus_pz)); if (_pz > 0) {_rap = - _rap;} } //if (this->E() != abs(this->pz())) { // _rap = 0.5*log((this->E() + this->pz())/(this->E() - this->pz())); // } else { // // Overlapping points can give problems. Let's lift the degeneracy // // in case of multiple 0-pT points (can be found at parton-level) // double MaxRapHere = MaxRap + abs(this->pz()); // if (this->pz() >= 0.0) {_rap = MaxRapHere;} else {_rap = -MaxRapHere;} //} }
void fastjet::PseudoJet::_reset_indices | ( | ) | [inline, private] |
set the indices to default values
Definition at line 299 of file PseudoJet.hh.
{ set_cluster_hist_index(-1); set_user_index(-1); }
double fastjet::PseudoJet::beam_distance | ( | ) | const [inline] |
returns distance between this jet and the beam
Definition at line 177 of file PseudoJet.hh.
{return _kt2;}
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 234 of file PseudoJet.cc.
References E(), m(), px(), py(), and pz().
{ 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; }
int fastjet::PseudoJet::cluster_hist_index | ( | ) | const [inline] |
return the cluster_hist_index, intended to be used by clustering routines.
Definition at line 135 of file PseudoJet.hh.
Referenced by fastjet::ClusterSequence::_potentially_valid(), fastjet::ClusterSequence::add_constituents(), 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::ClusterSequence::get_subhist_set(), fastjet::ClusterSequence::has_child(), fastjet::ClusterSequence::has_parents(), fastjet::ClusterSequence::has_partner(), fastjet::ClusterSequenceActiveAreaExplicitGhosts::is_pure_ghost(), fastjet::ClusterSequence::object_in_jet(), fastjet::SISConeBaseExtras::pass(), and fastjet::ClusterSequenceAreaBase::subtracted_jet().
{return _cluster_hist_index;}
int fastjet::PseudoJet::cluster_sequence_history_index | ( | ) | const [inline] |
alternative name for cluster_hist_index() [perhaps more meaningful]
Definition at line 140 of file PseudoJet.hh.
{ return cluster_hist_index();}
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 325 of file PseudoJet.cc.
References fastjet::d0::inline_maths::phi(), phi(), fastjet::pi, and twopi.
double fastjet::PseudoJet::e | ( | ) | const [inline] |
Definition at line 67 of file PseudoJet.hh.
{return _E;} // like CLHEP
double fastjet::PseudoJet::E | ( | ) | const [inline] |
Definition at line 66 of file PseudoJet.hh.
Referenced by fastjet::ClusterSequenceActiveArea::_throw_unless_jets_have_same_perp_or_E(), boost(), fastjet::dot_product(), fastjet::have_same_momentum(), fastjet::JadeBriefJet::init(), fastjet::operator+(), fastjet::operator-(), fastjet::JetDefinition::DefaultRecombiner::preprocess(), fastjet::SISConeSphericalPlugin::run_clustering(), fastjet::SISConePlugin::run_clustering(), fastjet::ATLASConePlugin::run_clustering(), and unboost().
{return _E;}
double fastjet::PseudoJet::Et | ( | ) | const [inline] |
return the transverse energy
Definition at line 112 of file PseudoJet.hh.
Referenced by fastjet::CMSIterativeConePlugin::run_clustering().
double fastjet::PseudoJet::Et2 | ( | ) | const [inline] |
double fastjet::PseudoJet::eta | ( | ) | const [inline] |
Definition at line 92 of file PseudoJet.hh.
Referenced by fastjet::CMSIterativeConePlugin::run_clustering().
{return pseudorapidity();}
valarray< double > fastjet::PseudoJet::four_mom | ( | ) | const |
double fastjet::PseudoJet::kt2 | ( | ) | const [inline] |
returns the squared transverse momentum
Definition at line 95 of file PseudoJet.hh.
Referenced by fastjet::ClusterSequence::jet_scale_for_algorithm().
{return _kt2;}
double fastjet::PseudoJet::kt_distance | ( | const PseudoJet & | other | ) | const |
returns kt distance (R=1) between this jet and another
Definition at line 302 of file PseudoJet.cc.
References _kt2, _phi, _rap, fastjet::d0::inline_maths::min(), fastjet::pi, and twopi.
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 321 of file PseudoJet.hh.
Referenced by boost(), main(), operator<<(), and unboost().
{ double mm = m2(); return mm < 0.0 ? -std::sqrt(-mm) : std::sqrt(mm); }
double fastjet::PseudoJet::m2 | ( | ) | const [inline] |
double fastjet::PseudoJet::modp2 | ( | ) | const [inline] |
return px^2+py^2+pz^2
Definition at line 110 of file PseudoJet.hh.
Referenced by fastjet::JadeBriefJet::init(), and fastjet::EECamBriefJet::init().
double fastjet::PseudoJet::mperp | ( | ) | const [inline] |
returns the transverse mass = sqrt(kt^2+m^2)
Definition at line 105 of file PseudoJet.hh.
{return sqrt(std::abs(mperp2()));}
double fastjet::PseudoJet::mperp2 | ( | ) | const [inline] |
double fastjet::PseudoJet::operator() | ( | int | i | ) | const |
returns component i, where X==0, Y==1, Z==2, E==3
Definition at line 117 of file PseudoJet.cc.
void fastjet::PseudoJet::operator*= | ( | double | coeff | ) |
void fastjet::PseudoJet::operator+= | ( | const PseudoJet & | other_jet | ) |
void fastjet::PseudoJet::operator-= | ( | const PseudoJet & | other_jet | ) |
void fastjet::PseudoJet::operator/= | ( | double | coeff | ) |
divide the jet's momentum by the coefficient
Definition at line 201 of file PseudoJet.cc.
{ (*this) *= 1.0/coeff; }
double fastjet::PseudoJet::operator[] | ( | int | i | ) | const [inline] |
returns component i, where X==0, Y==1, Z==2, E==3
Definition at line 119 of file PseudoJet.hh.
{ return (*this)(i); }; // this too
double fastjet::PseudoJet::perp | ( | ) | const [inline] |
returns the scalar transverse momentum
Definition at line 99 of file PseudoJet.hh.
Referenced by operator<<(), print_jet(), fastjet::JetDefinition::DefaultRecombiner::recombine(), fastjet::ClusterSequenceAreaBase::subtracted_jet(), and fastjet::ClusterSequenceAreaBase::subtracted_pt().
{return sqrt(_kt2);} // like CLHEP
double fastjet::PseudoJet::perp2 | ( | ) | const [inline] |
returns the squared transverse momentum
Definition at line 97 of file PseudoJet.hh.
Referenced by fastjet::ClusterSequenceActiveArea::_throw_unless_jets_have_same_perp_or_E(), fastjet::ClusterSequence::has_parents(), fastjet::ClusterSequence::inclusive_jets(), fastjet::JetDefinition::DefaultRecombiner::preprocess(), and fastjet::JetDefinition::DefaultRecombiner::recombine().
{return _kt2;} // like CLHEP
double fastjet::PseudoJet::phi | ( | ) | const [inline] |
returns phi (in the range 0..2pi)
Definition at line 73 of file PseudoJet.hh.
Referenced by fastjet::CircularRange::CircularRange(), delta_phi_to(), fastjet::RangeDefinition::is_in_range(), operator<<(), print_jet(), fastjet::JetDefinition::DefaultRecombiner::recombine(), fastjet::CMSIterativeConePlugin::run_clustering(), and fastjet::RangeDefinition::set_position().
{return phi_02pi();}
double fastjet::PseudoJet::phi_02pi | ( | ) | const [inline] |
double fastjet::PseudoJet::phi_std | ( | ) | const [inline] |
returns phi in the range -pi..pi
Definition at line 76 of file PseudoJet.hh.
References fastjet::pi, and twopi.
double fastjet::PseudoJet::plain_distance | ( | const PseudoJet & | other | ) | const |
returns squared cylinder (rap-phi) distance between this jet and another
Definition at line 315 of file PseudoJet.cc.
References _phi, _rap, fastjet::pi, and twopi.
Referenced by fastjet::TrackJetPlugin::run_clustering().
double fastjet::PseudoJet::pseudorapidity | ( | ) | const |
returns the pseudo-rapidity or some large value when the rapidity is infinite
Definition at line 137 of file PseudoJet.cc.
References fastjet::MaxRap, and fastjet::pi.
double fastjet::PseudoJet::px | ( | ) | const [inline] |
Definition at line 68 of file PseudoJet.hh.
Referenced by fastjet::ClusterSequenceActiveArea::_throw_unless_jets_have_same_perp_or_E(), boost(), fastjet::dot_product(), fastjet::have_same_momentum(), fastjet::JadeBriefJet::init(), fastjet::EECamBriefJet::init(), fastjet::operator+(), fastjet::operator-(), fastjet::JetDefinition::DefaultRecombiner::preprocess(), fastjet::SISConeSphericalPlugin::run_clustering(), fastjet::SISConePlugin::run_clustering(), fastjet::ATLASConePlugin::run_clustering(), and unboost().
{return _px;}
double fastjet::PseudoJet::py | ( | ) | const [inline] |
Definition at line 69 of file PseudoJet.hh.
Referenced by fastjet::ClusterSequenceActiveArea::_throw_unless_jets_have_same_perp_or_E(), boost(), fastjet::dot_product(), fastjet::have_same_momentum(), fastjet::JadeBriefJet::init(), fastjet::EECamBriefJet::init(), fastjet::operator+(), fastjet::operator-(), fastjet::JetDefinition::DefaultRecombiner::preprocess(), fastjet::SISConeSphericalPlugin::run_clustering(), fastjet::SISConePlugin::run_clustering(), fastjet::ATLASConePlugin::run_clustering(), and unboost().
{return _py;}
double fastjet::PseudoJet::pz | ( | ) | const [inline] |
Definition at line 70 of file PseudoJet.hh.
Referenced by fastjet::ClusterSequenceActiveArea::_throw_unless_jets_have_same_perp_or_E(), boost(), fastjet::dot_product(), fastjet::have_same_momentum(), fastjet::JadeBriefJet::init(), fastjet::EECamBriefJet::init(), fastjet::operator+(), fastjet::operator-(), fastjet::JetDefinition::DefaultRecombiner::preprocess(), fastjet::SISConeSphericalPlugin::run_clustering(), fastjet::SISConePlugin::run_clustering(), fastjet::ATLASConePlugin::run_clustering(), and unboost().
{return _pz;}
double fastjet::PseudoJet::rap | ( | ) | const [inline] |
returns the rapidity or some large value when the rapidity is infinite
Definition at line 84 of file PseudoJet.hh.
Referenced by fastjet::ClusterSequenceActiveArea::_transfer_areas(), fastjet::CircularRange::CircularRange(), fastjet::RangeDefinition::is_in_range(), main(), operator<<(), print_jet(), fastjet::JetDefinition::DefaultRecombiner::recombine(), and fastjet::RangeDefinition::set_position().
{return _rap;}
double fastjet::PseudoJet::rapidity | ( | ) | const [inline] |
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 200 of file PseudoJet.hh.
{ reset(some_four_vector[0], some_four_vector[1], some_four_vector[2], some_four_vector[3]); }
void fastjet::PseudoJet::reset | ( | double | px, |
double | py, | ||
double | pz, | ||
double | E | ||
) | [inline] |
reset the 4-momentum according to the supplied components and put the user and history indices back to their default values
Definition at line 327 of file PseudoJet.hh.
{ _px = px; _py = py; _pz = pz; _E = E; _finish_init(); _reset_indices(); }
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 (which does not know about indices...)
Definition at line 193 of file PseudoJet.hh.
{ (*this) = psjet; }
void fastjet::PseudoJet::set_cluster_hist_index | ( | const int | index | ) | [inline] |
set the cluster_hist_index, intended to be used by clustering routines.
Definition at line 137 of file PseudoJet.hh.
Referenced by fastjet::ClusterSequence::_do_ij_recombination_step(), fastjet::ClusterSequence::plugin_record_ij_recombination(), and fastjet::ClusterSequenceAreaBase::subtracted_jet().
{_cluster_hist_index = index;}
void fastjet::PseudoJet::set_cluster_sequence_history_index | ( | const int | index | ) | [inline] |
alternative name for set_cluster_hist_index(...) [perhaps more meaningful]
Definition at line 144 of file PseudoJet.hh.
{ set_cluster_hist_index(index);}
void fastjet::PseudoJet::set_user_index | ( | const int | index | ) | [inline] |
set the user_index, intended to allow the user to "add" information
Definition at line 151 of file PseudoJet.hh.
Referenced by main(), fastjet::JetDefinition::DefaultRecombiner::preprocess(), fastjet::JetDefinition::DefaultRecombiner::recombine(), FlavourRecombiner::recombine(), fastjet::SISConeSphericalPlugin::run_clustering(), fastjet::SISConePlugin::run_clustering(), and fastjet::ClusterSequenceAreaBase::subtracted_jet().
{_user_index = index;}
double fastjet::PseudoJet::squared_distance | ( | const PseudoJet & | other | ) | const [inline] |
returns squared cylinder (rap-phi) distance between this jet and another
Definition at line 164 of file PseudoJet.hh.
{ return plain_distance(other);}
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 261 of file PseudoJet.cc.
References E(), m(), px(), py(), and pz().
{ 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; }
int fastjet::PseudoJet::user_index | ( | ) | const [inline] |
return the user_index, intended to allow the user to "add" information
Definition at line 149 of file PseudoJet.hh.
Referenced by operator<<(), fastjet::JetDefinition::DefaultRecombiner::preprocess(), FlavourRecombiner::recombine(), and fastjet::ClusterSequenceAreaBase::subtracted_jet().
{return _user_index;}
int fastjet::PseudoJet::_cluster_hist_index [private] |
Definition at line 209 of file PseudoJet.hh.
double fastjet::PseudoJet::_E [private] |
Definition at line 207 of file PseudoJet.hh.
Referenced by operator+=(), and operator-=().
double fastjet::PseudoJet::_kt2 [private] |
Definition at line 208 of file PseudoJet.hh.
Referenced by kt_distance().
double fastjet::PseudoJet::_phi [private] |
Definition at line 208 of file PseudoJet.hh.
Referenced by kt_distance(), and plain_distance().
double fastjet::PseudoJet::_px [private] |
Definition at line 207 of file PseudoJet.hh.
Referenced by operator+=(), and operator-=().
double fastjet::PseudoJet::_py [private] |
Definition at line 207 of file PseudoJet.hh.
Referenced by operator+=(), and operator-=().
double fastjet::PseudoJet::_pz [private] |
Definition at line 207 of file PseudoJet.hh.
Referenced by operator+=(), and operator-=().
double fastjet::PseudoJet::_rap [private] |
Definition at line 208 of file PseudoJet.hh.
Referenced by kt_distance(), and plain_distance().
int fastjet::PseudoJet::_user_index [private] |
Definition at line 209 of file PseudoJet.hh.