32 #ifndef __FASTJET_PSEUDOJET_HH__
33 #define __FASTJET_PSEUDOJET_HH__
40 #include "fastjet/config.h"
41 #include "fastjet/internal/numconsts.hh"
42 #include "fastjet/internal/IsBase.hh"
43 #include "fastjet/SharedPtr.hh"
44 #include "fastjet/Error.hh"
45 #include "fastjet/PseudoJetStructureBase.hh"
47 FASTJET_BEGIN_NAMESPACE
57 const double pseudojet_invalid_rap = -1e200;
82 PseudoJet() : _px(0), _py(0), _pz(0), _E(0) {_finish_init(); _reset_indices();}
85 PseudoJet(
const double px,
const double py,
const double pz,
const double E);
89 template <
class L>
PseudoJet(
const L & some_four_vector);
98 #ifdef FASTJET_HAVE_THREAD_SAFETY
99 PseudoJet(
const PseudoJet &other){ (*this)=other; }
100 PseudoJet& operator=(
const PseudoJet& other);
116 inline double E()
const {
return _E;}
117 inline double e()
const {
return _E;}
118 inline double px()
const {
return _px;}
119 inline double py()
const {
return _py;}
120 inline double pz()
const {
return _pz;}
123 inline double phi()
const {
return phi_02pi();}
127 _ensure_valid_rap_phi();
128 return _phi > pi ? _phi-twopi : _phi;}
132 _ensure_valid_rap_phi();
138 inline double rap()
const {
139 _ensure_valid_rap_phi();
148 double pseudorapidity()
const;
149 double eta()
const {
return pseudorapidity();}
152 inline double pt2()
const {
return _kt2;}
154 inline double pt()
const {
return sqrt(_kt2);}
156 inline double perp2()
const {
return _kt2;}
158 inline double perp()
const {
return sqrt(_kt2);}
160 inline double kt2()
const {
return _kt2;}
163 inline double m2()
const {
return (_E+_pz)*(_E-_pz)-_kt2;}
166 inline double m()
const;
169 inline double mperp2()
const {
return (_E+_pz)*(_E-_pz);}
171 inline double mperp()
const {
return sqrt(std::abs(mperp2()));}
173 inline double mt2()
const {
return (_E+_pz)*(_E-_pz);}
175 inline double mt()
const {
return sqrt(std::abs(mperp2()));}
178 inline double modp2()
const {
return _kt2+_pz*_pz;}
180 inline double modp()
const {
return sqrt(_kt2+_pz*_pz);}
183 inline double Et()
const {
return (_kt2==0) ? 0.0 : _E/sqrt(1.0+_pz*_pz/_kt2);}
185 inline double Et2()
const {
return (_kt2==0) ? 0.0 : _E*_E/(1.0+_pz*_pz/_kt2);}
190 return std::min(1.0, std::max(-1.0, _pz/sqrt(modp2())));
196 double operator () (
int i)
const ;
198 inline double operator [] (
int i)
const {
return (*
this)(i); };
203 double kt_distance(
const PseudoJet & other)
const;
206 double plain_distance(
const PseudoJet & other)
const;
210 return plain_distance(other);}
215 return sqrt(squared_distance(other));
220 double delta_phi_to(
const PseudoJet & other)
const;
232 std::valarray<double> four_mom()
const;
237 enum { X=0, Y=1, Z=2, T=3, NUM_COORDINATES=4, SIZE=NUM_COORDINATES };
246 PseudoJet & boost(
const PseudoJet & prest);
249 PseudoJet & unboost(
const PseudoJet & prest);
251 PseudoJet & operator*=(
double);
252 PseudoJet & operator/=(
double);
253 PseudoJet & operator+=(
const PseudoJet &);
254 PseudoJet & operator-=(
const PseudoJet &);
290 inline void reset(
double px,
double py,
double pz,
double E);
307 template <
class L>
inline void reset(
const L & some_four_vector) {
319 const PseudoJet * pj = fastjet::cast_if_derived<const PseudoJet>(&some_four_vector);
324 reset(some_four_vector[0], some_four_vector[1],
325 some_four_vector[2], some_four_vector[3]);
333 inline void reset_PtYPhiM(
double pt_in,
double y_in,
double phi_in,
double m_in=0.0) {
334 reset_momentum_PtYPhiM(pt_in, y_in, phi_in, m_in);
342 inline void reset_momentum(
double px,
double py,
double pz,
double E);
347 inline void reset_momentum(
const PseudoJet & pj);
351 void reset_momentum_PtYPhiM(
double pt,
double y,
double phi,
double m=0.0);
357 reset_momentum(some_four_vector[0], some_four_vector[1],
358 some_four_vector[2], some_four_vector[3]);
373 void set_cached_rap_phi(
double rap,
double phi);
446 _user_info.reset(user_info_in);
480 return dynamic_cast<const L &
>(* _user_info.get());
485 return _user_info.get();
492 return _user_info.get() &&
dynamic_cast<const L *
>(_user_info.get());
500 return _user_info.get();
538 std::string description()
const;
552 bool has_associated_cluster_sequence()
const;
558 bool has_valid_cluster_sequence()
const;
566 const ClusterSequence* associated_cs()
const {
return associated_cluster_sequence();}
571 return validated_cs();
580 return validated_csab();
603 bool has_structure()
const;
637 template<
typename StructureType>
638 const StructureType & structure()
const;
643 template<
typename TransformerType>
644 bool has_structure_of()
const;
651 template<
typename TransformerType>
652 const typename TransformerType::StructureType & structure_of()
const;
671 virtual bool has_partner(
PseudoJet &partner)
const;
679 virtual bool has_child(
PseudoJet &child)
const;
694 virtual bool contains(
const PseudoJet &constituent)
const;
701 virtual bool is_inside(
const PseudoJet &jet)
const;
705 virtual bool has_constituents()
const;
711 virtual std::vector<PseudoJet> constituents()
const;
715 virtual bool has_exclusive_subjets()
const;
728 std::vector<PseudoJet> exclusive_subjets (
const double dcut)
const;
736 int n_exclusive_subjets(
const double dcut)
const;
746 std::vector<PseudoJet> exclusive_subjets (
int nsub)
const;
756 std::vector<PseudoJet> exclusive_subjets_up_to (
int nsub)
const;
765 double exclusive_subdmerge(
int nsub)
const;
775 double exclusive_subdmerge_max(
int nsub)
const;
785 virtual bool has_pieces()
const;
791 virtual std::vector<PseudoJet> pieces()
const;
800 virtual bool has_area()
const;
804 virtual double area()
const;
809 virtual double area_error()
const;
817 virtual bool is_pure_ghost()
const;
836 return cluster_hist_index();}
840 set_cluster_hist_index(index);}
852 double _px,_py,_pz,_E;
853 mutable double _phi, _rap;
855 int _cluster_hist_index, _user_index;
857 #ifdef FASTJET_HAVE_THREAD_SAFETY
864 mutable std::atomic<int> _init_status;
870 void _reset_indices();
891 #ifdef FASTJET_HAVE_THREAD_SAFETY
892 void _ensure_valid_rap_phi()
const;
894 inline void _ensure_valid_rap_phi()
const{
895 if (_phi == pseudojet_invalid_phi) _set_rap_phi();
900 void _set_rap_phi()
const;
903 friend PseudoJet
operator*(
double,
const PseudoJet &);
910 PseudoJet operator+(
const PseudoJet &,
const PseudoJet &);
911 PseudoJet operator-(
const PseudoJet &,
const PseudoJet &);
912 PseudoJet
operator*(
double,
const PseudoJet &);
913 PseudoJet
operator*(
const PseudoJet &,
double);
914 PseudoJet operator/(
const PseudoJet &,
double);
919 bool operator==(
const PseudoJet &,
const PseudoJet &);
926 bool operator==(
const PseudoJet & jet,
const double val);
927 inline bool operator==(
const double val,
const PseudoJet & jet) {
return jet == val;}
932 inline bool operator!=(
const double val,
const PseudoJet & a) {
return !(a==val);}
936 return a.E()*b.E() - a.px()*b.px() - a.py()*b.py() - a.pz()*b.pz();
941 double dot_3d = a.px()*b.px() + a.py()*b.py() + a.pz()*b.pz();
942 return std::min(1.0, std::max(-1.0, dot_3d/sqrt(a.
modp2()*b.
modp2())));
951 bool have_same_momentum(
const PseudoJet &,
const PseudoJet &);
955 PseudoJet PtYPhiM(
double pt,
double y,
double phi,
double m = 0.0);
961 std::vector<PseudoJet>
sorted_by_pt(
const std::vector<PseudoJet> & jets);
967 std::vector<PseudoJet>
sorted_by_E(
const std::vector<PseudoJet> & jets);
970 std::vector<PseudoJet>
sorted_by_pz(
const std::vector<PseudoJet> & jets);
977 void sort_indices(std::vector<int> & indices,
978 const std::vector<double> & values);
985 const std::vector<double> & values) {
987 if (objects.size() != values.size()){
988 throw Error(
"fastjet::objects_sorted_by_values(...): the size of the 'objects' vector must match the size of the 'values' vector");
992 std::vector<int> indices(values.size());
993 for (
size_t i = 0; i < indices.size(); i++) {indices[i] = i;}
996 sort_indices(indices, values);
999 std::vector<T> objects_sorted(objects.size());
1002 for (
size_t i = 0; i < indices.size(); i++) {
1003 objects_sorted[i] = objects[indices[i]];
1006 return objects_sorted;
1014 class IndexedSortHelper {
1016 inline IndexedSortHelper (
const std::vector<double> * reference_values) {
1017 _ref_values = reference_values;
1019 inline int operator() (
const int i1,
const int i2)
const {
1020 return (*_ref_values)[i1] < (*_ref_values)[i2];
1023 const std::vector<double> * _ref_values;
1032 template <
class L>
inline PseudoJet::PseudoJet(
const L & some_four_vector) {
1033 reset(some_four_vector);
1038 inline void PseudoJet::_reset_indices() {
1039 set_cluster_hist_index(-1);
1060 inline double PseudoJet::m()
const {
1062 return mm < 0.0 ? -std::sqrt(-mm) : std::sqrt(mm);
1066 inline void PseudoJet::reset(
double px_in,
double py_in,
double pz_in,
double E_in) {
1076 inline void PseudoJet::reset_momentum(
double px_in,
double py_in,
double pz_in,
double E_in) {
1084 inline void PseudoJet::reset_momentum(
const PseudoJet & pj) {
1103 template<
typename StructureType>
1104 const StructureType & PseudoJet::structure()
const{
1105 return dynamic_cast<const StructureType &
>(* validated_structure_ptr());
1111 template<
typename TransformerType>
1112 bool PseudoJet::has_structure_of()
const{
1113 if (!_structure)
return false;
1115 return dynamic_cast<const typename TransformerType::StructureType *
>(_structure.
get()) != 0;
1121 template<
typename TransformerType>
1122 const typename TransformerType::StructureType & PseudoJet::structure_of()
const{
1124 throw Error(
"Trying to access the structure of a PseudoJet without an associated structure");
1126 return dynamic_cast<const typename TransformerType::StructureType &
>(*_structure);
1143 PseudoJet join(
const std::vector<PseudoJet> & pieces);
1159 FASTJET_END_NAMESPACE
base class that sets interface for extensions of ClusterSequence that provide information about the a...
base class corresponding to errors that can be thrown by FastJet
Contains any information related to the clustering that should be directly accessible to PseudoJet.
error class to be thrown if accessing user info when it doesn't exist
a base class to hold extra user information in a PseudoJet
Class to contain pseudojets, including minimal information of use to jet-clustering routines.
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...
PseudoJet()
default constructor, which as of FJ3.0 provides an object for which all operations are now valid and ...
double cos_theta() const
cos of the polar angle should we have: min(1.0,max(-1.0,_pz/sqrt(modp2())));
double mt() const
returns the transverse mass = sqrt(kt^2+m^2)
int cluster_sequence_history_index() const
alternative name for cluster_hist_index() [perhaps more meaningful]
double rap() const
returns the rapidity or some large value when the rapidity is infinite
double Et() const
return the transverse energy
double delta_R(const PseudoJet &other) const
return the cylinder (rap-phi) distance between this jet and another, .
const UserInfoBase * user_info_ptr() const
retrieve a pointer to the (const) user information
double mperp() const
returns the transverse mass = sqrt(kt^2+m^2)
double modp() const
return the 3-vector modulus = sqrt(px^2+py^2+pz^2)
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...
void reset_momentum(const L &some_four_vector)
reset the 4-momentum according to the supplied generic 4-vector (accessible via indexing,...
bool has_user_info() const
returns true if the PseudoJet has user information
double squared_distance(const PseudoJet &other) const
returns squared cylinder (rap-phi) distance between this jet and another
bool has_valid_cs() const
shorthand for has_valid_cluster_sequence()
double phi() const
returns phi (in the range 0..2pi)
virtual ~PseudoJet()
default (virtual) destructor
double perp() const
returns the scalar transverse momentum
double perp2() const
returns the squared transverse momentum
SharedPtr< UserInfoBase > & user_info_shared_ptr()
retrieve a (non-const) shared pointer to the user information; you can use this, for example,...
double m2() const
returns the squared invariant mass // like CLHEP
void set_cluster_hist_index(const int index)
set the cluster_hist_index, intended to be used by clustering routines.
double kt2() const
returns the squared transverse momentum
bool has_associated_cs() const
shorthand for has_associated_cluster_sequence()
bool has_user_info() const
returns true if the PseudoJet has user information than can be cast to the template argument type.
double theta() const
polar angle
const L & user_info() const
returns a reference to the dynamic cast conversion of user_info to type L.
double phi_std() const
returns phi in the range -pi..pi
double phi_02pi() const
returns phi in the range 0..2pi
int cluster_hist_index() const
return the cluster_hist_index, intended to be used by clustering routines.
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...
double pt() const
returns the scalar transverse momentum
void set_cluster_sequence_history_index(const int index)
alternative name for set_cluster_hist_index(...) [perhaps more meaningful]
double mperp2() const
returns the squared transverse mass = kt^2+m^2
void reset(const L &some_four_vector)
reset the 4-momentum according to the supplied generic 4-vector (accessible via indexing,...
int user_index() const
return the user_index,
double modp2() const
return the squared 3-vector modulus = px^2+py^2+pz^2
const SharedPtr< UserInfoBase > & user_info_shared_ptr() const
retrieve a (const) shared pointer to the user information
double Et2() const
return the transverse energy squared
double rapidity() const
the same as rap()
double beam_distance() const
returns distance between this jet and the beam
double pt2() const
returns the squared transverse momentum
void set_user_info(UserInfoBase *user_info_in)
sets the internal shared pointer to the user information.
const ClusterSequenceAreaBase * validated_cluster_sequence_area_base() const
if the jet has valid area information then return a pointer to the associated ClusterSequenceAreaBase...
void reset(const PseudoJet &psjet)
reset the PseudoJet to be equal to psjet (including its indices); NB if the argument is derived from ...
double mt2() const
returns the squared transverse mass = kt^2+m^2
An implementation of shared pointers that is broadly similar to C++11 shared_ptr (https://en....
T * get() const
get the stored pointer
void reset()
reset the pointer to default value (NULL)
Selector operator*(const Selector &s1, const Selector &s2)
successive application of 2 selectors
bool operator!=(const PseudoJet &a, const PseudoJet &b)
inequality test which is exact opposite of operator==
vector< PseudoJet > sorted_by_pz(const vector< PseudoJet > &jets)
return a vector of jets sorted into increasing pz
vector< PseudoJet > sorted_by_pt(const vector< PseudoJet > &jets)
return a vector of jets sorted into decreasing kt2
const double MaxRap
Used to protect against parton-level events where pt can be zero for some partons,...
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,...
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...
vector< PseudoJet > sorted_by_rapidity(const vector< PseudoJet > &jets)
return a vector of jets sorted into increasing rapidity
double dot_product(const PseudoJet &a, const PseudoJet &b)
returns the 4-vector dot product of a and b
const double pseudojet_invalid_phi
default value for phi, meaning it (and rapidity) have yet to be calculated)
double theta(const PseudoJet &a, const PseudoJet &b)
returns the angle between a and b
double cos_theta(const PseudoJet &a, const PseudoJet &b)
returns the cosine of the angle between a and b
vector< PseudoJet > sorted_by_E(const vector< PseudoJet > &jets)
return a vector of jets sorted into decreasing energy