32 #include "fastjet/Error.hh" 33 #include "fastjet/PseudoJet.hh" 34 #include "fastjet/ClusterSequence.hh" 36 #include "fastjet/ClusterSequenceAreaBase.hh" 38 #include "fastjet/CompositeJetStructure.hh" 46 FASTJET_BEGIN_NAMESPACE
53 PseudoJet::PseudoJet(
const double px_in,
const double py_in,
const double pz_in,
const double E_in) {
70 void PseudoJet::_finish_init () {
71 _kt2 = this->px()*this->px() + this->py()*this->py();
80 _rap = pseudojet_invalid_rap;
84 void PseudoJet::_set_rap_phi()
const {
89 _phi = atan2(this->py(),this->px());
91 if (_phi < 0.0) {_phi += twopi;}
92 if (_phi >= twopi) {_phi -= twopi;}
93 if (this->E() == abs(this->pz()) && _kt2 == 0) {
98 double MaxRapHere =
MaxRap + abs(this->pz());
99 if (this->pz() >= 0.0) {_rap = MaxRapHere;}
else {_rap = -MaxRapHere;}
104 double effective_m2 = max(0.0,m2());
105 double E_plus_pz = _E + abs(_pz);
107 _rap = 0.5*log((_kt2 + effective_m2)/(E_plus_pz*E_plus_pz));
108 if (_pz > 0) {_rap = - _rap;}
116 valarray<double> PseudoJet::four_mom()
const {
117 valarray<double> mom(4);
128 double PseudoJet::operator () (
int i)
const {
140 err <<
"PseudoJet subscripting: bad index (" << i <<
")";
141 throw Error(err.str());
148 double PseudoJet::pseudorapidity()
const {
149 if (px() == 0.0 && py() ==0.0)
return MaxRap;
150 if (pz() == 0.0)
return 0.0;
152 double theta = atan(perp()/pz());
153 if (theta < 0) theta += pi;
154 return -log(tan(theta/2));
164 jet1.E() +jet2.E() );
174 jet1.E() -jet2.E() );
183 jet._ensure_valid_rap_phi();
187 coeff_times_jet *= coeff;
188 return coeff_times_jet;
200 return (1.0/coeff)*jet;
213 _ensure_valid_rap_phi();
226 (*this) *= 1.0/coeff;
234 _px += other_jet._px;
235 _py += other_jet._py;
236 _pz += other_jet._pz;
246 _px -= other_jet._px;
247 _py -= other_jet._py;
248 _pz -= other_jet._pz;
256 if (a.px() != b.px())
return false;
257 if (a.py() != b.py())
return false;
258 if (a.pz() != b.pz())
return false;
259 if (a.E () != b.E ())
return false;
273 throw Error(
"comparing a PseudoJet with a non-zero constant (double) is not allowed.");
274 return (jet.px() == 0 && jet.py() == 0 &&
275 jet.pz() == 0 && jet.E() == 0);
288 if (prest.px() == 0.0 && prest.py() == 0.0 && prest.pz() == 0.0)
291 double m_local = prest.
m();
292 assert(m_local != 0);
294 double pf4 = ( px()*prest.px() + py()*prest.py()
295 + pz()*prest.pz() + E()*prest.E() )/m_local;
296 double fn = (pf4 + E()) / (prest.E() + m_local);
297 _px += fn*prest.px();
298 _py += fn*prest.py();
299 _pz += fn*prest.pz();
315 if (prest.px() == 0.0 && prest.py() == 0.0 && prest.pz() == 0.0)
318 double m_local = prest.
m();
319 assert(m_local != 0);
321 double pf4 = ( -px()*prest.px() - py()*prest.py()
322 - pz()*prest.pz() + E()*prest.E() )/m_local;
323 double fn = (pf4 + E()) / (prest.E() + m_local);
324 _px -= fn*prest.px();
325 _py -= fn*prest.py();
326 _pz -= fn*prest.pz();
337 return jeta.px() == jetb.px()
338 && jeta.py() == jetb.py()
339 && jeta.pz() == jetb.pz()
340 && jeta.E() == jetb.E();
344 void PseudoJet::set_cached_rap_phi(
double rap_in,
double phi_in) {
345 _rap = rap_in; _phi = phi_in;
346 if (_phi >= twopi) _phi -= twopi;
347 if (_phi < 0) _phi += twopi;
351 void PseudoJet::reset_momentum_PtYPhiM(
double pt_in,
double y_in,
double phi_in,
double m_in) {
352 assert(phi_in < 2*twopi && phi_in > -twopi);
353 double ptm = (m_in == 0) ? pt_in : sqrt(pt_in*pt_in+m_in*m_in);
354 double exprap = exp(y_in);
355 double pminus = ptm/exprap;
356 double pplus = ptm*exprap;
357 double px_local = pt_in*cos(phi_in);
358 double py_local = pt_in*sin(phi_in);
359 reset_momentum(px_local,py_local,0.5*(pplus-pminus),0.5*(pplus+pminus));
360 set_cached_rap_phi(y_in,phi_in);
366 assert(phi < 2*twopi && phi > -twopi);
367 double ptm = (m == 0) ? pt : sqrt(pt*pt+m*m);
368 double exprap = exp(y);
369 double pminus = ptm/exprap;
370 double pplus = ptm*exprap;
371 double px = pt*cos(phi);
372 double py = pt*sin(phi);
373 PseudoJet mom(px,py,0.5*(pplus-pminus),0.5*(pplus+pminus));
382 double PseudoJet::kt_distance(
const PseudoJet & other)
const {
384 double distance = min(_kt2, other._kt2);
385 double dphi = abs(phi() - other.
phi());
386 if (dphi > pi) {dphi = twopi - dphi;}
387 double drap = rap() - other.
rap();
388 distance = distance * (dphi*dphi + drap*drap);
395 double PseudoJet::plain_distance(
const PseudoJet & other)
const {
396 double dphi = abs(phi() - other.
phi());
397 if (dphi > pi) {dphi = twopi - dphi;}
398 double drap = rap() - other.
rap();
399 return (dphi*dphi + drap*drap);
405 double PseudoJet::delta_phi_to(
const PseudoJet & other)
const {
406 double dphi = other.
phi() - phi();
407 if (dphi > pi) dphi -= twopi;
408 if (dphi < -pi) dphi += twopi;
413 string PseudoJet::description()
const{
416 return "standard PseudoJet (with no associated clustering information)";
419 return _structure->description();
434 bool PseudoJet::has_associated_cluster_sequence()
const{
435 return (_structure) && (_structure->has_associated_cluster_sequence());
442 if (! has_associated_cluster_sequence())
return NULL;
444 return _structure->associated_cluster_sequence();
451 bool PseudoJet::has_valid_cluster_sequence()
const{
452 return (_structure) && (_structure->has_valid_cluster_sequence());
462 return validated_structure_ptr()->validated_cs();
469 _structure = structure_in;
474 bool PseudoJet::has_structure()
const{
475 return bool(_structure);
485 return _structure.get();
500 return _structure.get();
510 throw Error(
"Trying to access the structure of a PseudoJet which has no associated structure");
511 return _structure.get();
530 return validated_structure_ptr()->has_partner(*
this, partner);
541 return validated_structure_ptr()->has_child(*
this, child);
552 return validated_structure_ptr()->has_parents(*
this, parent1, parent2);
561 bool PseudoJet::contains(
const PseudoJet &constituent)
const{
562 return validated_structure_ptr()->object_in_jet(constituent, *
this);
572 return validated_structure_ptr()->object_in_jet(*
this, jet);
578 bool PseudoJet::has_constituents()
const{
579 return (_structure) && (_structure->has_constituents());
584 vector<PseudoJet> PseudoJet::constituents()
const{
585 return validated_structure_ptr()->constituents(*
this);
591 bool PseudoJet::has_exclusive_subjets()
const{
592 return (_structure) && (_structure->has_exclusive_subjets());
607 std::vector<PseudoJet> PseudoJet::exclusive_subjets (
const double dcut)
const {
608 return validated_structure_ptr()->exclusive_subjets(*
this, dcut);
618 int PseudoJet::n_exclusive_subjets(
const double dcut)
const {
619 return validated_structure_ptr()->n_exclusive_subjets(*
this, dcut);
631 std::vector<PseudoJet> PseudoJet::exclusive_subjets_up_to (
int nsub)
const {
632 return validated_structure_ptr()->exclusive_subjets_up_to(*
this, nsub);
638 std::vector<PseudoJet> PseudoJet::exclusive_subjets (
int nsub)
const {
639 vector<PseudoJet> subjets = exclusive_subjets_up_to(nsub);
640 if (
int(subjets.size()) < nsub) {
642 err <<
"Requested " << nsub <<
" exclusive subjets, but there were only " 643 << subjets.size() <<
" particles in the jet";
644 throw Error(err.str());
655 double PseudoJet::exclusive_subdmerge(
int nsub)
const {
656 return validated_structure_ptr()->exclusive_subdmerge(*
this, nsub);
666 double PseudoJet::exclusive_subdmerge_max(
int nsub)
const {
667 return validated_structure_ptr()->exclusive_subdmerge_max(*
this, nsub);
675 bool PseudoJet::has_pieces()
const{
676 return ((_structure) && (_structure->has_pieces(*
this)));
684 std::vector<PseudoJet> PseudoJet::pieces()
const{
685 return validated_structure_ptr()->pieces(*
this);
705 if (csab == NULL)
throw Error(
"you requested jet-area related information, but the PseudoJet does not have associated area information.");
712 bool PseudoJet::has_area()
const{
714 if (! has_structure())
return false;
715 return (validated_structure_ptr()->has_area() != 0);
721 double PseudoJet::area()
const{
722 return validated_structure_ptr()->
area(*
this);
729 double PseudoJet::area_error()
const{
730 return validated_structure_ptr()->area_error(*
this);
743 bool PseudoJet::is_pure_ghost()
const{
758 PseudoJet::InexistentUserInfo::InexistentUserInfo() :
Error(
"you attempted to perform a dynamic cast of a PseudoJet's extra info, but the extra info pointer was null")
765 void sort_indices(vector<int> & indices,
766 const vector<double> & values) {
767 IndexedSortHelper index_sort_helper(&values);
768 sort(indices.begin(), indices.end(), index_sort_helper);
775 vector<double> minus_kt2(jets.size());
776 for (
size_t i = 0; i < jets.size(); i++) {minus_kt2[i] = -jets[i].kt2();}
783 vector<double> rapidities(jets.size());
784 for (
size_t i = 0; i < jets.size(); i++) {rapidities[i] = jets[i].rap();}
791 vector<double> energies(jets.size());
792 for (
size_t i = 0; i < jets.size(); i++) {energies[i] = -jets[i].E();}
799 vector<double> pz(jets.size());
800 for (
size_t i = 0; i < jets.size(); i++) {pz[i] = jets[i].pz();}
818 for (
unsigned int i=0; i<pieces.size(); i++)
832 return join(vector<PseudoJet>(1,j1));
839 pieces.push_back(j1);
840 pieces.push_back(j2);
848 pieces.push_back(j1);
849 pieces.push_back(j2);
850 pieces.push_back(j3);
858 pieces.push_back(j1);
859 pieces.push_back(j2);
860 pieces.push_back(j3);
861 pieces.push_back(j4);
868 FASTJET_END_NAMESPACE
PseudoJet PtYPhiM(double pt, double y, double phi, double m)
return a pseudojet with the given pt, y, phi and mass
double theta(const PseudoJet &a, const PseudoJet &b)
returns the angle between a and b
void set_cached_rap_phi(double rap, double phi)
in some cases when setting a 4-momentum, the user/program knows what rapidity and azimuth are associa...
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).
virtual PseudoJet area_4vector() const
return the jet 4-vector area.
const double MaxRap
Used to protect against parton-level events where pt can be zero for some partons, giving rapidity=infinity.
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
vector< PseudoJet > sorted_by_rapidity(const vector< PseudoJet > &jets)
return a vector of jets sorted into increasing rapidity
The structure for a jet made of pieces.
virtual double area(const PseudoJet &) const
return the area associated with the given jet; this base class returns 0.
Contains any information related to the clustering that should be directly accessible to PseudoJet...
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...
virtual bool is_pure_ghost() const
true if this jet is made exclusively of ghosts.
const PseudoJetStructureBase * structure_ptr() const
return a pointer to the structure (of type PseudoJetStructureBase*) associated with this PseudoJet...
base class that sets interface for extensions of ClusterSequence that provide information about the a...
vector< PseudoJet > sorted_by_E(const vector< PseudoJet > &jets)
return a vector of jets sorted into decreasing energy
virtual std::vector< PseudoJet > pieces() const
retrieve the pieces that make up the jet.
base class corresponding to errors that can be thrown by FastJet
double m() const
returns the invariant mass (If m2() is negative then -sqrt(-m2()) is returned, as in CLHEP) ...
int cluster_hist_index() const
return the cluster_hist_index, intended to be used by clustering routines.
an implementation of C++0x shared pointers (or boost's)
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
bool have_same_momentum(const PseudoJet &jeta, const PseudoJet &jetb)
returns true if the momenta of the two input jets are identical
void set_structure_shared_ptr(const SharedPtr< PseudoJetStructureBase > &structure)
set the associated structure
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...
int user_index() const
return the user_index,