32 #ifndef __FASTJET_PSEUDOJET_HH__ 
   33 #define __FASTJET_PSEUDOJET_HH__ 
   40 #include "fastjet/internal/numconsts.hh" 
   41 #include "fastjet/internal/IsBase.hh" 
   42 #include "fastjet/SharedPtr.hh" 
   43 #include "fastjet/Error.hh" 
   44 #include "fastjet/PseudoJetStructureBase.hh" 
   46 FASTJET_BEGIN_NAMESPACE      
 
   56 const double pseudojet_invalid_rap = -1e200;
 
   79   PseudoJet() : _px(0), _py(0), _pz(0), _E(0) {_finish_init(); _reset_indices();}
 
   82   PseudoJet(
const double px, 
const double py, 
const double pz, 
const double E);
 
   86   template <
class L> 
PseudoJet(
const L & some_four_vector);
 
  103   inline double E()
   const {
return _E;}
 
  104   inline double e()
   const {
return _E;} 
 
  105   inline double px()
  const {
return _px;}
 
  106   inline double py()
  const {
return _py;}
 
  107   inline double pz()
  const {
return _pz;}
 
  110   inline double phi()
 const {
return phi_02pi();}
 
  114     _ensure_valid_rap_phi();
 
  115     return _phi > pi ? _phi-twopi : _phi;}
 
  119     _ensure_valid_rap_phi();
 
  125   inline double rap()
 const {
 
  126     _ensure_valid_rap_phi();
 
  135   double pseudorapidity() 
const;
 
  136   double eta()
 const {
return pseudorapidity();}
 
  139   inline double pt2()
 const {
return _kt2;}
 
  141   inline double  pt()
 const {
return sqrt(_kt2);} 
 
  143   inline double perp2()
 const {
return _kt2;}  
 
  145   inline double  perp()
 const {
return sqrt(_kt2);}    
 
  147   inline double kt2()
 const {
return _kt2;} 
 
  150   inline double  m2()
 const {
return (_E+_pz)*(_E-_pz)-_kt2;}    
 
  153   inline double  m() 
const;    
 
  156   inline double mperp2()
 const {
return (_E+_pz)*(_E-_pz);}
 
  158   inline double mperp()
 const {
return sqrt(std::abs(mperp2()));}
 
  160   inline double mt2()
 const {
return (_E+_pz)*(_E-_pz);}
 
  162   inline double mt()
 const {
return sqrt(std::abs(mperp2()));}
 
  165   inline double modp2()
 const {
return _kt2+_pz*_pz;}
 
  167   inline double modp()
 const {
return sqrt(_kt2+_pz*_pz);}
 
  170   inline double Et()
 const {
return (_kt2==0) ? 0.0 : _E/sqrt(1.0+_pz*_pz/_kt2);}
 
  172   inline double Et2()
 const {
return (_kt2==0) ? 0.0 : _E*_E/(1.0+_pz*_pz/_kt2);}
 
  177     return std::min(1.0, std::max(-1.0, _pz/sqrt(modp2())));
 
  183   double operator () (
int i) 
const ; 
 
  185   inline double operator [] (
int i)
 const { 
return (*
this)(i); }; 
 
  190   double kt_distance(
const PseudoJet & other) 
const;
 
  193   double plain_distance(
const PseudoJet & other) 
const;
 
  197     return plain_distance(other);}
 
  202     return sqrt(squared_distance(other));
 
  207   double delta_phi_to(
const PseudoJet & other) 
const;
 
  219   std::valarray<double> four_mom() 
const;
 
  224   enum { X=0, Y=1, Z=2, T=3, NUM_COORDINATES=4, SIZE=NUM_COORDINATES };
 
  233   PseudoJet & boost(
const PseudoJet & prest);
 
  236   PseudoJet & unboost(
const PseudoJet & prest);
 
  238   PseudoJet & operator*=(
double);
 
  239   PseudoJet & operator/=(
double);
 
  240   PseudoJet & operator+=(
const PseudoJet &);
 
  241   PseudoJet & operator-=(
const PseudoJet &);
 
  245   inline void reset(
double px, 
double py, 
double pz, 
double E);
 
  262   template <
class L> 
inline void reset(
const L & some_four_vector) {
 
  274     const PseudoJet * pj = fastjet::cast_if_derived<const PseudoJet>(&some_four_vector);
 
  279       reset(some_four_vector[0], some_four_vector[1],
 
  280             some_four_vector[2], some_four_vector[3]);
 
  288   inline void reset_PtYPhiM(
double pt_in, 
double y_in, 
double phi_in, 
double m_in=0.0) {
 
  289     reset_momentum_PtYPhiM(pt_in, y_in, phi_in, m_in);
 
  296   inline void reset_momentum(
double px, 
double py, 
double pz, 
double E);
 
  301   inline void reset_momentum(
const PseudoJet & pj);
 
  305   void reset_momentum_PtYPhiM(
double pt, 
double y, 
double phi, 
double m=0.0);
 
  311     reset_momentum(some_four_vector[0], some_four_vector[1],
 
  312                    some_four_vector[2], some_four_vector[3]);
 
  327   void set_cached_rap_phi(
double rap, 
double phi);
 
  400     _user_info.reset(user_info_in);
 
  434     return dynamic_cast<const L &>(* _user_info.get());
 
  439     return _user_info.get();
 
  446     return _user_info.get() && dynamic_cast<const L *>(_user_info.get());
 
  454     return _user_info.get();
 
  492   std::string description() 
const;
 
  506   bool has_associated_cluster_sequence() 
const;
 
  512   bool has_valid_cluster_sequence() 
const;
 
  520   const ClusterSequence* associated_cs()
 const {
return associated_cluster_sequence();}
 
  525     return validated_cs();
 
  534     return validated_csab();
 
  557   bool has_structure() 
const;
 
  591   template<
typename StructureType>
 
  592   const StructureType & structure() 
const;
 
  597   template<
typename TransformerType>
 
  598   bool has_structure_of() 
const;
 
  605   template<
typename TransformerType>
 
  606   const typename TransformerType::StructureType & structure_of() 
const;
 
  625   virtual bool has_partner(
PseudoJet &partner) 
const;
 
  633   virtual bool has_child(
PseudoJet &child) 
const;
 
  648   virtual bool contains(
const PseudoJet &constituent) 
const;
 
  655   virtual bool is_inside(
const PseudoJet &jet) 
const;
 
  659   virtual bool has_constituents() 
const;
 
  665   virtual std::vector<PseudoJet> constituents() 
const;
 
  669   virtual bool has_exclusive_subjets() 
const;
 
  682   std::vector<PseudoJet> exclusive_subjets (
const double dcut) 
const;
 
  690   int n_exclusive_subjets(
const double dcut) 
const;
 
  700   std::vector<PseudoJet> exclusive_subjets (
int nsub) 
const;
 
  710   std::vector<PseudoJet> exclusive_subjets_up_to (
int nsub) 
const;
 
  719   double exclusive_subdmerge(
int nsub) 
const;
 
  729   double exclusive_subdmerge_max(
int nsub) 
const;
 
  739   virtual bool has_pieces() 
const;
 
  745   virtual std::vector<PseudoJet> pieces() 
const;
 
  754   virtual bool has_area() 
const;
 
  758   virtual double area() 
const;
 
  763   virtual double area_error() 
const;
 
  771   virtual bool is_pure_ghost() 
const;
 
  790     return cluster_hist_index();}
 
  794     set_cluster_hist_index(index);}
 
  806   double _px,_py,_pz,_E;
 
  807   mutable double _phi, _rap;
 
  809   int    _cluster_hist_index, _user_index;
 
  814   void _reset_indices();
 
  818   inline void _ensure_valid_rap_phi()
 const {
 
  819     if (_phi == pseudojet_invalid_phi) _set_rap_phi();
 
  823   void _set_rap_phi() 
const;
 
  826   friend PseudoJet 
operator*(
double, 
const PseudoJet &);
 
  833 PseudoJet operator+(
const PseudoJet &, 
const PseudoJet &);
 
  834 PseudoJet operator-(
const PseudoJet &, 
const PseudoJet &);
 
  835 PseudoJet 
operator*(
double, 
const PseudoJet &);
 
  836 PseudoJet 
operator*(
const PseudoJet &, 
double);
 
  837 PseudoJet operator/(
const PseudoJet &, 
double);
 
  842 bool operator==(
const PseudoJet &, 
const PseudoJet &);
 
  849 bool operator==(
const PseudoJet & jet, 
const double val);
 
  850 inline bool operator==(
const double val, 
const PseudoJet & jet) {
return jet == val;}
 
  855 inline bool operator!=( 
const double val, 
const PseudoJet & a) {
return !(a==val);}
 
  859   return a.E()*b.E() - a.px()*b.px() - a.py()*b.py() - a.pz()*b.pz();
 
  864   double dot_3d = a.px()*b.px() + a.py()*b.py() + a.pz()*b.pz();
 
  865   return std::min(1.0, std::max(-1.0, dot_3d/sqrt(a.
modp2()*b.
modp2())));
 
  874 bool have_same_momentum(
const PseudoJet &, 
const PseudoJet &);
 
  878 PseudoJet PtYPhiM(
double pt, 
double y, 
double phi, 
double m = 0.0);
 
  884 std::vector<PseudoJet> sorted_by_pt(
const std::vector<PseudoJet> & jets);
 
  887 std::vector<PseudoJet> sorted_by_rapidity(
const std::vector<PseudoJet> & jets);
 
  890 std::vector<PseudoJet> sorted_by_E(
const std::vector<PseudoJet> & jets);
 
  893 std::vector<PseudoJet> sorted_by_pz(
const std::vector<PseudoJet> & jets);
 
  900 void sort_indices(std::vector<int> & indices, 
 
  901                   const std::vector<double> & values);
 
  908                                               const std::vector<double> & values) {
 
  910   if (objects.size() != values.size()){
 
  911     throw Error(
"fastjet::objects_sorted_by_values(...): the size of the 'objects' vector must match the size of the 'values' vector");
 
  915   std::vector<int> indices(values.size());
 
  916   for (
size_t i = 0; i < indices.size(); i++) {indices[i] = i;}
 
  919   sort_indices(indices, values);
 
  922   std::vector<T> objects_sorted(objects.size());
 
  925   for (
size_t i = 0; i < indices.size(); i++) {
 
  926     objects_sorted[i] = objects[indices[i]];
 
  929   return objects_sorted;
 
  937 class IndexedSortHelper {
 
  939   inline IndexedSortHelper (
const std::vector<double> * reference_values) {
 
  940     _ref_values = reference_values;
 
  942   inline int operator() (
const int i1, 
const int i2)
 const {
 
  943     return  (*_ref_values)[i1] < (*_ref_values)[i2];
 
  946   const std::vector<double> * _ref_values;
 
  955 template <
class L> 
inline  PseudoJet::PseudoJet(
const L & some_four_vector) {
 
  956   reset(some_four_vector);
 
  961 inline void PseudoJet::_reset_indices() { 
 
  962   set_cluster_hist_index(-1);
 
  970 inline double PseudoJet::m()
 const {
 
  972   return mm < 0.0 ? -std::sqrt(-mm) : std::sqrt(mm);
 
  976 inline void PseudoJet::reset(
double px_in, 
double py_in, 
double pz_in, 
double E_in) {
 
  985 inline void PseudoJet::reset_momentum(
double px_in, 
double py_in, 
double pz_in, 
double E_in) {
 
  993 inline void PseudoJet::reset_momentum(
const PseudoJet & pj) {
 
 1012 template<
typename StructureType>
 
 1013 const StructureType & PseudoJet::structure()
 const{
 
 1014   return dynamic_cast<const StructureType &>(* validated_structure_ptr());
 
 1020 template<
typename TransformerType>
 
 1021 bool PseudoJet::has_structure_of()
 const{
 
 1022   if (!_structure) 
return false;
 
 1024   return dynamic_cast<const typename TransformerType::StructureType *>(_structure.
get()) != 0;
 
 1030 template<
typename TransformerType>
 
 1031 const typename TransformerType::StructureType & PseudoJet::structure_of()
 const{
 
 1033     throw Error(
"Trying to access the structure of a PseudoJet without an associated structure");
 
 1035   return dynamic_cast<const typename TransformerType::StructureType &>(*_structure);
 
 1052 PseudoJet join(
const std::vector<PseudoJet> & pieces);
 
 1068 FASTJET_END_NAMESPACE
 
 1070 #endif // __FASTJET_PSEUDOJET_HH__