PseudoJet.hh

00001 //STARTHEADER
00002 // $Id: PseudoJet.hh 1987 2011-03-10 10:01:56Z salam $
00003 //
00004 // Copyright (c) 2005-2011, Matteo Cacciari, Gavin Salam and Gregory Soyez
00005 //
00006 //----------------------------------------------------------------------
00007 // This file is part of FastJet.
00008 //
00009 //  FastJet is free software; you can redistribute it and/or modify
00010 //  it under the terms of the GNU General Public License as published by
00011 //  the Free Software Foundation; either version 2 of the License, or
00012 //  (at your option) any later version.
00013 //
00014 //  The algorithms that underlie FastJet have required considerable
00015 //  development and are described in hep-ph/0512210. If you use
00016 //  FastJet as part of work towards a scientific publication, please
00017 //  include a citation to the FastJet paper.
00018 //
00019 //  FastJet is distributed in the hope that it will be useful,
00020 //  but WITHOUT ANY WARRANTY; without even the implied warranty of
00021 //  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
00022 //  GNU General Public License for more details.
00023 //
00024 //  You should have received a copy of the GNU General Public License
00025 //  along with FastJet; if not, write to the Free Software
00026 //  Foundation, Inc.:
00027 //      59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
00028 //----------------------------------------------------------------------
00029 //ENDHEADER
00030 
00031 
00032 #ifndef __FASTJET_PSEUDOJET_HH__
00033 #define __FASTJET_PSEUDOJET_HH__
00034 
00035 #include<valarray>
00036 #include<vector>
00037 #include<cassert>
00038 #include<cmath>
00039 #include<iostream>
00040 #include "fastjet/internal/numconsts.hh"
00041 #include "fastjet/internal/IsBase.hh"
00042 #include "fastjet/internal/DerivedPseudoJetHelper.hh"
00043 #include "fastjet/SharedPtr.hh"
00044 #include "fastjet/Error.hh"
00045 #include "fastjet/PseudoJetStructureBase.hh"
00046 
00047 FASTJET_BEGIN_NAMESPACE      // defined in fastjet/internal/base.hh
00048 
00049 //using namespace std;
00050 
00051 /// Used to protect against parton-level events where pt can be zero
00052 /// for some partons, giving rapidity=infinity. KtJet fails in those cases.
00053 const double MaxRap = 1e5;
00054 
00055 /// default value for phi, meaning it (and rapidity) have yet to be calculated) 
00056 const double pseudojet_invalid_phi = -100.0;
00057 
00058 // forward definition
00059 class ClusterSequenceAreaBase;
00060 
00061 /// @ingroup basic_classes
00062 /// \class PseudoJet
00063 /// Class to contain pseudojets, including minimal information of use to
00064 /// to jet-clustering routines.
00065 class PseudoJet {
00066 
00067  public:
00068   //----------------------------------------------------------------------
00069   /// @name Constructors and destructor
00070   //\{
00071   /// default constructor leaves PseudoJet unusable
00072   PseudoJet() {};
00073   /// construct a pseudojet from explicit components
00074   PseudoJet(const double px, const double py, const double pz, const double E);
00075   /// constructor from any object that has px,py,pz,E = some_four_vector[0--3],
00076   template <class L> PseudoJet(const L & some_four_vector) ;
00077 
00078   /// default (virtual) destructor
00079   virtual ~PseudoJet(){};
00080   //\} ---- end of constructors and destructors --------------------------
00081 
00082   //----------------------------------------------------------------------
00083   /// @name Kinematic access functions
00084   //\{
00085   //----------------------------------------------------------------------
00086   inline double E()   const {return _E;}
00087   inline double e()   const {return _E;} // like CLHEP
00088   inline double px()  const {return _px;}
00089   inline double py()  const {return _py;}
00090   inline double pz()  const {return _pz;}
00091 
00092   /// returns phi (in the range 0..2pi)
00093   inline double phi() const {return phi_02pi();}
00094 
00095   /// returns phi in the range -pi..pi
00096   inline double phi_std()  const {
00097     _ensure_valid_rap_phi();
00098     return _phi > pi ? _phi-twopi : _phi;}
00099 
00100   /// returns phi in the range 0..2pi
00101   inline double phi_02pi() const {
00102     _ensure_valid_rap_phi();
00103     return _phi;
00104   }
00105 
00106   /// returns the rapidity or some large value when the rapidity
00107   /// is infinite
00108   inline double rap() const {
00109     _ensure_valid_rap_phi();
00110     return _rap;
00111   }
00112 
00113   /// the same as rap()
00114   inline double rapidity() const {return rap();} // like CLHEP
00115 
00116   /// returns the pseudo-rapidity or some large value when the
00117   /// rapidity is infinite
00118   double pseudorapidity() const;
00119   double eta() const {return pseudorapidity();}
00120 
00121   /// returns the squared transverse momentum
00122   inline double kt2() const {return _kt2;}
00123   /// returns the squared transverse momentum
00124   inline double perp2() const {return _kt2;}  // like CLHEP
00125   /// returns the scalar transverse momentum
00126   inline double  perp() const {return sqrt(_kt2);}    // like CLHEP
00127   /// returns the squared invariant mass // like CLHEP
00128   inline double  m2() const {return (_E+_pz)*(_E-_pz)-_kt2;}    
00129   /// returns the squared transverse mass = kt^2+m^2
00130   inline double mperp2() const {return (_E+_pz)*(_E-_pz);}
00131   /// returns the transverse mass = sqrt(kt^2+m^2)
00132   inline double mperp() const {return sqrt(std::abs(mperp2()));}
00133   /// returns the invariant mass 
00134   /// (If m2() is negative then -sqrt(-m2()) is returned, as in CLHEP)
00135   inline double  m() const;    
00136   /// return px^2+py^2+pz^2
00137   inline double modp2() const {return _kt2+_pz*_pz;}
00138   /// return the transverse energy
00139   inline double Et() const {return (_kt2==0) ? 0.0 : _E/sqrt(1.0+_pz*_pz/_kt2);}
00140   /// return the transverse energy squared
00141   inline double Et2() const {return (_kt2==0) ? 0.0 : _E*_E/(1.0+_pz*_pz/_kt2);}
00142 
00143   /// returns component i, where X==0, Y==1, Z==2, E==3
00144   double operator () (int i) const ; 
00145   /// returns component i, where X==0, Y==1, Z==2, E==3
00146   inline double operator [] (int i) const { return (*this)(i); }; // this too
00147 
00148 
00149 
00150   /// returns kt distance (R=1) between this jet and another
00151   double kt_distance(const PseudoJet & other) const;
00152 
00153   /// returns squared cylinder (rap-phi) distance between this jet and another
00154   double plain_distance(const PseudoJet & other) const;
00155   /// returns squared cylinder (rap-phi) distance between this jet and
00156   /// another
00157   inline double squared_distance(const PseudoJet & other) const {
00158     return plain_distance(other);}
00159 
00160   /// returns other.phi() - this.phi(), constrained to be in 
00161   /// range -pi .. pi
00162   double delta_phi_to(const PseudoJet & other) const;
00163 
00164   //// this seemed to compile except if it was used
00165   //friend inline double 
00166   //  kt_distance(const PseudoJet & jet1, const PseudoJet & jet2) { 
00167   //                                      return jet1.kt_distance(jet2);}
00168 
00169   /// returns distance between this jet and the beam
00170   inline double beam_distance() const {return _kt2;}
00171 
00172   /// return a valarray containing the four-momentum (components 0-2
00173   /// are 3-mom, component 3 is energy).
00174   std::valarray<double> four_mom() const;
00175 
00176   //\}  ------- end of kinematic access functions
00177 
00178   // taken from CLHEP
00179   enum { X=0, Y=1, Z=2, T=3, NUM_COORDINATES=4, SIZE=NUM_COORDINATES };
00180 
00181 
00182   //----------------------------------------------------------------------
00183   /// @name Kinematic modification functions
00184   //\{
00185   //----------------------------------------------------------------------
00186   /// transform this jet (given in the rest frame of prest) into a jet
00187   /// in the lab frame [NOT FULLY TESTED]
00188   PseudoJet & boost(const PseudoJet & prest);
00189   /// transform this jet (given in lab) into a jet in the rest
00190   /// frame of prest  [NOT FULLY TESTED]
00191   PseudoJet & unboost(const PseudoJet & prest);
00192 
00193   void operator*=(double);
00194   void operator/=(double);
00195   void operator+=(const PseudoJet &);
00196   void operator-=(const PseudoJet &);
00197 
00198   /// reset the 4-momentum according to the supplied components and
00199   /// put the user and history indices back to their default values
00200   inline void reset(double px, double py, double pz, double E);
00201   
00202   /// reset the PseudoJet to be equal to psjet (including its
00203   /// indices); NB if the argument is derived from a PseudoJet then
00204   /// the "reset" used will be the templated version (which does not
00205   /// know about indices...)
00206   inline void reset(const PseudoJet & psjet) {
00207     (*this) = psjet;
00208   }
00209 
00210   /// reset the 4-momentum according to the supplied generic 4-vector
00211   /// (accessible via indexing, [0]==px,...[3]==E) and put the user
00212   /// and history indices back to their default values.
00213   template <class L> inline void reset(const L & some_four_vector) {
00214     reset(some_four_vector[0], some_four_vector[1],
00215           some_four_vector[2], some_four_vector[3]);
00216   }
00217 
00218   //\} --- end of kin mod functions ------------------------------------
00219 
00220   //----------------------------------------------------------------------
00221   /// @name User index functions
00222   ///
00223   /// To allow the user to set and access an integer index which can
00224   /// be exploited by the user to associate extra information with a
00225   /// particle/jet (for example pdg id, or an indication of a
00226   /// particle's origin within the user's analysis)
00227   //
00228   //\{
00229 
00230   /// return the user_index, 
00231   inline int user_index() const {return _user_index;}
00232   /// set the user_index, intended to allow the user to add simple
00233   /// identifying information to a particle/jet
00234   inline void set_user_index(const int index) {_user_index = index;}
00235 
00236   //\} ----- end of use index functions ---------------------------------
00237 
00238   //----------------------------------------------------------------------
00239   /// @name User information types and functions
00240   ///
00241   /// Allows PseudoJet to carry extra user info (as an object derived from
00242   /// UserInfoBase).
00243   //\{
00244 
00245   /// @ingroup user_info
00246   /// \class UserInfoBase
00247   /// a base class to hold extra user information in a PseudoJet
00248   ///
00249   /// This is a base class to help associate extra user information
00250   /// with a jet. The user should store their information in a class
00251   /// derived from this. This allows information of arbitrary
00252   /// complexity to be easily associated with a PseudoJet (in contrast
00253   /// to the user index). For example, in a Monte Carlo simulation,
00254   /// the user information might include the PDG ID, and the position
00255   /// of the production vertex for the particle.
00256   ///
00257   /// The PseudoJet is able to store a shared pointer to any object
00258   /// derived from UserInfo. The use of a shared pointer frees the
00259   /// user of the need to handle the memory management associated with
00260   /// the information.
00261   ///
00262   /// Having the user information derive from a common base class also
00263   /// facilitates dynamic casting, etc.
00264   ///
00265   class UserInfoBase{
00266   public:
00267     // dummy ctor
00268     UserInfoBase(){};
00269 
00270     // dummy virtual dtor
00271     // makes it polymorphic to allow for dynamic_cast
00272     virtual ~UserInfoBase(){}; 
00273   };
00274 
00275   /// error class to be thrown if accessing user info when it doesn't
00276   /// exist
00277   class InexistentUserInfo : public Error {
00278   public:
00279     InexistentUserInfo();
00280   };
00281 
00282   /// sets the internal shared pointer to the user information.
00283   ///
00284   /// Note that the PseudoJet will now _own_ the pointer, and delete
00285   /// the corresponding object when it (the jet, and any copies of the jet)
00286   /// goes out of scope. 
00287   void set_user_info(UserInfoBase * user_info_in) {
00288     _user_info.reset(user_info_in);
00289   }
00290 
00291   /// returns a reference to the dynamic cast conversion of user_info
00292   /// to type L.
00293   ///
00294   /// Usage: suppose you have previously set the user info with a pointer
00295   /// to an object of type MyInfo, 
00296   ///
00297   ///   class MyInfo: public PseudoJet::UserInfoBase {
00298   ///      MyInfo(int id) : _pdg_id(id);
00299   ///      int pdg_id() const {return _pdg_id;}
00300   ///      int _pdg_id;
00301   ///   };
00302   ///
00303   ///   PseudoJet particle(...);
00304   ///   particle.set_user_info(new MyInfo(its_pdg_id));
00305   ///
00306   /// Then you would access that pdg_id() as
00307   ///
00308   ///   particle.user_info<MyInfo>().pdg_id();
00309   ///
00310   /// It's overkill for just a single integer, but scales easily to
00311   /// more extensive information.
00312   ///
00313   /// Note that user_info() throws an InexistentUserInfo() error if
00314   /// there is no user info; throws a std::bad_cast if the conversion
00315   /// doesn't work
00316   ///
00317   /// If this behaviour does not fit your needs, use instead the the
00318   /// user_info_ptr() or user_info_shared_ptr() member functions.
00319   template<class L>
00320   const L & user_info() const{
00321     if (_user_info.get() == 0) throw InexistentUserInfo();
00322     return dynamic_cast<const L &>(* _user_info.get());
00323   }
00324 
00325   /// retrieve a pointer to the (const) user information
00326   const UserInfoBase * user_info_ptr() const{
00327     if (!_user_info()) return NULL;
00328     return _user_info.get();
00329   }
00330 
00331 
00332   /// retrieve a (const) shared pointer to the user information
00333   const SharedPtr<UserInfoBase> & user_info_shared_ptr() const{
00334     return _user_info;
00335   }
00336 
00337   /// retrieve a (non-const) shared pointer to the user information
00338   SharedPtr<UserInfoBase> & user_info_shared_ptr(){
00339     return _user_info;
00340   }
00341 
00342   // \} --- end of extra info functions ---------------------------------
00343 
00344   //----------------------------------------------------------------------
00345   /// @name Description
00346   ///
00347   /// Since a PseudoJet can have a structure that contains a variety
00348   /// of information, we provide a description that allows one to check
00349   /// exactly what kind of PseudoJet we are dealing with
00350   //
00351   //\{
00352 
00353   /// return a string describing what kind of PseudoJet we are dealing with 
00354   std::string description() const;
00355 
00356   //\} ----- end of description functions ---------------------------------
00357 
00358   //-------------------------------------------------------------
00359   /// @name Access to the associated ClusterSequence object.
00360   ///
00361   /// In addition to having kinematic information, jets may contain a
00362   /// reference to an associated ClusterSequence (this is the case,
00363   /// for example, if the jet has been returned by a ClusterSequence
00364   /// member function).
00365   //\{
00366   //-------------------------------------------------------------
00367   /// returns true if this PseudoJet has an associated (and still
00368   /// valid) ClusterSequence.
00369   bool has_associated_cluster_sequence() const;
00370 
00371   /// get a (const) pointer to the parent ClusterSequence (NULL if
00372   /// inexistent)
00373   const ClusterSequence* associated_cluster_sequence() const;
00374 
00375   /// if the jet has a valid associated cluster sequence then return a
00376   /// pointer to it; otherwise throw an error
00377   const ClusterSequence * validated_cs() const;
00378 
00379   /// if the jet has valid area information then return a pointer to
00380   /// the associated ClusterSequenceAreaBase object; otherwise throw an error
00381   const ClusterSequenceAreaBase * validated_csab() const;
00382 
00383   //\}
00384 
00385   //-------------------------------------------------------------
00386   /// @name Access to the associated PseudoJetStructureBase object.
00387   ///
00388   /// In addition to having kinematic information, jets may contain a
00389   /// reference to an associated ClusterSequence (this is the case,
00390   /// for example, if the jet has been returned by a ClusterSequence
00391   /// member function).
00392   //\{
00393   //-------------------------------------------------------------
00394 
00395   /// set the associated structure
00396   void set_structure_shared_ptr(const SharedPtr<PseudoJetStructureBase> &structure);
00397 
00398   /// return true if there is some strusture associated with this PseudoJet
00399   bool has_structure() const;
00400 
00401   /// return a pointer to the structure (of type
00402   /// PseudoJetStructureBase*) associated wioth this PseudoJet.
00403   ///
00404   /// return NULL if there is no associated structure
00405   const PseudoJetStructureBase* structure_ptr() const;
00406   
00407   /// return a pointer to the structure (of type
00408   /// PseudoJetStructureBase*) associated wioth this PseudoJet.
00409   ///
00410   /// throw an error if there is no associated structure
00411   const PseudoJetStructureBase* validated_structure_ptr() const;
00412   
00413   /// return a reference to the shared pointer to the
00414   /// PseudoJetStructureBase associated wioth this PseudoJet
00415   const SharedPtr<PseudoJetStructureBase> & structure_shared_ptr() const;
00416 
00417   /// returns a reference to the structure casted to the requested
00418   /// structure type
00419   ///
00420   /// If there is no sructure associated, an Error is thrown.
00421   /// If the type is not met, a std::bad_cast error is thrown.
00422   template<typename StructureType>
00423   const StructureType & structure() const;
00424 
00425 
00426 
00427   /// check if the PseudoJet has the properties of the result of a Transformer 
00428   /// (that is, its structure is compatible with a Transformer::StructureType)
00429   /// if there is no structure, false is returned
00430   template<typename TransformerType>
00431   bool has_properties_of() const;
00432 
00433   /// this is a helper to access an structuree created by a Transformer 
00434   /// (that is, of type Transformer::StructureType)
00435   /// NULL is returned if the corresponding type is not met
00436   /// if there is no structure, an error is thrown
00437   template<typename TransformerType>
00438   const typename TransformerType::StructureType * extra_properties() const;
00439 
00440   //\}
00441 
00442   //-------------------------------------------------------------
00443   /// @name Methods for access to information about jet structure
00444   ///
00445   /// These allow access to jet constituents, and other jet
00446   /// subtructure information. They only work if the jet is associated
00447   /// with a ClusterSequence.
00448   //-------------------------------------------------------------
00449   //\{
00450 
00451   /// check if it has been recombined with another PseudoJet in which
00452   /// case, return its partner through the argument. Otherwise,
00453   /// 'partner' is set to 0.
00454   ///
00455   /// an Error is thrown if this PseudoJet has no currently valid
00456   /// associated ClusterSequence
00457   virtual bool has_partner(PseudoJet &partner) const;
00458 
00459   /// check if it has been recombined with another PseudoJet in which
00460   /// case, return its child through the argument. Otherwise, 'child'
00461   /// is set to 0.
00462   /// 
00463   /// an Error is thrown if this PseudoJet has no currently valid
00464   /// associated ClusterSequence
00465   virtual bool has_child(PseudoJet &child) const;
00466 
00467   /// check if it is the product of a recombination, in which case
00468   /// return the 2 parents through the 'parent1' and 'parent2'
00469   /// arguments. Otherwise, set these to 0.
00470   ///
00471   /// an Error is thrown if this PseudoJet has no currently valid
00472   /// associated ClusterSequence
00473   virtual bool has_parents(PseudoJet &parent1, PseudoJet &parent2) const;
00474 
00475   /// check if the current PseudoJet contains the one passed as
00476   /// argument.
00477   ///
00478   /// an Error is thrown if this PseudoJet has no currently valid
00479   /// associated ClusterSequence
00480   virtual bool contains(const PseudoJet &constituent) const;
00481 
00482   /// check if the current PseudoJet is contained the one passed as
00483   /// argument.
00484   ///
00485   /// an Error is thrown if this PseudoJet has no currently valid
00486   /// associated ClusterSequence
00487   virtual bool is_inside(const PseudoJet &jet) const;
00488 
00489 
00490   /// returns true if the PseudoJet has constituents
00491   virtual bool has_constituents() const;
00492 
00493   /// retrieve the constituents. 
00494   ///
00495   /// an Error is thrown if this PseudoJet has no currently valid
00496   /// associated ClusterSequence or other substructure information
00497   virtual std::vector<PseudoJet> constituents() const;
00498 
00499 
00500   /// returns true if the PseudoJet has support for exclusive subjets
00501   virtual bool has_exclusive_subjets() const;
00502 
00503   /// return a vector of all subjets of the current jet (in the sense
00504   /// of the exclusive algorithm) that would be obtained when running
00505   /// the algorithm with the given dcut. 
00506   ///
00507   /// Time taken is O(m ln m), where m is the number of subjets that
00508   /// are found. If m gets to be of order of the total number of
00509   /// constituents in the jet, this could be substantially slower than
00510   /// just getting that list of constituents.
00511   ///
00512   /// an Error is thrown if this PseudoJet has no currently valid
00513   /// associated ClusterSequence
00514   std::vector<PseudoJet> exclusive_subjets (const double & dcut) const;
00515 
00516   /// return the size of exclusive_subjets(...); still n ln n with same
00517   /// coefficient, but marginally more efficient than manually taking
00518   /// exclusive_subjets.size()
00519   ///
00520   /// an Error is thrown if this PseudoJet has no currently valid
00521   /// associated ClusterSequence
00522   int n_exclusive_subjets(const double & dcut) const;
00523 
00524   /// return the list of subjets obtained by unclustering the supplied
00525   /// jet down to n subjets (or all constituents if there are fewer
00526   /// than n).
00527   ///
00528   /// requires n ln n time
00529   ///
00530   /// an Error is thrown if this PseudoJet has no currently valid
00531   /// associated ClusterSequence
00532   std::vector<PseudoJet> exclusive_subjets (int nsub) const;
00533 
00534   /// return the dij that was present in the merging nsub+1 -> nsub 
00535   /// subjets inside this jet.
00536   ///
00537   /// an Error is thrown if this PseudoJet has no currently valid
00538   /// associated ClusterSequence
00539   double exclusive_subdmerge(int nsub) const;
00540 
00541   /// return the maximum dij that occurred in the whole event at the
00542   /// stage that the nsub+1 -> nsub merge of subjets occurred inside 
00543   /// this jet.
00544   ///
00545   /// an Error is thrown if this PseudoJet has no currently valid
00546   /// associated ClusterSequence
00547   double exclusive_subdmerge_max(int nsub) const;
00548 
00549 
00550   /// returns true if a jet has pieces
00551   ///
00552   /// By default a single particle or a jet coming from a
00553   /// ClusterSequence have no pieces and this methos will return false.
00554   ///
00555   /// In practice, this is equivalent to have an structure of type
00556   /// CompositeJetStructure.
00557   virtual bool has_pieces() const;
00558 
00559 
00560   /// retrieve the pieces that make up the jet. 
00561   ///
00562   /// If the jet does not support pieces, an error is throw
00563   virtual std::vector<PseudoJet> pieces() const;
00564 
00565 
00566   // the following ones require a computation of the area in the
00567   // parent ClusterSequence (See ClusterSequenceAreaBase for details)
00568   //------------------------------------------------------------------
00569 
00570   /// check if it has a defined area
00571   virtual bool has_area() const;
00572 
00573   /// return the jet (scalar) area.
00574   /// throws an Error if there is no support for area in the parent CS
00575   virtual double area() const;
00576 
00577   /// return the error (uncertainty) associated with the determination
00578   /// of the area of this jet.
00579   /// throws an Error if there is no support for area in the parent CS
00580   virtual double area_error() const;
00581 
00582   /// return the jet 4-vector area.
00583   /// throws an Error if there is no support for area in the parent CS
00584   virtual PseudoJet area_4vector() const;
00585 
00586   /// true if this jet is made exclusively of ghosts.
00587   /// throws an Error if there is no support for area in the parent CS
00588   virtual bool is_pure_ghost() const;
00589 
00590   //\} --- end of jet structure -------------------------------------
00591 
00592 
00593 
00594   //----------------------------------------------------------------------
00595   /// @name Members mainly intended for internal use
00596   //----------------------------------------------------------------------
00597   //\{
00598   /// return the cluster_hist_index, intended to be used by clustering
00599   /// routines.
00600   inline int cluster_hist_index() const {return _cluster_hist_index;}
00601   /// set the cluster_hist_index, intended to be used by clustering routines.
00602   inline void set_cluster_hist_index(const int index) {_cluster_hist_index = index;}
00603 
00604   /// alternative name for cluster_hist_index() [perhaps more meaningful]
00605   inline int cluster_sequence_history_index() const {
00606     return cluster_hist_index();}
00607   /// alternative name for set_cluster_hist_index(...) [perhaps more
00608   /// meaningful]
00609   inline void set_cluster_sequence_history_index(const int index) {
00610     set_cluster_hist_index(index);}
00611 
00612   //\} ---- end of internal use functions ---------------------------
00613 
00614  protected:  
00615 
00616   SharedPtr<PseudoJetStructureBase> _structure;
00617   SharedPtr<UserInfoBase> _user_info;
00618 
00619 
00620  private: 
00621   // NB: following order must be kept for things to behave sensibly...
00622   double _px,_py,_pz,_E;
00623   mutable double _phi, _rap;
00624   double _kt2; 
00625   int    _cluster_hist_index, _user_index;
00626 
00627   /// calculate phi, rap, kt2 based on the 4-momentum components
00628   void _finish_init();
00629   /// set the indices to default values
00630   void _reset_indices();
00631 
00632   /// ensure that the internal values for rapidity and phi 
00633   /// correspond to 4-momentum structure
00634   inline void _ensure_valid_rap_phi() const {
00635     if (_phi == pseudojet_invalid_phi) _set_rap_phi();
00636   }
00637 
00638   /// set cached rapidity and phi values
00639   void _set_rap_phi() const;
00640 };
00641 
00642 
00643 //----------------------------------------------------------------------
00644 // routines for basic binary operations
00645 
00646 PseudoJet operator+(const PseudoJet &, const PseudoJet &);
00647 PseudoJet operator-(const PseudoJet &, const PseudoJet &);
00648 PseudoJet operator*(double, const PseudoJet &);
00649 PseudoJet operator*(const PseudoJet &, double);
00650 PseudoJet operator/(const PseudoJet &, double);
00651 
00652 inline double dot_product(const PseudoJet & a, const PseudoJet & b) {
00653   return a.E()*b.E() - a.px()*b.px() - a.py()*b.py() - a.pz()*b.pz();
00654 }
00655 
00656 /// returns true if the momenta of the two input jets are identical
00657 bool have_same_momentum(const PseudoJet &, const PseudoJet &);
00658 
00659 /// return a pseudojet with the given pt, y, phi and mass
00660 PseudoJet PtYPhiM(double pt, double y, double phi, double m = 0.0);
00661 
00662 //----------------------------------------------------------------------
00663 // Routines to do with providing sorted arrays of vectors.
00664 
00665 /// return a vector of jets sorted into decreasing transverse momentum
00666 std::vector<PseudoJet> sorted_by_pt(const std::vector<PseudoJet> & jets);
00667 
00668 /// return a vector of jets sorted into increasing rapidity
00669 std::vector<PseudoJet> sorted_by_rapidity(const std::vector<PseudoJet> & jets);
00670 
00671 /// return a vector of jets sorted into decreasing energy
00672 std::vector<PseudoJet> sorted_by_E(const std::vector<PseudoJet> & jets);
00673 
00674 /// return a vector of jets sorted into increasing pz
00675 std::vector<PseudoJet> sorted_by_pz(const std::vector<PseudoJet> & jets);
00676 
00677 //----------------------------------------------------------------------
00678 // some code to help sorting
00679 
00680 /// sort the indices so that values[indices[0->n-1]] is sorted
00681 /// into increasing order 
00682 void sort_indices(std::vector<int> & indices, 
00683                   const std::vector<double> & values);
00684 
00685 /// given a vector of values with a one-to-one correspondence with the
00686 /// vector of objects, sort objects into an order such that the
00687 /// associated values would be in increasing order (but don't actually
00688 /// touch the values vector in the process).
00689 template<class T> std::vector<T> objects_sorted_by_values(const std::vector<T> & objects, 
00690                                               const std::vector<double> & values);
00691 
00692 /// \if internal_doc
00693 /// @ingroup internal
00694 /// \class IndexedSortHelper
00695 /// a class that helps us carry out indexed sorting.
00696 /// \endif
00697 class IndexedSortHelper {
00698 public:
00699   inline IndexedSortHelper (const std::vector<double> * reference_values) {
00700     _ref_values = reference_values;
00701   };
00702   inline int operator() (const int & i1, const int & i2) const {
00703     return  (*_ref_values)[i1] < (*_ref_values)[i2];
00704   };
00705 private:
00706   const std::vector<double> * _ref_values;
00707 };
00708 
00709 
00710 //----------------------------------------------------------------------
00711 /// constructor from any object that has px,py,pz,E = some_four_vector[0--3],
00712 // NB: do not know if it really needs to be inline, but when it wasn't
00713 //     linking failed with g++ (who knows what was wrong...)
00714 template <class L> inline  PseudoJet::PseudoJet(const L & some_four_vector) {
00715   // now check whether L is simply a class that implements
00716   // some_fuor_vector[0--3] or actually is derived from PseudoJet and
00717   // has extra information
00718   DerivedPseudoJetHelper<L, IsBaseAndDerived<PseudoJet,L>::value> dpj_helper(some_four_vector);
00719 
00720   if (dpj_helper() != NULL){
00721     const PseudoJet *pj = dpj_helper();
00722     reset(*pj);
00723   } else {
00724     reset(some_four_vector);
00725   }
00726 }
00727 
00728 
00729 //----------------------------------------------------------------------
00730 inline void PseudoJet::_reset_indices() { 
00731   set_cluster_hist_index(-1);
00732   set_user_index(-1);
00733   _structure.reset();
00734 }
00735 
00736 //----------------------------------------------------------------------
00737 /// specialization of the "reset" template for case where something
00738 /// is reset to a pseudojet -- it then takes the user and history
00739 /// indices from the psjet
00740 // template<> inline void PseudoJet::reset<PseudoJet>(const PseudoJet & psjet) {
00741 //   (*this) = psjet;
00742 // }
00743 
00744 ////// fun and games...
00745 ////template<class L> class FJVector : public L {
00746 //////  /** Default Constructor: create jet with no constituents */
00747 //////  Vector<L>();
00748 ////
00749 ////};
00750 ////
00751 
00752 // taken literally from CLHEP
00753 inline double PseudoJet::m() const {
00754   double mm = m2();
00755   return mm < 0.0 ? -std::sqrt(-mm) : std::sqrt(mm);
00756 }
00757 
00758 
00759 inline void PseudoJet::reset(double px, double py, double pz, double E) {
00760   _px = px;
00761   _py = py;
00762   _pz = pz;
00763   _E  = E;
00764   _finish_init();
00765   _reset_indices();
00766 }
00767 
00768 
00769 //-------------------------------------------------------------------------------
00770 // implementation of the templated accesses to the underlying structyre
00771 //-------------------------------------------------------------------------------
00772 
00773 // returns a reference to the structure casted to the requested
00774 // structure type
00775 //
00776 // If there is no sructure associated, an Error is thrown.
00777 // If the type is not met, a std::bad_cast error is thrown.
00778 template<typename StructureType>
00779 const StructureType & PseudoJet::structure() const{
00780   return dynamic_cast<const StructureType &>(* validated_structure_ptr());
00781   
00782 }
00783 
00784 // check if the PseudoJet has the properties of the result of a Transformer 
00785 // (that is, its structure is compatible with a Transformer::StructureType)
00786 template<typename TransformerType>
00787 bool PseudoJet::has_properties_of() const{
00788   if (!_structure()) return false;
00789 
00790   return dynamic_cast<const typename TransformerType::StructureType *>(_structure.get()) != 0;
00791 }
00792 
00793 // this is a helper to access a structure created by a Transformer 
00794 // (that is, of type Transformer::StructureType)
00795 // NULL is returned if the corresponding type is not met
00796 template<typename TransformerType>
00797 const typename TransformerType::StructureType * PseudoJet::extra_properties() const{
00798   if (!_structure()) 
00799     throw Error("Trying to access the structure of a PseudoJet without an associated structure");
00800 
00801   return dynamic_cast<const typename TransformerType::StructureType *>(_structure.get());
00802 }
00803 
00804 
00805 
00806 //-------------------------------------------------------------------------------
00807 // helper functions to build a jet made of pieces
00808 //-------------------------------------------------------------------------------
00809 
00810 /// build a "CompositeJet" from the vector of its pieces
00811 ///
00812 /// In this case, E-scheme recombination is assumed to compute the
00813 /// total momentum
00814 PseudoJet join(const std::vector<PseudoJet> & pieces);
00815 
00816 /// build a MergedJet from a single PseudoJet
00817 PseudoJet join(const PseudoJet & j1);
00818 
00819 /// build a MergedJet from 2 PseudoJet
00820 PseudoJet join(const PseudoJet & j1, const PseudoJet & j2);
00821 
00822 /// build a MergedJet from 3 PseudoJet
00823 PseudoJet join(const PseudoJet & j1, const PseudoJet & j2, const PseudoJet & j3);
00824 
00825 /// build a MergedJet from 4 PseudoJet
00826 PseudoJet join(const PseudoJet & j1, const PseudoJet & j2, const PseudoJet & j3, const PseudoJet & j4);
00827 
00828 
00829 
00830 
00831 FASTJET_END_NAMESPACE
00832 
00833 #endif // __FASTJET_PSEUDOJET_HH__