PseudoJet Class Reference

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

#include <PseudoJet.hh>

Inheritance diagram for PseudoJet:
Inheritance graph
[legend]

List of all members.

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
PseudoJetboost (const PseudoJet &prest)
 transform this jet (given in the rest frame of prest) into a jet in the lab frame [NOT FULLY TESTED]
PseudoJetunboost (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(.
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,.
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

Detailed Description

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

Definition at line 52 of file PseudoJet.hh.


Member Enumeration Documentation

anonymous enum
Enumerator:
X 
Y 
Z 
T 
NUM_COORDINATES 
SIZE 

Definition at line 123 of file PseudoJet.hh.

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


Constructor & Destructor Documentation

PseudoJet::PseudoJet (  )  [inline]

Definition at line 55 of file PseudoJet.hh.

Referenced by PseudoJet().

00055 {};

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.

References _E, _finish_init(), _px, _py, _pz, and _reset_indices().

00047                                                                                       {
00048   
00049   _E  = E ;
00050   _px = px;
00051   _py = py;
00052   _pz = pz;
00053 
00054   this->_finish_init();
00055 
00056   // some default values for the history and user indices
00057   _reset_indices();
00058 
00059 }

template<class L >
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.

References _E, _finish_init(), _px, _py, _pz, and _reset_indices().

00286                                                                           {
00287 
00288   _px = some_four_vector[0];
00289   _py = some_four_vector[1];
00290   _pz = some_four_vector[2];
00291   _E  = some_four_vector[3];
00292   _finish_init();
00293   // some default values for these two indices
00294   _reset_indices();
00295 }

template<>
PseudoJet::PseudoJet ( const siscone::Cmomentum &  four_vector  )  [inline]

shortcut for converting siscone Cmomentum into PseudoJet

Definition at line 19 of file SISConePlugin.cc.

References PseudoJet().

00019                                                                     {
00020   (*this) = PseudoJet(four_vector.px,four_vector.py,four_vector.pz,
00021                       four_vector.E);
00022 }

template<>
PseudoJet::PseudoJet ( const siscone_spherical::CSphmomentum &  four_vector  )  [inline]

shortcut for converting siscone CSphmomentum into PseudoJet

Definition at line 19 of file SISConeSphericalPlugin.cc.

References PseudoJet().

00019                                                                                  {
00020   (*this) = PseudoJet(four_vector.px,four_vector.py,four_vector.pz,
00021                       four_vector.E);
00022 }


Member Function Documentation

void 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 _E, _kt2, _phi, _pz, _rap, E(), m2(), MaxRap, px(), py(), pz(), and twopi.

Referenced by boost(), operator+=(), operator-=(), PseudoJet(), reset(), and unboost().

00064                               {
00065   _kt2 = this->px()*this->px() + this->py()*this->py();
00066   if (_kt2 == 0.0) {
00067     _phi = 0.0; } 
00068   else {
00069     _phi = atan2(this->py(),this->px());
00070   }
00071   if (_phi < 0.0) {_phi += twopi;}
00072   if (_phi >= twopi) {_phi -= twopi;} // can happen if phi=-|eps<1e-15|?
00073   if (this->E() == abs(this->pz()) && _kt2 == 0) {
00074     // Point has infinite rapidity -- convert that into a very large
00075     // number, but in such a way that different 0-pt momenta will have
00076     // different rapidities (so as to lift the degeneracy between
00077     // them) [this can be relevant at parton-level]
00078     double MaxRapHere = MaxRap + abs(this->pz());
00079     if (this->pz() >= 0.0) {_rap = MaxRapHere;} else {_rap = -MaxRapHere;}
00080   } else {
00081     // get the rapidity in a way that's modestly insensitive to roundoff
00082     // error when things pz,E are large (actually the best we can do without
00083     // explicit knowledge of mass)
00084     double effective_m2 = max(0.0,m2()); // force non tachyonic mass
00085     double E_plus_pz    = _E + abs(_pz); // the safer of p+, p-
00086     // p+/p- = (p+ p-) / (p-)^2 = (kt^2+m^2)/(p-)^2
00087     _rap = 0.5*log((_kt2 + effective_m2)/(E_plus_pz*E_plus_pz));
00088     if (_pz > 0) {_rap = - _rap;}
00089   }
00090 
00092   //if (this->E() != abs(this->pz())) {
00093   //  _rap = 0.5*log((this->E() + this->pz())/(this->E() - this->pz()));
00094   //    } else {
00095   //  // Overlapping points can give problems. Let's lift the degeneracy
00096   //  // in case of multiple 0-pT points (can be found at parton-level)
00097   //  double MaxRapHere = MaxRap + abs(this->pz());
00098   //  if (this->pz() >= 0.0) {_rap = MaxRapHere;} else {_rap = -MaxRapHere;}
00099   //}
00100 }

void PseudoJet::_reset_indices (  )  [inline, private]

set the indices to default values

Definition at line 299 of file PseudoJet.hh.

References set_cluster_hist_index(), and set_user_index().

Referenced by PseudoJet(), and reset().

00299                                       { 
00300   set_cluster_hist_index(-1);
00301   set_user_index(-1);
00302 }

double PseudoJet::beam_distance (  )  const [inline]

returns distance between this jet and the beam

Definition at line 177 of file PseudoJet.hh.

References _kt2.

00177 {return _kt2;}

PseudoJet & 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 lab) into a jet in the rest frame of prest

Definition at line 234 of file PseudoJet.cc.

References _E, _finish_init(), _px, _py, _pz, E(), m(), px(), py(), and pz().

00234                                                     {
00235   
00236   if (prest.px() == 0.0 && prest.py() == 0.0 && prest.pz() == 0.0) 
00237     return *this;
00238 
00239   double m = prest.m();
00240   assert(m != 0);
00241 
00242   double pf4  = (  px()*prest.px() + py()*prest.py()
00243                  + pz()*prest.pz() + E()*prest.E() )/m;
00244   double fn   = (pf4 + E()) / (prest.E() + m);
00245   _px +=  fn*prest.px();
00246   _py +=  fn*prest.py();
00247   _pz +=  fn*prest.pz();
00248   _E = pf4;
00249 
00250   _finish_init(); // we need to recalculate phi,rap,kt2
00251   return *this;
00252 }

int PseudoJet::cluster_hist_index (  )  const [inline]
int PseudoJet::cluster_sequence_history_index (  )  const [inline]

alternative name for cluster_hist_index() [perhaps more meaningful]

Definition at line 140 of file PseudoJet.hh.

References cluster_hist_index().

00140                                                     {
00141     return cluster_hist_index();}

double 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 _phi, fastjet::pi, and twopi.

00325                                                             {
00326   double dphi = abs(other._phi - _phi);
00327   if (dphi >  pi) dphi -= twopi;
00328   if (dphi < -pi) dphi += twopi;
00329   return dphi;
00330 }

double PseudoJet::e (  )  const [inline]

Definition at line 67 of file PseudoJet.hh.

References _E.

Referenced by operator()().

00067 {return _E;} // like CLHEP

double PseudoJet::E (  )  const [inline]
double PseudoJet::Et (  )  const [inline]

return the transverse energy

Definition at line 112 of file PseudoJet.hh.

References _E, _kt2, and _pz.

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

00112 {return (_kt2==0) ? 0.0 : _E/sqrt(1.0+_pz*_pz/_kt2);}

double PseudoJet::Et2 (  )  const [inline]

return the transverse energy squared

Definition at line 114 of file PseudoJet.hh.

References _E, _kt2, and _pz.

00114 {return (_kt2==0) ? 0.0 : _E*_E/(1.0+_pz*_pz/_kt2);}

double PseudoJet::eta (  )  const [inline]

Definition at line 92 of file PseudoJet.hh.

References pseudorapidity().

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

00092 {return pseudorapidity();}

valarray< double > PseudoJet::four_mom (  )  const

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

Definition at line 105 of file PseudoJet.cc.

References _E, _px, _py, and _pz.

00105                                            {
00106   valarray<double> mom(4);
00107   mom[0] = _px;
00108   mom[1] = _py;
00109   mom[2] = _pz;
00110   mom[3] = _E ;
00111   return mom;
00112 }

double PseudoJet::kt2 (  )  const [inline]

returns the squared transverse momentum

Definition at line 95 of file PseudoJet.hh.

References _kt2.

Referenced by ClusterSequence::jet_scale_for_algorithm().

00095 {return _kt2;}

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

00302                                                            {
00303   //double distance = min(this->kt2(), other.kt2());
00304   double distance = min(_kt2, other._kt2);
00305   double dphi = abs(_phi - other._phi);
00306   if (dphi > pi) {dphi = twopi - dphi;}
00307   double drap = _rap - other._rap;
00308   distance = distance * (dphi*dphi + drap*drap);
00309   return distance;
00310 }

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

References m2().

Referenced by boost(), and unboost().

00321                                  {
00322   double mm = m2();
00323   return mm < 0.0 ? -std::sqrt(-mm) : std::sqrt(mm);
00324 }

double PseudoJet::m2 (  )  const [inline]

returns the squared invariant mass // like CLHEP

Definition at line 101 of file PseudoJet.hh.

References _E, _kt2, and _pz.

Referenced by _finish_init(), and m().

00101 {return (_E+_pz)*(_E-_pz)-_kt2;}    

double PseudoJet::modp2 (  )  const [inline]

return px^2+py^2+pz^2

Definition at line 110 of file PseudoJet.hh.

References _kt2, and _pz.

Referenced by fastjet::JadeBriefJet::init(), and fastjet::EECamBriefJet::init().

00110 {return _kt2+_pz*_pz;}

double PseudoJet::mperp (  )  const [inline]

returns the transverse mass = sqrt(kt^2+m^2)

Definition at line 105 of file PseudoJet.hh.

References mperp2().

00105 {return sqrt(std::abs(mperp2()));}

double PseudoJet::mperp2 (  )  const [inline]

returns the squared transverse mass = kt^2+m^2

Definition at line 103 of file PseudoJet.hh.

References _E, and _pz.

Referenced by mperp().

00103 {return (_E+_pz)*(_E-_pz);}

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

References e(), px(), py(), pz(), T, X, Y, and Z.

00117                                           {
00118   switch(i) {
00119   case X:
00120     return px();
00121   case Y:
00122     return py();
00123   case Z:
00124     return pz();
00125   case T:
00126     return e();
00127   default:
00128     ostringstream err;
00129     err << "PseudoJet subscripting: bad index (" << i << ")";
00130     throw Error(err.str());
00131   }
00132   return 0.;
00133 }  

void PseudoJet::operator*= ( double  coeff  ) 

multiply the jet's momentum by the coefficient

Definition at line 190 of file PseudoJet.cc.

References _E, _kt2, _px, _py, and _pz.

00190                                        {
00191   _px *= coeff;
00192   _py *= coeff;
00193   _pz *= coeff;
00194   _E  *= coeff;
00195   _kt2*= coeff*coeff;
00196   // phi and rap are unchanged
00197 }

void PseudoJet::operator+= ( const PseudoJet other_jet  ) 

add the other jet's momentum to this jet

Definition at line 208 of file PseudoJet.cc.

References _E, _finish_init(), _px, _py, and _pz.

00208                                                       {
00209   _px += other_jet._px;
00210   _py += other_jet._py;
00211   _pz += other_jet._pz;
00212   _E  += other_jet._E ;
00213   _finish_init(); // we need to recalculate phi,rap,kt2
00214 }

void PseudoJet::operator-= ( const PseudoJet other_jet  ) 

subtract the other jet's momentum from this jet

Definition at line 219 of file PseudoJet.cc.

References _E, _finish_init(), _px, _py, and _pz.

00219                                                       {
00220   _px -= other_jet._px;
00221   _py -= other_jet._py;
00222   _pz -= other_jet._pz;
00223   _E  -= other_jet._E ;
00224   _finish_init(); // we need to recalculate phi,rap,kt2
00225 }

void PseudoJet::operator/= ( double  coeff  ) 

divide the jet's momentum by the coefficient

Definition at line 201 of file PseudoJet.cc.

00201                                        {
00202   (*this) *= 1.0/coeff;
00203 }

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

00119 { return (*this)(i); }; // this too

double PseudoJet::perp (  )  const [inline]
double PseudoJet::perp2 (  )  const [inline]
double PseudoJet::phi (  )  const [inline]
double PseudoJet::phi_02pi (  )  const [inline]

returns phi in the range 0..2pi

Definition at line 80 of file PseudoJet.hh.

References _phi.

Referenced by phi().

00080 {return _phi;}

double PseudoJet::phi_std (  )  const [inline]

returns phi in the range -pi..pi

Definition at line 76 of file PseudoJet.hh.

References _phi, fastjet::pi, and twopi.

00076                                  {
00077     return _phi > pi ? _phi-twopi : _phi;}

double 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(), and squared_distance().

00315                                                               {
00316   double dphi = abs(_phi - other._phi);
00317   if (dphi > pi) {dphi = twopi - dphi;}
00318   double drap = _rap - other._rap;
00319   return (dphi*dphi + drap*drap);
00320 }

double 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 MaxRap, perp(), fastjet::pi, px(), py(), and pz().

Referenced by eta().

00137                                        {
00138   if (px() == 0.0 && py() ==0.0) return MaxRap;
00139   if (pz() == 0.0) return 0.0;
00140 
00141   double theta = atan(perp()/pz());
00142   if (theta < 0) theta += pi;
00143   return -log(tan(theta/2));
00144 }

double PseudoJet::px (  )  const [inline]
double PseudoJet::py (  )  const [inline]
double PseudoJet::pz (  )  const [inline]
double PseudoJet::rap (  )  const [inline]

returns the rapidity or some large value when the rapidity is infinite

Definition at line 84 of file PseudoJet.hh.

References _rap.

Referenced by RangeDefinition::is_in_range(), JetDefinition::DefaultRecombiner::recombine(), and RangeDefinition::set_position().

00084 {return _rap;}

double PseudoJet::rapidity (  )  const [inline]

the same as rap()

Definition at line 87 of file PseudoJet.hh.

References _rap.

00087 {return _rap;} // like CLHEP

template<class L >
void 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.

References reset().

00200                                                                    {
00201     reset(some_four_vector[0], some_four_vector[1],
00202           some_four_vector[2], some_four_vector[3]);
00203   }

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

00193                                              {
00194     (*this) = psjet;
00195   }

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

References _E, _finish_init(), _px, _py, _pz, and _reset_indices().

Referenced by reset().

00327                                                                       {
00328   _px = px;
00329   _py = py;
00330   _pz = pz;
00331   _E  = E;
00332   _finish_init();
00333   _reset_indices();
00334 }

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

References _cluster_hist_index.

Referenced by _reset_indices(), set_cluster_sequence_history_index(), and ClusterSequenceAreaBase::subtracted_jet().

00137 {_cluster_hist_index = index;}

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

References set_cluster_hist_index().

00144                                                                   {
00145     set_cluster_hist_index(index);}

void PseudoJet::set_user_index ( const int  index  )  [inline]
double 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.

References plain_distance().

00164                                                                 {
00165     return plain_distance(other);}

PseudoJet & 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 the rest frame of prest) into a jet in the lab frame;

Definition at line 261 of file PseudoJet.cc.

References _E, _finish_init(), _px, _py, _pz, E(), m(), px(), py(), and pz().

00261                                                       {
00262   
00263   if (prest.px() == 0.0 && prest.py() == 0.0 && prest.pz() == 0.0) 
00264     return *this;
00265 
00266   double m = prest.m();
00267   assert(m != 0);
00268 
00269   double pf4  = ( -px()*prest.px() - py()*prest.py()
00270                  - pz()*prest.pz() + E()*prest.E() )/m;
00271   double fn   = (pf4 + E()) / (prest.E() + m);
00272   _px -=  fn*prest.px();
00273   _py -=  fn*prest.py();
00274   _pz -=  fn*prest.pz();
00275   _E = pf4;
00276 
00277   _finish_init(); // we need to recalculate phi,rap,kt2
00278   return *this;
00279 }

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

References _user_index.

Referenced by JetDefinition::DefaultRecombiner::preprocess(), and ClusterSequenceAreaBase::subtracted_jet().

00149 {return _user_index;}


Member Data Documentation

Definition at line 209 of file PseudoJet.hh.

Referenced by cluster_hist_index(), and set_cluster_hist_index().

double PseudoJet::_E [private]
double PseudoJet::_kt2 [private]
double PseudoJet::_phi [private]

Definition at line 208 of file PseudoJet.hh.

Referenced by _finish_init(), delta_phi_to(), kt_distance(), phi_02pi(), phi_std(), and plain_distance().

double PseudoJet::_px [private]

Definition at line 207 of file PseudoJet.hh.

Referenced by boost(), four_mom(), operator*=(), operator+=(), operator-=(), PseudoJet(), px(), reset(), and unboost().

double PseudoJet::_py [private]

Definition at line 207 of file PseudoJet.hh.

Referenced by boost(), four_mom(), operator*=(), operator+=(), operator-=(), PseudoJet(), py(), reset(), and unboost().

double PseudoJet::_pz [private]
double PseudoJet::_rap [private]

Definition at line 208 of file PseudoJet.hh.

Referenced by _finish_init(), kt_distance(), plain_distance(), rap(), and rapidity().

int PseudoJet::_user_index [private]

Definition at line 209 of file PseudoJet.hh.

Referenced by set_user_index(), and user_index().


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

Generated on 26 Feb 2010 for fastjet by  doxygen 1.6.1