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()  );
 
  169 PseudoJet operator- (
const PseudoJet & jet1, 
const PseudoJet & jet2) {
 
  171   return PseudoJet(jet1.px()-jet2.px(),
 
  174                    jet1.E() -jet2.E()  );
 
  179 PseudoJet 
operator* (
double coeff, 
const PseudoJet & jet) {
 
  183   jet._ensure_valid_rap_phi(); 
 
  186   PseudoJet coeff_times_jet(jet);
 
  187   coeff_times_jet *= coeff;
 
  188   return coeff_times_jet;
 
  193 PseudoJet 
operator* (
const PseudoJet & jet, 
double coeff) {
 
  199 PseudoJet operator/ (
const PseudoJet & jet, 
double coeff) {
 
  200   return (1.0/coeff)*jet;
 
  205 void PseudoJet::operator*=(
double coeff) {
 
  213   _ensure_valid_rap_phi(); 
 
  224 void PseudoJet::operator/=(
double coeff) {
 
  225   (*this) *= 1.0/coeff;
 
  231 void PseudoJet::operator+=(
const PseudoJet & other_jet) {
 
  232   _px += other_jet._px;
 
  233   _py += other_jet._py;
 
  234   _pz += other_jet._pz;
 
  242 void PseudoJet::operator-=(
const PseudoJet & other_jet) {
 
  243   _px -= other_jet._px;
 
  244   _py -= other_jet._py;
 
  245   _pz -= other_jet._pz;
 
  252   if (a.px() != b.px()) 
return false;
 
  253   if (a.py() != b.py()) 
return false;
 
  254   if (a.pz() != b.pz()) 
return false;
 
  255   if (a.E () != b.E ()) 
return false;
 
  269     throw Error(
"comparing a PseudoJet with a non-zero constant (double) is not allowed.");
 
  270   return (jet.px() == 0 && jet.py() == 0 && 
 
  271           jet.pz() == 0 && jet.E() == 0);
 
  284   if (prest.px() == 0.0 && prest.py() == 0.0 && prest.pz() == 0.0) 
 
  287   double m_local = prest.
m();
 
  288   assert(m_local != 0);
 
  290   double pf4  = (  px()*prest.px() + py()*prest.py()
 
  291                  + pz()*prest.pz() + E()*prest.E() )/m_local;
 
  292   double fn   = (pf4 + E()) / (prest.E() + m_local);
 
  293   _px +=  fn*prest.px();
 
  294   _py +=  fn*prest.py();
 
  295   _pz +=  fn*prest.pz();
 
  311   if (prest.px() == 0.0 && prest.py() == 0.0 && prest.pz() == 0.0) 
 
  314   double m_local = prest.
m();
 
  315   assert(m_local != 0);
 
  317   double pf4  = ( -px()*prest.px() - py()*prest.py()
 
  318                  - pz()*prest.pz() + E()*prest.E() )/m_local;
 
  319   double fn   = (pf4 + E()) / (prest.E() + m_local);
 
  320   _px -=  fn*prest.px();
 
  321   _py -=  fn*prest.py();
 
  322   _pz -=  fn*prest.pz();
 
  333   return jeta.px() == jetb.px()
 
  334     &&   jeta.py() == jetb.py()
 
  335     &&   jeta.pz() == jetb.pz()
 
  336     &&   jeta.E()  == jetb.E();
 
  340 void PseudoJet::set_cached_rap_phi(
double rap_in, 
double phi_in) {
 
  341   _rap = rap_in; _phi = phi_in;
 
  342   if (_phi >= twopi) _phi -= twopi;
 
  343   if (_phi < 0)      _phi += twopi;
 
  347 void PseudoJet::reset_momentum_PtYPhiM(
double pt_in, 
double y_in, 
double phi_in, 
double m_in) {
 
  348   assert(phi_in < 2*twopi && phi_in > -twopi);
 
  349   double ptm = (m_in == 0) ? pt_in : sqrt(pt_in*pt_in+m_in*m_in);
 
  350   double exprap = exp(y_in);
 
  351   double pminus = ptm/exprap;
 
  352   double pplus  = ptm*exprap;
 
  353   double px_local = pt_in*cos(phi_in);
 
  354   double py_local = pt_in*sin(phi_in);
 
  355   reset_momentum(px_local,py_local,0.5*(pplus-pminus),0.5*(pplus+pminus));
 
  356   set_cached_rap_phi(y_in,phi_in);
 
  362   assert(phi < 2*twopi && phi > -twopi);
 
  363   double ptm = (m == 0) ? pt : sqrt(pt*pt+m*m);
 
  364   double exprap = exp(y);
 
  365   double pminus = ptm/exprap;
 
  366   double pplus  = ptm*exprap;
 
  367   double px = pt*cos(phi);
 
  368   double py = pt*sin(phi);
 
  369   PseudoJet mom(px,py,0.5*(pplus-pminus),0.5*(pplus+pminus));
 
  378 double PseudoJet::kt_distance(
const PseudoJet & other)
 const {
 
  380   double distance = min(_kt2, other._kt2);
 
  381   double dphi = abs(phi() - other.
phi());
 
  382   if (dphi > pi) {dphi = twopi - dphi;}
 
  383   double drap = rap() - other.
rap();
 
  384   distance = distance * (dphi*dphi + drap*drap);
 
  391 double PseudoJet::plain_distance(
const PseudoJet & other)
 const {
 
  392   double dphi = abs(phi() - other.
phi());
 
  393   if (dphi > pi) {dphi = twopi - dphi;}
 
  394   double drap = rap() - other.
rap();
 
  395   return (dphi*dphi + drap*drap);
 
  401 double PseudoJet::delta_phi_to(
const PseudoJet & other)
 const {
 
  402   double dphi = other.
phi() - phi();
 
  403   if (dphi >  pi) dphi -= twopi;
 
  404   if (dphi < -pi) dphi += twopi;
 
  409 string PseudoJet::description()
 const{
 
  412     return "standard PseudoJet (with no associated clustering information)";
 
  415   return _structure()->description();
 
  430 bool PseudoJet::has_associated_cluster_sequence()
 const{
 
  431   return (_structure()) && (_structure->has_associated_cluster_sequence());
 
  438   if (! has_associated_cluster_sequence()) 
return NULL;
 
  440   return _structure->associated_cluster_sequence();
 
  447 bool PseudoJet::has_valid_cluster_sequence()
 const{
 
  448   return (_structure()) && (_structure->has_valid_cluster_sequence());
 
  458   return validated_structure_ptr()->validated_cs();
 
  465   _structure = structure_in;
 
  470 bool PseudoJet::has_structure()
 const{
 
  480   if (!_structure()) 
return NULL;
 
  495   if (!_structure()) 
return NULL;
 
  506     throw Error(
"Trying to access the structure of a PseudoJet which has no associated structure");
 
  526   return validated_structure_ptr()->has_partner(*
this, partner);
 
  537   return validated_structure_ptr()->has_child(*
this, child);
 
  548   return validated_structure_ptr()->has_parents(*
this, parent1, parent2);
 
  557 bool PseudoJet::contains(
const PseudoJet &constituent)
 const{
 
  558   return validated_structure_ptr()->object_in_jet(constituent, *
this);
 
  568   return validated_structure_ptr()->object_in_jet(*
this, jet);
 
  574 bool PseudoJet::has_constituents()
 const{
 
  575   return (_structure()) && (_structure->has_constituents());
 
  580 vector<PseudoJet> PseudoJet::constituents()
 const{
 
  581   return validated_structure_ptr()->constituents(*
this);
 
  587 bool PseudoJet::has_exclusive_subjets()
 const{
 
  588   return (_structure()) && (_structure->has_exclusive_subjets());
 
  603 std::vector<PseudoJet> PseudoJet::exclusive_subjets (
const double dcut)
 const {
 
  604   return validated_structure_ptr()->exclusive_subjets(*
this, dcut);
 
  614 int PseudoJet::n_exclusive_subjets(
const double dcut)
 const {
 
  615   return validated_structure_ptr()->n_exclusive_subjets(*
this, dcut);
 
  627 std::vector<PseudoJet> PseudoJet::exclusive_subjets_up_to (
int nsub)
 const {
 
  628   return validated_structure_ptr()->exclusive_subjets_up_to(*
this, nsub);
 
  634 std::vector<PseudoJet> PseudoJet::exclusive_subjets (
int nsub)
 const {
 
  635   vector<PseudoJet> subjets = exclusive_subjets_up_to(nsub);
 
  636   if (
int(subjets.size()) < nsub) {
 
  638     err << 
"Requested " << nsub << 
" exclusive subjets, but there were only "  
  639         << subjets.size() << 
" particles in the jet";
 
  640     throw Error(err.str());
 
  651 double PseudoJet::exclusive_subdmerge(
int nsub)
 const {
 
  652   return validated_structure_ptr()->exclusive_subdmerge(*
this, nsub);
 
  662 double PseudoJet::exclusive_subdmerge_max(
int nsub)
 const {
 
  663   return validated_structure_ptr()->exclusive_subdmerge_max(*
this, nsub);
 
  671 bool PseudoJet::has_pieces()
 const{
 
  672   return ((_structure()) && (_structure->has_pieces(*
this)));
 
  680 std::vector<PseudoJet> PseudoJet::pieces()
 const{
 
  681   return validated_structure_ptr()->pieces(*
this);
 
  701   if (csab == NULL) 
throw Error(
"you requested jet-area related information, but the PseudoJet does not have associated area information.");
 
  708 bool PseudoJet::has_area()
 const{
 
  710   if (! has_structure()) 
return false;
 
  711   return (validated_structure_ptr()->has_area() != 0);
 
  717 double PseudoJet::area()
 const{
 
  718   return validated_structure_ptr()->
area(*
this);
 
  725 double PseudoJet::area_error()
 const{
 
  726   return validated_structure_ptr()->area_error(*
this);
 
  739 bool PseudoJet::is_pure_ghost()
 const{
 
  754 PseudoJet::InexistentUserInfo::InexistentUserInfo() : 
Error(
"you attempted to perform a dynamic cast of a PseudoJet's extra info, but the extra info pointer was null")
 
  761 void sort_indices(vector<int> & indices, 
 
  762                          const vector<double> & values) {
 
  763   IndexedSortHelper index_sort_helper(&values);
 
  764   sort(indices.begin(), indices.end(), index_sort_helper);
 
  774                        const vector<T> & objects, 
 
  775                        const vector<double> & values) {
 
  777   assert(objects.size() == values.size());
 
  780   vector<int> indices(values.size());
 
  781   for (
size_t i = 0; i < indices.size(); i++) {indices[i] = i;}
 
  784   sort_indices(indices, values);
 
  787   vector<T> objects_sorted(objects.size());
 
  790   for (
size_t i = 0; i < indices.size(); i++) {
 
  791     objects_sorted[i] = objects[indices[i]];
 
  794   return objects_sorted;
 
  800   vector<double> minus_kt2(jets.size());
 
  801   for (
size_t i = 0; i < jets.size(); i++) {minus_kt2[i] = -jets[i].kt2();}
 
  808   vector<double> rapidities(jets.size());
 
  809   for (
size_t i = 0; i < jets.size(); i++) {rapidities[i] = jets[i].rap();}
 
  816   vector<double> energies(jets.size());
 
  817   for (
size_t i = 0; i < jets.size(); i++) {energies[i] = -jets[i].E();}
 
  824   vector<double> pz(jets.size());
 
  825   for (
size_t i = 0; i < jets.size(); i++) {pz[i] = jets[i].pz();}
 
  843   for (
unsigned int i=0; i<pieces.size(); i++)
 
  848   CompositeJetStructure *cj_struct = 
new CompositeJetStructure(pieces);
 
  850   result.set_structure_shared_ptr(SharedPtr<PseudoJetStructureBase>(cj_struct));
 
  857   return join(vector<PseudoJet>(1,j1));
 
  864   pieces.push_back(j1);
 
  865   pieces.push_back(j2);
 
  873   pieces.push_back(j1);
 
  874   pieces.push_back(j2);
 
  875   pieces.push_back(j3);
 
  883   pieces.push_back(j1);
 
  884   pieces.push_back(j2);
 
  885   pieces.push_back(j3);
 
  886   pieces.push_back(j4);
 
  893 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 rap() const 
returns the rapidity or some large value when the rapidity is infinite 
 
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 
 
const UserInfoBase * user_info_ptr() const 
retrieve a pointer to the (const) user information 
 
const double MaxRap
Used to protect against parton-level events where pt can be zero for some partons, giving rapidity=infinity. 
 
int cluster_hist_index() const 
return the cluster_hist_index, intended to be used by clustering routines. 
 
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 
 
Contains any information related to the clustering that should be directly accessible to PseudoJet...
 
virtual std::vector< PseudoJet > pieces() const 
retrieve the pieces that make up the jet. 
 
vector< T > objects_sorted_by_values(const vector< T > &objects, const 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 
 
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 double area(const PseudoJet &) const 
return the area associated with the given jet; this base class returns 0. 
 
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 
 
base class corresponding to errors that can be thrown by FastJet 
 
an implementation of C++0x shared pointers (or boost's) 
 
virtual bool is_pure_ghost() const 
true if this jet is made exclusively of ghosts. 
 
bool have_same_momentum(const PseudoJet &jeta, const PseudoJet &jetb)
returns true if the momenta of the two input jets are identical 
 
double phi() const 
returns phi (in the range 0..2pi) 
 
int user_index() const 
return the user_index, 
 
PseudoJet()
default constructor, which as of FJ3.0 provides an object for which all operations are now valid and ...
 
double m() const 
returns the invariant mass (If m2() is negative then -sqrt(-m2()) is returned, as in CLHEP) ...
 
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 PseudoJetStructureBase * structure_ptr() const 
return a pointer to the structure (of type PseudoJetStructureBase*) associated with this PseudoJet...
 
virtual PseudoJet area_4vector() const 
return the jet 4-vector area.