PseudoJet.cc

00001 //STARTHEADER
00002 // $Id: PseudoJet.cc 1985 2011-03-09 22:26:40Z soyez $
00003 //
00004 // Copyright (c) 2005-2010, 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 #include "fastjet/Error.hh"
00033 #include "fastjet/PseudoJet.hh"
00034 #include "fastjet/ClusterSequence.hh"
00035 #include "fastjet/ClusterSequenceAreaBase.hh"
00036 #include "fastjet/CompositeJetStructure.hh"
00037 #include<valarray>
00038 #include<iostream>
00039 #include<sstream>
00040 #include<cmath>
00041 #include<algorithm>
00042 #include <cstdarg>
00043 
00044 FASTJET_BEGIN_NAMESPACE      // defined in fastjet/internal/base.hh
00045 
00046 using namespace std;
00047 
00048 
00049 //----------------------------------------------------------------------
00050 // another constructor...
00051 PseudoJet::PseudoJet(const double px, const double py, const double pz, const double E) {
00052   
00053   _E  = E ;
00054   _px = px;
00055   _py = py;
00056   _pz = pz;
00057 
00058   this->_finish_init();
00059 
00060   // some default values for the history and user indices
00061   _reset_indices();
00062 
00063 }
00064 
00065 
00066 //----------------------------------------------------------------------
00067 /// do standard end of initialisation
00068 void PseudoJet::_finish_init () {
00069   _kt2 = this->px()*this->px() + this->py()*this->py();
00070   _phi = pseudojet_invalid_phi;
00071 }
00072 
00073 //----------------------------------------------------------------------
00074 void PseudoJet::_set_rap_phi() const {
00075 
00076   if (_kt2 == 0.0) {
00077     _phi = 0.0; } 
00078   else {
00079     _phi = atan2(this->py(),this->px());
00080   }
00081   if (_phi < 0.0) {_phi += twopi;}
00082   if (_phi >= twopi) {_phi -= twopi;} // can happen if phi=-|eps<1e-15|?
00083   if (this->E() == abs(this->pz()) && _kt2 == 0) {
00084     // Point has infinite rapidity -- convert that into a very large
00085     // number, but in such a way that different 0-pt momenta will have
00086     // different rapidities (so as to lift the degeneracy between
00087     // them) [this can be relevant at parton-level]
00088     double MaxRapHere = MaxRap + abs(this->pz());
00089     if (this->pz() >= 0.0) {_rap = MaxRapHere;} else {_rap = -MaxRapHere;}
00090   } else {
00091     // get the rapidity in a way that's modestly insensitive to roundoff
00092     // error when things pz,E are large (actually the best we can do without
00093     // explicit knowledge of mass)
00094     double effective_m2 = max(0.0,m2()); // force non tachyonic mass
00095     double E_plus_pz    = _E + abs(_pz); // the safer of p+, p-
00096     // p+/p- = (p+ p-) / (p-)^2 = (kt^2+m^2)/(p-)^2
00097     _rap = 0.5*log((_kt2 + effective_m2)/(E_plus_pz*E_plus_pz));
00098     if (_pz > 0) {_rap = - _rap;}
00099   }
00100 
00101 }
00102 
00103 
00104 //----------------------------------------------------------------------
00105 // return a valarray four-momentum
00106 valarray<double> PseudoJet::four_mom() const {
00107   valarray<double> mom(4);
00108   mom[0] = _px;
00109   mom[1] = _py;
00110   mom[2] = _pz;
00111   mom[3] = _E ;
00112   return mom;
00113 }
00114 
00115 //----------------------------------------------------------------------
00116 // Return the component corresponding to the specified index.
00117 // taken from CLHEP
00118 double PseudoJet::operator () (int i) const {
00119   switch(i) {
00120   case X:
00121     return px();
00122   case Y:
00123     return py();
00124   case Z:
00125     return pz();
00126   case T:
00127     return e();
00128   default:
00129     ostringstream err;
00130     err << "PseudoJet subscripting: bad index (" << i << ")";
00131     throw Error(err.str());
00132   }
00133   return 0.;
00134 }  
00135 
00136 //----------------------------------------------------------------------
00137 // return the pseudorapidity
00138 double PseudoJet::pseudorapidity() const {
00139   if (px() == 0.0 && py() ==0.0) return MaxRap;
00140   if (pz() == 0.0) return 0.0;
00141 
00142   double theta = atan(perp()/pz());
00143   if (theta < 0) theta += pi;
00144   return -log(tan(theta/2));
00145 }
00146 
00147 //----------------------------------------------------------------------
00148 // return "sum" of two pseudojets
00149 PseudoJet operator+ (const PseudoJet & jet1, const PseudoJet & jet2) {
00150   //return PseudoJet(jet1.four_mom()+jet2.four_mom());
00151   return PseudoJet(jet1.px()+jet2.px(),
00152                    jet1.py()+jet2.py(),
00153                    jet1.pz()+jet2.pz(),
00154                    jet1.E() +jet2.E()  );
00155 } 
00156 
00157 //----------------------------------------------------------------------
00158 // return difference of two pseudojets
00159 PseudoJet operator- (const PseudoJet & jet1, const PseudoJet & jet2) {
00160   //return PseudoJet(jet1.four_mom()-jet2.four_mom());
00161   return PseudoJet(jet1.px()-jet2.px(),
00162                    jet1.py()-jet2.py(),
00163                    jet1.pz()-jet2.pz(),
00164                    jet1.E() -jet2.E()  );
00165 } 
00166 
00167 //----------------------------------------------------------------------
00168 // return the product, coeff * jet
00169 PseudoJet operator* (double coeff, const PseudoJet & jet) {
00170   //return PseudoJet(coeff*jet.four_mom());
00171   // the following code is hopefully more efficient
00172   PseudoJet coeff_times_jet(jet);
00173   coeff_times_jet *= coeff;
00174   return coeff_times_jet;
00175 } 
00176 
00177 //----------------------------------------------------------------------
00178 // return the product, coeff * jet
00179 PseudoJet operator* (const PseudoJet & jet, double coeff) {
00180   return coeff*jet;
00181 } 
00182 
00183 //----------------------------------------------------------------------
00184 // return the ratio, jet / coeff
00185 PseudoJet operator/ (const PseudoJet & jet, double coeff) {
00186   return (1.0/coeff)*jet;
00187 } 
00188 
00189 //----------------------------------------------------------------------
00190 /// multiply the jet's momentum by the coefficient
00191 void PseudoJet::operator*=(double coeff) {
00192   _px *= coeff;
00193   _py *= coeff;
00194   _pz *= coeff;
00195   _E  *= coeff;
00196   _kt2*= coeff*coeff;
00197   // phi and rap are unchanged
00198 }
00199 
00200 //----------------------------------------------------------------------
00201 /// divide the jet's momentum by the coefficient
00202 void PseudoJet::operator/=(double coeff) {
00203   (*this) *= 1.0/coeff;
00204 }
00205 
00206 
00207 //----------------------------------------------------------------------
00208 /// add the other jet's momentum to this jet
00209 void PseudoJet::operator+=(const PseudoJet & other_jet) {
00210   _px += other_jet._px;
00211   _py += other_jet._py;
00212   _pz += other_jet._pz;
00213   _E  += other_jet._E ;
00214   _finish_init(); // we need to recalculate phi,rap,kt2
00215 }
00216 
00217 
00218 //----------------------------------------------------------------------
00219 /// subtract the other jet's momentum from this jet
00220 void PseudoJet::operator-=(const PseudoJet & other_jet) {
00221   _px -= other_jet._px;
00222   _py -= other_jet._py;
00223   _pz -= other_jet._pz;
00224   _E  -= other_jet._E ;
00225   _finish_init(); // we need to recalculate phi,rap,kt2
00226 }
00227 
00228 
00229 //----------------------------------------------------------------------
00230 /// transform this jet (given in lab) into a jet in the rest
00231 /// frame of prest
00232 //
00233 // NB: code adapted from that in herwig f77 (checked how it worked
00234 // long ago)
00235 PseudoJet & PseudoJet::boost(const PseudoJet & prest) {
00236   
00237   if (prest.px() == 0.0 && prest.py() == 0.0 && prest.pz() == 0.0) 
00238     return *this;
00239 
00240   double m = prest.m();
00241   assert(m != 0);
00242 
00243   double pf4  = (  px()*prest.px() + py()*prest.py()
00244                  + pz()*prest.pz() + E()*prest.E() )/m;
00245   double fn   = (pf4 + E()) / (prest.E() + m);
00246   _px +=  fn*prest.px();
00247   _py +=  fn*prest.py();
00248   _pz +=  fn*prest.pz();
00249   _E = pf4;
00250 
00251   _finish_init(); // we need to recalculate phi,rap,kt2
00252   return *this;
00253 }
00254 
00255 
00256 //----------------------------------------------------------------------
00257 /// transform this jet (given in the rest frame of prest) into a jet
00258 /// in the lab frame;
00259 //
00260 // NB: code adapted from that in herwig f77 (checked how it worked
00261 // long ago)
00262 PseudoJet & PseudoJet::unboost(const PseudoJet & prest) {
00263   
00264   if (prest.px() == 0.0 && prest.py() == 0.0 && prest.pz() == 0.0) 
00265     return *this;
00266 
00267   double m = prest.m();
00268   assert(m != 0);
00269 
00270   double pf4  = ( -px()*prest.px() - py()*prest.py()
00271                  - pz()*prest.pz() + E()*prest.E() )/m;
00272   double fn   = (pf4 + E()) / (prest.E() + m);
00273   _px -=  fn*prest.px();
00274   _py -=  fn*prest.py();
00275   _pz -=  fn*prest.pz();
00276   _E = pf4;
00277 
00278   _finish_init(); // we need to recalculate phi,rap,kt2
00279   return *this;
00280 }
00281 
00282 
00283 //----------------------------------------------------------------------
00284 /// returns true if the momenta of the two input jets are identical
00285 bool have_same_momentum(const PseudoJet & jeta, const PseudoJet & jetb) {
00286   return jeta.px() == jetb.px()
00287     &&   jeta.py() == jetb.py()
00288     &&   jeta.pz() == jetb.pz()
00289     &&   jeta.E()  == jetb.E();
00290 }
00291 
00292 
00293 //----------------------------------------------------------------------
00294 /// return a pseudojet with the given pt, y, phi and mass
00295 PseudoJet PtYPhiM(double pt, double y, double phi, double m) {
00296   double ptm = sqrt(pt*pt+m*m);
00297   return PseudoJet(pt*cos(phi), pt*sin(phi), ptm*sinh(y), ptm*cosh(y));
00298 }
00299 
00300 
00301 //----------------------------------------------------------------------
00302 // return kt-distance between this jet and another one
00303 double PseudoJet::kt_distance(const PseudoJet & other) const {
00304   //double distance = min(this->kt2(), other.kt2());
00305   double distance = min(_kt2, other._kt2);
00306   double dphi = abs(phi() - other.phi());
00307   if (dphi > pi) {dphi = twopi - dphi;}
00308   double drap = rap() - other.rap();
00309   distance = distance * (dphi*dphi + drap*drap);
00310   return distance;
00311 }
00312 
00313 
00314 //----------------------------------------------------------------------
00315 // return squared cylinder (eta-phi) distance between this jet and another one
00316 double PseudoJet::plain_distance(const PseudoJet & other) const {
00317   double dphi = abs(phi() - other.phi());
00318   if (dphi > pi) {dphi = twopi - dphi;}
00319   double drap = rap() - other.rap();
00320   return (dphi*dphi + drap*drap);
00321 }
00322 
00323 //----------------------------------------------------------------------
00324 /// returns other.phi() - this.phi(), i.e. the phi distance to
00325 /// other, constrained to be in range -pi .. pi
00326 double PseudoJet::delta_phi_to(const PseudoJet & other) const {
00327   double dphi = other.phi() - phi();
00328   if (dphi >  pi) dphi -= twopi;
00329   if (dphi < -pi) dphi += twopi;
00330   return dphi;
00331 }
00332 
00333 
00334 string PseudoJet::description() const{
00335   // the "default" case of a PJ which does not belong to any cluster sequence
00336   if (!_structure())
00337     return "standard PseudoJet (with no associated Clustering information)";
00338   
00339   // for all the other cases, the descition comes from the structure
00340   return _structure()->description();
00341 }
00342 
00343 
00344 
00345 //----------------------------------------------------------------------
00346 //
00347 // The following methods access the associated jet structure (if any)
00348 //
00349 //----------------------------------------------------------------------
00350 
00351 
00352 //----------------------------------------------------------------------
00353 // check whether this PseudoJet has an associated parent
00354 // ClusterSequence
00355 bool PseudoJet::has_associated_cluster_sequence() const{
00356   return (_structure()) && (_structure->has_associated_cluster_sequence());
00357 }
00358 
00359 //----------------------------------------------------------------------
00360 // get a (const) pointer to the associated ClusterSequence (NULL if
00361 // inexistent)
00362 const ClusterSequence* PseudoJet::associated_cluster_sequence() const{
00363   if (! has_associated_cluster_sequence()) return NULL;
00364 
00365   return _structure->associated_cluster_sequence();
00366 }
00367 
00368 
00369 //----------------------------------------------------------------------
00370 // If there is a valid cluster sequence associated with this jet,
00371 // returns a pointer to it; otherwise throws an Error.
00372 //
00373 // Open question: should these errors be upgraded to classes of their
00374 // own so that they can be caught? [Maybe, but later]
00375 const ClusterSequence * PseudoJet::validated_cs() const {
00376   return validated_structure_ptr()->validated_cs();
00377 }
00378 
00379 
00380 //----------------------------------------------------------------------
00381 // set the associated structure
00382 void PseudoJet::set_structure_shared_ptr(const SharedPtr<PseudoJetStructureBase> &structure){
00383   _structure = structure;
00384 }
00385 
00386 //----------------------------------------------------------------------
00387 // return true if there is some strusture associated with this PseudoJet
00388 bool PseudoJet::has_structure() const{
00389   return _structure();
00390 }
00391 
00392 //----------------------------------------------------------------------
00393 // return a pointer to the structure (of type
00394 // PseudoJetStructureBase*) associated wioth this PseudoJet.
00395 //
00396 // return NULL if there is no associated structure
00397 const PseudoJetStructureBase* PseudoJet::structure_ptr() const {
00398   if (!_structure()) return NULL;
00399   return _structure();
00400 }
00401   
00402 //----------------------------------------------------------------------
00403 // return a pointer to the structure (of type
00404 // PseudoJetStructureBase*) associated wioth this PseudoJet.
00405 //
00406 // throw an error if there is no associated structure
00407 const PseudoJetStructureBase* PseudoJet::validated_structure_ptr() const {
00408   if (!_structure()) 
00409     throw Error("Trying to access the structure of a PseudoJet which has no associated structure");
00410   return _structure();
00411 }
00412   
00413 //----------------------------------------------------------------------
00414 // return a reference to the shared pointer to the
00415 // PseudoJetStructureBase associated wioth this PseudoJet
00416 const SharedPtr<PseudoJetStructureBase> & PseudoJet::structure_shared_ptr() const {
00417   return _structure;
00418 }
00419 
00420 
00421 //----------------------------------------------------------------------
00422 // check if it has been recombined with another PseudoJet in which
00423 // case, return its partner through the argument. Otherwise,
00424 // 'partner' is set to 0.
00425 //
00426 // false is also returned if this PseudoJet has no associated
00427 // ClusterSequence
00428 bool PseudoJet::has_partner(PseudoJet &partner) const{
00429   return validated_structure_ptr()->has_partner(*this, partner);
00430 }
00431 
00432 //----------------------------------------------------------------------
00433 // check if it has been recombined with another PseudoJet in which
00434 // case, return its child through the argument. Otherwise, 'child'
00435 // is set to 0.
00436 // 
00437 // false is also returned if this PseudoJet has no associated
00438 // ClusterSequence, with the child set to 0
00439 bool PseudoJet::has_child(PseudoJet &child) const{
00440   return validated_structure_ptr()->has_child(*this, child);
00441 }
00442 
00443 //----------------------------------------------------------------------
00444 // check if it is the product of a recombination, in which case
00445 // return the 2 parents through the 'parent1' and 'parent2'
00446 // arguments. Otherwise, set these to 0.
00447 //
00448 // false is also returned if this PseudoJet has no parent
00449 // ClusterSequence
00450 bool PseudoJet::has_parents(PseudoJet &parent1, PseudoJet &parent2) const{
00451   return validated_structure_ptr()->has_parents(*this, parent1, parent2);
00452 }
00453 
00454 //----------------------------------------------------------------------
00455 // check if the current PseudoJet contains the one passed as
00456 // argument
00457 //
00458 // false is also returned if this PseudoJet has no associated
00459 // ClusterSequence.
00460 bool PseudoJet::contains(const PseudoJet &constituent) const{
00461   return validated_structure_ptr()->object_in_jet(constituent, *this);
00462 }
00463 
00464 //----------------------------------------------------------------------
00465 // check if the current PseudoJet is contained the one passed as
00466 // argument
00467 //
00468 // false is also returned if this PseudoJet has no associated
00469 // ClusterSequence
00470 bool PseudoJet::is_inside(const PseudoJet &jet) const{
00471   return validated_structure_ptr()->object_in_jet(*this, jet);
00472 }
00473 
00474 
00475 //----------------------------------------------------------------------
00476 // returns true if the PseudoJet has constituents
00477 bool PseudoJet::has_constituents() const{
00478   return (_structure()) && (_structure->has_constituents());
00479 }
00480 
00481 //----------------------------------------------------------------------
00482 // retrieve the constituents.
00483 vector<PseudoJet> PseudoJet::constituents() const{
00484   return validated_structure_ptr()->constituents(*this);
00485 }
00486 
00487 
00488 //----------------------------------------------------------------------
00489 // returns true if the PseudoJet has support for exclusive subjets
00490 bool PseudoJet::has_exclusive_subjets() const{
00491   return (_structure()) && (_structure->has_exclusive_subjets());
00492 }
00493 
00494 //----------------------------------------------------------------------
00495 // return a vector of all subjets of the current jet (in the sense
00496 // of the exclusive algorithm) that would be obtained when running
00497 // the algorithm with the given dcut. 
00498 //
00499 // Time taken is O(m ln m), where m is the number of subjets that
00500 // are found. If m gets to be of order of the total number of
00501 // constituents in the jet, this could be substantially slower than
00502 // just getting that list of constituents.
00503 //
00504 // an Error is thrown if this PseudoJet has no currently valid
00505 // associated ClusterSequence
00506 std::vector<PseudoJet> PseudoJet::exclusive_subjets (const double & dcut) const {
00507   return validated_structure_ptr()->exclusive_subjets(*this, dcut);
00508 }
00509 
00510 //----------------------------------------------------------------------
00511 // return the size of exclusive_subjets(...); still n ln n with same
00512 // coefficient, but marginally more efficient than manually taking
00513 // exclusive_subjets.size()
00514 //
00515 // an Error is thrown if this PseudoJet has no currently valid
00516 // associated ClusterSequence
00517 int PseudoJet::n_exclusive_subjets(const double & dcut) const {
00518   return validated_structure_ptr()->n_exclusive_subjets(*this, dcut);
00519 }
00520 
00521 //----------------------------------------------------------------------
00522 // return the list of subjets obtained by unclustering the supplied
00523 // jet down to n subjets (or all constituents if there are fewer
00524 // than n).
00525 //
00526 // requires n ln n time
00527 //
00528 // an Error is thrown if this PseudoJet has no currently valid
00529 // associated ClusterSequence
00530 std::vector<PseudoJet> PseudoJet::exclusive_subjets (int nsub) const {
00531   return validated_structure_ptr()->exclusive_subjets(*this, nsub);
00532 }
00533 
00534 //----------------------------------------------------------------------
00535 // return the dij that was present in the merging nsub+1 -> nsub 
00536 // subjets inside this jet.
00537 //
00538 // an Error is thrown if this PseudoJet has no currently valid
00539 // associated ClusterSequence
00540 double PseudoJet::exclusive_subdmerge(int nsub) const {
00541   return validated_structure_ptr()->exclusive_subdmerge(*this, nsub);
00542 }
00543 
00544 //----------------------------------------------------------------------
00545 // return the maximum dij that occurred in the whole event at the
00546 // stage that the nsub+1 -> nsub merge of subjets occurred inside 
00547 // this jet.
00548 //
00549 // an Error is thrown if this PseudoJet has no currently valid
00550 // associated ClusterSequence
00551 double PseudoJet::exclusive_subdmerge_max(int nsub) const {
00552   return validated_structure_ptr()->exclusive_subdmerge_max(*this, nsub);
00553 }
00554 
00555 
00556 // returns true if a jet has pieces
00557 //
00558 // By default a single particle or a jet coming from a
00559 // ClusterSequence have no pieces and this methos will return false.
00560 bool PseudoJet::has_pieces() const{
00561   return ((_structure()) && (_structure->has_pieces()));
00562 }
00563 
00564 // retrieve the pieces that make up the jet. 
00565 //
00566 // By default a jet does not have pieces.
00567 // If the underlying interface supports "pieces" retrieve the
00568 // pieces from there.
00569 std::vector<PseudoJet> PseudoJet::pieces() const{
00570   if (!has_pieces())
00571     throw Error("Trying to retrieve the pieces of a PseudoJet that has no support for pieces.");
00572 
00573   return _structure->pieces(*this);
00574 }
00575 
00576 
00577 //----------------------------------------------------------------------
00578 // the following ones require a computation of the area in the
00579 // associated ClusterSequence (See ClusterSequenceAreaBase for details)
00580 //----------------------------------------------------------------------
00581 
00582 //----------------------------------------------------------------------
00583 // if possible, return a valid ClusterSequenceAreaBase pointer; otherwise
00584 // throw an error
00585 const ClusterSequenceAreaBase * PseudoJet::validated_csab() const {
00586   const ClusterSequenceAreaBase *csab = dynamic_cast<const ClusterSequenceAreaBase*>(validated_cs());
00587   if (csab == NULL) throw Error("you requested jet-area related information, but the PseudoJet does not have associated area information.");
00588   return csab;
00589 }
00590 
00591 
00592 //----------------------------------------------------------------------
00593 // check if it has a defined area
00594 bool PseudoJet::has_area() const{
00595   if (! has_associated_cluster_sequence()) return false;
00596   return (validated_structure_ptr()->has_area() != 0);
00597 }
00598 
00599 //----------------------------------------------------------------------
00600 // return the jet (scalar) area.
00601 // throw an Error if there is no support for area in the associated CS
00602 double PseudoJet::area() const{
00603   return validated_structure_ptr()->area(*this);
00604 }
00605 
00606 //----------------------------------------------------------------------
00607 // return the error (uncertainty) associated with the determination
00608 // of the area of this jet.
00609 // throws an Error if there is no support for area in the associated CS
00610 double PseudoJet::area_error() const{
00611   return validated_structure_ptr()->area_error(*this);
00612 }
00613 
00614 //----------------------------------------------------------------------
00615 // return the jet 4-vector area
00616 // throws an Error if there is no support for area in the associated CS
00617 PseudoJet PseudoJet::area_4vector() const{
00618   return validated_structure_ptr()->area_4vector(*this);
00619 }
00620 
00621 //----------------------------------------------------------------------
00622 // true if this jet is made exclusively of ghosts
00623 // throws an Error if there is no support for area in the associated CS
00624 bool PseudoJet::is_pure_ghost() const{
00625   return validated_structure_ptr()->is_pure_ghost(*this);
00626 }
00627 
00628 
00629 //----------------------------------------------------------------------
00630 //
00631 // end of the methods accessing the information in the associated
00632 // Cluster Sequence
00633 //
00634 //----------------------------------------------------------------------
00635 
00636 //----------------------------------------------------------------------
00637 /// provide a meaningful error message for InexistentUserInfo
00638 PseudoJet::InexistentUserInfo::InexistentUserInfo() : Error("you attempted to perform a dynamic cast of a PseudoJet's extra info, but the extra info pointer was null")
00639 {}
00640 
00641 
00642 //----------------------------------------------------------------------
00643 // sort the indices so that values[indices[0..n-1]] is sorted
00644 // into increasing order 
00645 void sort_indices(vector<int> & indices, 
00646                          const vector<double> & values) {
00647   IndexedSortHelper index_sort_helper(&values);
00648   sort(indices.begin(), indices.end(), index_sort_helper);
00649 }
00650 
00651 
00652 
00653 //----------------------------------------------------------------------
00654 /// given a vector of values with a one-to-one correspondence with the
00655 /// vector of objects, sort objects into an order such that the
00656 /// associated values would be in increasing order
00657 template<class T> vector<T>  objects_sorted_by_values(
00658                        const vector<T> & objects, 
00659                        const vector<double> & values) {
00660 
00661   assert(objects.size() == values.size());
00662 
00663   // get a vector of indices
00664   vector<int> indices(values.size());
00665   for (size_t i = 0; i < indices.size(); i++) {indices[i] = i;}
00666   
00667   // sort the indices
00668   sort_indices(indices, values);
00669   
00670   // copy the objects 
00671   vector<T> objects_sorted(objects.size());
00672   
00673   // place the objects in the correct order
00674   for (size_t i = 0; i < indices.size(); i++) {
00675     objects_sorted[i] = objects[indices[i]];
00676   }
00677 
00678   return objects_sorted;
00679 }
00680 
00681 //----------------------------------------------------------------------
00682 /// return a vector of jets sorted into decreasing kt2
00683 vector<PseudoJet> sorted_by_pt(const vector<PseudoJet> & jets) {
00684   vector<double> minus_kt2(jets.size());
00685   for (size_t i = 0; i < jets.size(); i++) {minus_kt2[i] = -jets[i].kt2();}
00686   return objects_sorted_by_values(jets, minus_kt2);
00687 }
00688 
00689 //----------------------------------------------------------------------
00690 /// return a vector of jets sorted into increasing rapidity
00691 vector<PseudoJet> sorted_by_rapidity(const vector<PseudoJet> & jets) {
00692   vector<double> rapidities(jets.size());
00693   for (size_t i = 0; i < jets.size(); i++) {rapidities[i] = jets[i].rap();}
00694   return objects_sorted_by_values(jets, rapidities);
00695 }
00696 
00697 //----------------------------------------------------------------------
00698 /// return a vector of jets sorted into decreasing energy
00699 vector<PseudoJet> sorted_by_E(const vector<PseudoJet> & jets) {
00700   vector<double> energies(jets.size());
00701   for (size_t i = 0; i < jets.size(); i++) {energies[i] = -jets[i].E();}
00702   return objects_sorted_by_values(jets, energies);
00703 }
00704 
00705 //----------------------------------------------------------------------
00706 /// return a vector of jets sorted into increasing pz
00707 vector<PseudoJet> sorted_by_pz(const vector<PseudoJet> & jets) {
00708   vector<double> pz(jets.size());
00709   for (size_t i = 0; i < jets.size(); i++) {pz[i] = jets[i].pz();}
00710   return objects_sorted_by_values(jets, pz);
00711 }
00712 
00713 
00714 //-------------------------------------------------------------------------------
00715 // helper functions to build a jet made of pieces
00716 //-------------------------------------------------------------------------------
00717 
00718 // build a "CompositeJet" from the vector of its pieces
00719 //
00720 // In this case, E-scheme recombination is assumed to compute the
00721 // total momentum
00722 PseudoJet join(const vector<PseudoJet> & pieces){
00723   PseudoJet result(0.0,0.0,0.0,0.0);
00724   for (unsigned int i=0; i<pieces.size(); i++){
00725     const PseudoJet it = pieces[i];
00726     result += it;
00727   }
00728 
00729   CompositeJetStructure *cj_struct = new CompositeJetStructure(pieces);
00730   result.set_structure_shared_ptr(SharedPtr<PseudoJetStructureBase>(cj_struct));
00731 
00732   return result;
00733 }
00734 
00735 // build a "CompositeJet" from a single PseudoJet
00736 PseudoJet join(const PseudoJet & j1){
00737   return join(vector<PseudoJet>(1,j1));
00738 }
00739 
00740 // build a "CompositeJet" from two PseudoJet
00741 PseudoJet join(const PseudoJet & j1, const PseudoJet & j2){
00742   vector<PseudoJet> pieces;
00743   pieces.push_back(j1);
00744   pieces.push_back(j2);
00745   return join(pieces);
00746 }
00747 
00748 // build a "CompositeJet" from 3 PseudoJet
00749 PseudoJet join(const PseudoJet & j1, const PseudoJet & j2, const PseudoJet & j3){
00750   vector<PseudoJet> pieces;
00751   pieces.push_back(j1);
00752   pieces.push_back(j2);
00753   pieces.push_back(j3);
00754   return join(pieces);
00755 }
00756 
00757 // build a "CompositeJet" from 4 PseudoJet
00758 PseudoJet join(const PseudoJet & j1, const PseudoJet & j2, const PseudoJet & j3, const PseudoJet & j4){
00759   vector<PseudoJet> pieces;
00760   pieces.push_back(j1);
00761   pieces.push_back(j2);
00762   pieces.push_back(j3);
00763   pieces.push_back(j4);
00764   return join(pieces);
00765 }
00766 
00767 
00768 
00769 FASTJET_END_NAMESPACE
00770