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);}
175 double operator () (
int i)
const ;
177 inline double operator [] (
int i)
const {
return (*
this)(i); };
182 double kt_distance(
const PseudoJet & other)
const;
185 double plain_distance(
const PseudoJet & other)
const;
189 return plain_distance(other);}
194 return sqrt(squared_distance(other));
199 double delta_phi_to(
const PseudoJet & other)
const;
211 std::valarray<double> four_mom()
const;
216 enum { X=0, Y=1, Z=2, T=3, NUM_COORDINATES=4, SIZE=NUM_COORDINATES };
237 inline void reset(
double px,
double py,
double pz,
double E);
254 template <
class L>
inline void reset(
const L & some_four_vector) {
266 const PseudoJet * pj = fastjet::cast_if_derived<const PseudoJet>(&some_four_vector);
271 reset(some_four_vector[0], some_four_vector[1],
272 some_four_vector[2], some_four_vector[3]);
280 inline void reset_PtYPhiM(
double pt_in,
double y_in,
double phi_in,
double m_in=0.0) {
281 reset_momentum_PtYPhiM(pt_in, y_in, phi_in, m_in);
288 inline void reset_momentum(
double px,
double py,
double pz,
double E);
293 inline void reset_momentum(
const PseudoJet & pj);
297 void reset_momentum_PtYPhiM(
double pt,
double y,
double phi,
double m=0.0);
303 reset_momentum(some_four_vector[0], some_four_vector[1],
304 some_four_vector[2], some_four_vector[3]);
319 void set_cached_rap_phi(
double rap,
double phi);
392 _user_info.reset(user_info_in);
426 return dynamic_cast<const L &
>(* _user_info.get());
431 return _user_info.get();
438 return _user_info.get() &&
dynamic_cast<const L *
>(_user_info.get());
446 return _user_info.get();
484 std::string description()
const;
498 bool has_associated_cluster_sequence()
const;
504 bool has_valid_cluster_sequence()
const;
512 const ClusterSequence* associated_cs()
const {
return associated_cluster_sequence();}
517 return validated_cs();
526 return validated_csab();
549 bool has_structure()
const;
583 template<
typename StructureType>
584 const StructureType & structure()
const;
589 template<
typename TransformerType>
590 bool has_structure_of()
const;
597 template<
typename TransformerType>
598 const typename TransformerType::StructureType & structure_of()
const;
617 virtual bool has_partner(
PseudoJet &partner)
const;
625 virtual bool has_child(
PseudoJet &child)
const;
640 virtual bool contains(
const PseudoJet &constituent)
const;
647 virtual bool is_inside(
const PseudoJet &jet)
const;
651 virtual bool has_constituents()
const;
657 virtual std::vector<PseudoJet> constituents()
const;
661 virtual bool has_exclusive_subjets()
const;
674 std::vector<PseudoJet> exclusive_subjets (
const double dcut)
const;
682 int n_exclusive_subjets(
const double dcut)
const;
692 std::vector<PseudoJet> exclusive_subjets (
int nsub)
const;
702 std::vector<PseudoJet> exclusive_subjets_up_to (
int nsub)
const;
711 double exclusive_subdmerge(
int nsub)
const;
721 double exclusive_subdmerge_max(
int nsub)
const;
731 virtual bool has_pieces()
const;
737 virtual std::vector<PseudoJet> pieces()
const;
746 virtual bool has_area()
const;
750 virtual double area()
const;
755 virtual double area_error()
const;
763 virtual bool is_pure_ghost()
const;
782 return cluster_hist_index();}
786 set_cluster_hist_index(index);}
798 double _px,_py,_pz,_E;
799 mutable double _phi, _rap;
801 int _cluster_hist_index, _user_index;
806 void _reset_indices();
810 inline void _ensure_valid_rap_phi()
const {
815 void _set_rap_phi()
const;
851 return a.E()*b.E() - a.px()*b.px() - a.py()*b.py() - a.pz()*b.pz();
865 std::vector<PseudoJet>
sorted_by_pt(
const std::vector<PseudoJet> & jets);
871 std::vector<PseudoJet>
sorted_by_E(
const std::vector<PseudoJet> & jets);
874 std::vector<PseudoJet>
sorted_by_pz(
const std::vector<PseudoJet> & jets);
881 void sort_indices(std::vector<int> & indices,
882 const std::vector<double> & values);
889 const std::vector<double> & values) {
891 if (objects.size() != values.size()){
892 throw Error(
"fastjet::objects_sorted_by_values(...): the size of the 'objects' vector must match the size of the 'values' vector");
896 std::vector<int> indices(values.size());
897 for (
size_t i = 0; i < indices.size(); i++) {indices[i] = i;}
900 sort_indices(indices, values);
903 std::vector<T> objects_sorted(objects.size());
906 for (
size_t i = 0; i < indices.size(); i++) {
907 objects_sorted[i] = objects[indices[i]];
910 return objects_sorted;
918 class IndexedSortHelper {
920 inline IndexedSortHelper (
const std::vector<double> * reference_values) {
921 _ref_values = reference_values;
923 inline int operator() (
const int i1,
const int i2)
const {
924 return (*_ref_values)[i1] < (*_ref_values)[i2];
927 const std::vector<double> * _ref_values;
936 template <
class L>
inline PseudoJet::PseudoJet(
const L & some_four_vector) {
937 reset(some_four_vector);
942 inline void PseudoJet::_reset_indices() {
943 set_cluster_hist_index(-1);
951 inline double PseudoJet::m()
const {
953 return mm < 0.0 ? -std::sqrt(-mm) : std::sqrt(mm);
957 inline void PseudoJet::reset(
double px_in,
double py_in,
double pz_in,
double E_in) {
966 inline void PseudoJet::reset_momentum(
double px_in,
double py_in,
double pz_in,
double E_in) {
974 inline void PseudoJet::reset_momentum(
const PseudoJet & pj) {
993 template<
typename StructureType>
994 const StructureType & PseudoJet::structure()
const{
995 return dynamic_cast<const StructureType &
>(* validated_structure_ptr());
1001 template<
typename TransformerType>
1002 bool PseudoJet::has_structure_of()
const{
1003 if (!_structure)
return false;
1005 return dynamic_cast<const typename TransformerType::StructureType *
>(_structure.get()) != 0;
1011 template<
typename TransformerType>
1012 const typename TransformerType::StructureType & PseudoJet::structure_of()
const{
1014 throw Error(
"Trying to access the structure of a PseudoJet without an associated structure");
1016 return dynamic_cast<const typename TransformerType::StructureType &
>(*_structure);
1033 PseudoJet join(
const std::vector<PseudoJet> & pieces);
1049 FASTJET_END_NAMESPACE
1051 #endif // __FASTJET_PSEUDOJET_HH__ PseudoJet PtYPhiM(double pt, double y, double phi, double m)
return a pseudojet with the given pt, y, phi and mass
void reset_momentum(const L &some_four_vector)
reset the 4-momentum according to the supplied generic 4-vector (accessible via indexing, [0]==px,...[3]==E), but leave all other information (indices, user info, etc.) untouched
double Et2() const
return the transverse energy squared
double m2() const
returns the squared invariant mass // like CLHEP
double mt2() const
returns the squared transverse mass = kt^2+m^2
vector< PseudoJet > sorted_by_pt(const vector< PseudoJet > &jets)
return a vector of jets sorted into decreasing kt2
std::vector< T > objects_sorted_by_values(const std::vector< T > &objects, const std::vector< double > &values)
given a vector of values with a one-to-one correspondence with the vector of objects, sort objects into an order such that the associated values would be in increasing order (but don't actually touch the values vector in the process).
const ClusterSequence * validated_cluster_sequence() const
if the jet has a valid associated cluster sequence then return a pointer to it; otherwise throw an er...
void set_user_index(const int index)
set the user_index, intended to allow the user to add simple identifying information to a particle/je...
double Et() const
return the transverse energy
const double MaxRap
Used to protect against parton-level events where pt can be zero for some partons, giving rapidity=infinity.
void reset(const PseudoJet &psjet)
reset the PseudoJet to be equal to psjet (including its indices); NB if the argument is derived from ...
int cluster_sequence_history_index() const
alternative name for cluster_hist_index() [perhaps more meaningful]
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.
bool has_valid_cs() const
shorthand for has_valid_cluster_sequence()
const UserInfoBase * user_info_ptr() const
retrieve a pointer to the (const) user information
Selector operator*(const Selector &s1, const Selector &s2)
successive application of 2 selectors
vector< PseudoJet > sorted_by_pz(const vector< PseudoJet > &jets)
return a vector of jets sorted into increasing pz
double phi_02pi() const
returns phi in the range 0..2pi
double squared_distance(const PseudoJet &other) const
returns squared cylinder (rap-phi) distance between this jet and another
void reset_PtYPhiM(double pt_in, double y_in, double phi_in, double m_in=0.0)
reset the PseudoJet according to the specified pt, rapidity, azimuth and mass (also resetting indices...
void set_user_info(UserInfoBase *user_info_in)
sets the internal shared pointer to the user information.
vector< PseudoJet > sorted_by_rapidity(const vector< PseudoJet > &jets)
return a vector of jets sorted into increasing rapidity
double pt() const
returns the scalar transverse momentum
double rapidity() const
the same as rap()
Contains any information related to the clustering that should be directly accessible to PseudoJet...
void set_cluster_sequence_history_index(const int index)
alternative name for set_cluster_hist_index(...) [perhaps more meaningful]
bool has_user_info() const
returns true if the PseudoJet has user information than can be cast to the template argument type...
double phi_std() const
returns phi in the range -pi..pi
double modp2() const
return the squared 3-vector modulus = px^2+py^2+pz^2
bool operator==(const PseudoJet &a, const PseudoJet &b)
returns true if the 4 momentum components of the two PseudoJets are identical and all the internal in...
const SharedPtr< UserInfoBase > & user_info_shared_ptr() const
retrieve a (const) shared pointer to the user information
virtual ~PseudoJet()
default (virtual) destructor
double dot_product(const PseudoJet &a, const PseudoJet &b)
returns the 4-vector dot product of a and b
double perp() const
returns the scalar transverse momentum
void set_cluster_hist_index(const int index)
set the cluster_hist_index, intended to be used by clustering routines.
bool has_associated_cs() const
shorthand for has_associated_cluster_sequence()
base class that sets interface for extensions of ClusterSequence that provide information about the a...
error class to be thrown if accessing user info when it doesn't exist
vector< PseudoJet > sorted_by_E(const vector< PseudoJet > &jets)
return a vector of jets sorted into decreasing energy
bool has_user_info() const
returns true if the PseudoJet has user information
base class corresponding to errors that can be thrown by FastJet
double kt2() const
returns the squared transverse momentum
int cluster_hist_index() const
return the cluster_hist_index, intended to be used by clustering routines.
double perp2() const
returns the squared transverse momentum
double delta_R(const PseudoJet &other) const
return the cylinder (rap-phi) distance between this jet and another, .
an implementation of C++0x shared pointers (or boost's)
bool operator!=(const PseudoJet &a, const PseudoJet &b)
inequality test which is exact opposite of operator==
double phi() const
returns phi (in the range 0..2pi)
double rap() const
returns the rapidity or some large value when the rapidity is infinite
double modp() const
return the 3-vector modulus = sqrt(px^2+py^2+pz^2)
const L & user_info() const
returns a reference to the dynamic cast conversion of user_info to type L.
bool have_same_momentum(const PseudoJet &jeta, const PseudoJet &jetb)
returns true if the momenta of the two input jets are identical
double mt() const
returns the transverse mass = sqrt(kt^2+m^2)
a base class to hold extra user information in a PseudoJet
double mperp() const
returns the transverse mass = sqrt(kt^2+m^2)
double pt2() const
returns the squared transverse momentum
PseudoJet()
default constructor, which as of FJ3.0 provides an object for which all operations are now valid and ...
const double pseudojet_invalid_phi
default value for phi, meaning it (and rapidity) have yet to be calculated)
Class to contain pseudojets, including minimal information of use to jet-clustering routines...
const ClusterSequenceAreaBase * validated_cluster_sequence_area_base() const
if the jet has valid area information then return a pointer to the associated ClusterSequenceAreaBase...
int user_index() const
return the user_index,
double beam_distance() const
returns distance between this jet and the beam
SharedPtr< UserInfoBase > & user_info_shared_ptr()
retrieve a (non-const) shared pointer to the user information; you can use this, for example...
double mperp2() const
returns the squared transverse mass = kt^2+m^2