FastJet 3.0beta1
|
00001 //STARTHEADER 00002 // $Id: PseudoJet.cc 2330 2011-06-30 15:15:09Z salam $ 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 bool operator==(const PseudoJet & a, const PseudoJet & b) { 00230 if (a.px() != b.px()) return false; 00231 if (a.py() != b.py()) return false; 00232 if (a.pz() != b.pz()) return false; 00233 if (a.E () != b.E ()) return false; 00234 00235 if (a.user_index() != b.user_index()) return false; 00236 if (a.cluster_hist_index() != b.cluster_hist_index()) return false; 00237 if (a.user_info_ptr() != b.user_info_ptr()) return false; 00238 if (a.structure_ptr() != b.structure_ptr()) return false; 00239 00240 return true; 00241 } 00242 00243 //---------------------------------------------------------------------- 00244 // check if the jet has zero momentum 00245 bool operator==(const PseudoJet & jet, const double val) { 00246 if (val != 0) 00247 throw Error("comparing a PseudoJet with a non-zero constant (double) is not allowed."); 00248 return (jet.px() == 0 && jet.py() == 0 && 00249 jet.pz() == 0 && jet.E() == 0); 00250 } 00251 00252 00253 00254 //---------------------------------------------------------------------- 00255 /// transform this jet (given in lab) into a jet in the rest 00256 /// frame of prest 00257 // 00258 // NB: code adapted from that in herwig f77 (checked how it worked 00259 // long ago) 00260 PseudoJet & PseudoJet::boost(const PseudoJet & prest) { 00261 00262 if (prest.px() == 0.0 && prest.py() == 0.0 && prest.pz() == 0.0) 00263 return *this; 00264 00265 double m = prest.m(); 00266 assert(m != 0); 00267 00268 double pf4 = ( px()*prest.px() + py()*prest.py() 00269 + pz()*prest.pz() + E()*prest.E() )/m; 00270 double fn = (pf4 + E()) / (prest.E() + m); 00271 _px += fn*prest.px(); 00272 _py += fn*prest.py(); 00273 _pz += fn*prest.pz(); 00274 _E = pf4; 00275 00276 _finish_init(); // we need to recalculate phi,rap,kt2 00277 return *this; 00278 } 00279 00280 00281 //---------------------------------------------------------------------- 00282 /// transform this jet (given in the rest frame of prest) into a jet 00283 /// in the lab frame; 00284 // 00285 // NB: code adapted from that in herwig f77 (checked how it worked 00286 // long ago) 00287 PseudoJet & PseudoJet::unboost(const PseudoJet & prest) { 00288 00289 if (prest.px() == 0.0 && prest.py() == 0.0 && prest.pz() == 0.0) 00290 return *this; 00291 00292 double m = prest.m(); 00293 assert(m != 0); 00294 00295 double pf4 = ( -px()*prest.px() - py()*prest.py() 00296 - pz()*prest.pz() + E()*prest.E() )/m; 00297 double fn = (pf4 + E()) / (prest.E() + m); 00298 _px -= fn*prest.px(); 00299 _py -= fn*prest.py(); 00300 _pz -= fn*prest.pz(); 00301 _E = pf4; 00302 00303 _finish_init(); // we need to recalculate phi,rap,kt2 00304 return *this; 00305 } 00306 00307 00308 //---------------------------------------------------------------------- 00309 /// returns true if the momenta of the two input jets are identical 00310 bool have_same_momentum(const PseudoJet & jeta, const PseudoJet & jetb) { 00311 return jeta.px() == jetb.px() 00312 && jeta.py() == jetb.py() 00313 && jeta.pz() == jetb.pz() 00314 && jeta.E() == jetb.E(); 00315 } 00316 00317 //---------------------------------------------------------------------- 00318 void PseudoJet::set_cached_rap_phi(double rap, double phi) { 00319 _rap = rap; _phi = phi; 00320 if (_phi >= twopi) _phi -= twopi; 00321 if (_phi < 0) _phi += twopi; 00322 } 00323 00324 //---------------------------------------------------------------------- 00325 void PseudoJet::reset_momentum_PtYPhiM(double pt, double y, double phi, double m) { 00326 assert(phi < 2*twopi && phi > -twopi); 00327 double ptm = (m == 0) ? pt : sqrt(pt*pt+m*m); 00328 double exprap = exp(y); 00329 double pminus = ptm/exprap; 00330 double pplus = ptm*exprap; 00331 double px = pt*cos(phi); 00332 double py = pt*sin(phi); 00333 reset_momentum(px,py,0.5*(pplus-pminus),0.5*(pplus+pminus)); 00334 set_cached_rap_phi(y,phi); 00335 } 00336 00337 //---------------------------------------------------------------------- 00338 /// return a pseudojet with the given pt, y, phi and mass 00339 PseudoJet PtYPhiM(double pt, double y, double phi, double m) { 00340 assert(phi < 2*twopi && phi > -twopi); 00341 double ptm = (m == 0) ? pt : sqrt(pt*pt+m*m); 00342 double exprap = exp(y); 00343 double pminus = ptm/exprap; 00344 double pplus = ptm*exprap; 00345 double px = pt*cos(phi); 00346 double py = pt*sin(phi); 00347 PseudoJet mom(px,py,0.5*(pplus-pminus),0.5*(pplus+pminus)); 00348 mom.set_cached_rap_phi(y,phi); 00349 return mom; 00350 //return PseudoJet(pt*cos(phi), pt*sin(phi), ptm*sinh(y), ptm*cosh(y)); 00351 } 00352 00353 00354 //---------------------------------------------------------------------- 00355 // return kt-distance between this jet and another one 00356 double PseudoJet::kt_distance(const PseudoJet & other) const { 00357 //double distance = min(this->kt2(), other.kt2()); 00358 double distance = min(_kt2, other._kt2); 00359 double dphi = abs(phi() - other.phi()); 00360 if (dphi > pi) {dphi = twopi - dphi;} 00361 double drap = rap() - other.rap(); 00362 distance = distance * (dphi*dphi + drap*drap); 00363 return distance; 00364 } 00365 00366 00367 //---------------------------------------------------------------------- 00368 // return squared cylinder (eta-phi) distance between this jet and another one 00369 double PseudoJet::plain_distance(const PseudoJet & other) const { 00370 double dphi = abs(phi() - other.phi()); 00371 if (dphi > pi) {dphi = twopi - dphi;} 00372 double drap = rap() - other.rap(); 00373 return (dphi*dphi + drap*drap); 00374 } 00375 00376 //---------------------------------------------------------------------- 00377 /// returns other.phi() - this.phi(), i.e. the phi distance to 00378 /// other, constrained to be in range -pi .. pi 00379 double PseudoJet::delta_phi_to(const PseudoJet & other) const { 00380 double dphi = other.phi() - phi(); 00381 if (dphi > pi) dphi -= twopi; 00382 if (dphi < -pi) dphi += twopi; 00383 return dphi; 00384 } 00385 00386 00387 string PseudoJet::description() const{ 00388 // the "default" case of a PJ which does not belong to any cluster sequence 00389 if (!_structure()) 00390 return "standard PseudoJet (with no associated Clustering information)"; 00391 00392 // for all the other cases, the descition comes from the structure 00393 return _structure()->description(); 00394 } 00395 00396 00397 00398 //---------------------------------------------------------------------- 00399 // 00400 // The following methods access the associated jet structure (if any) 00401 // 00402 //---------------------------------------------------------------------- 00403 00404 00405 //---------------------------------------------------------------------- 00406 // check whether this PseudoJet has an associated parent 00407 // ClusterSequence 00408 bool PseudoJet::has_associated_cluster_sequence() const{ 00409 return (_structure()) && (_structure->has_associated_cluster_sequence()); 00410 } 00411 00412 //---------------------------------------------------------------------- 00413 // get a (const) pointer to the associated ClusterSequence (NULL if 00414 // inexistent) 00415 const ClusterSequence* PseudoJet::associated_cluster_sequence() const{ 00416 if (! has_associated_cluster_sequence()) return NULL; 00417 00418 return _structure->associated_cluster_sequence(); 00419 } 00420 00421 00422 //---------------------------------------------------------------------- 00423 // check whether this PseudoJet has an associated parent 00424 // ClusterSequence that is still valid 00425 bool PseudoJet::has_valid_cluster_sequence() const{ 00426 return (_structure()) && (_structure->has_valid_cluster_sequence()); 00427 } 00428 00429 //---------------------------------------------------------------------- 00430 // If there is a valid cluster sequence associated with this jet, 00431 // returns a pointer to it; otherwise throws an Error. 00432 // 00433 // Open question: should these errors be upgraded to classes of their 00434 // own so that they can be caught? [Maybe, but later] 00435 const ClusterSequence * PseudoJet::validated_cs() const { 00436 return validated_structure_ptr()->validated_cs(); 00437 } 00438 00439 00440 //---------------------------------------------------------------------- 00441 // set the associated structure 00442 void PseudoJet::set_structure_shared_ptr(const SharedPtr<PseudoJetStructureBase> &structure){ 00443 _structure = structure; 00444 } 00445 00446 //---------------------------------------------------------------------- 00447 // return true if there is some strusture associated with this PseudoJet 00448 bool PseudoJet::has_structure() const{ 00449 return _structure(); 00450 } 00451 00452 //---------------------------------------------------------------------- 00453 // return a pointer to the structure (of type 00454 // PseudoJetStructureBase*) associated with this PseudoJet. 00455 // 00456 // return NULL if there is no associated structure 00457 const PseudoJetStructureBase* PseudoJet::structure_ptr() const { 00458 if (!_structure()) return NULL; 00459 return _structure(); 00460 } 00461 00462 //---------------------------------------------------------------------- 00463 // return a non-const pointer to the structure (of type 00464 // PseudoJetStructureBase*) associated with this PseudoJet. 00465 // 00466 // return NULL if there is no associated structure 00467 // 00468 // Only use this if you know what you are doing. In any case, 00469 // prefer the 'structure_ptr()' (the const version) to this method, 00470 // unless you really need a write access to the PseudoJet's 00471 // underlying structure. 00472 PseudoJetStructureBase* PseudoJet::structure_non_const_ptr(){ 00473 if (!_structure()) return NULL; 00474 return _structure(); 00475 } 00476 00477 //---------------------------------------------------------------------- 00478 // return a pointer to the structure (of type 00479 // PseudoJetStructureBase*) associated with this PseudoJet. 00480 // 00481 // throw an error if there is no associated structure 00482 const PseudoJetStructureBase* PseudoJet::validated_structure_ptr() const { 00483 if (!_structure()) 00484 throw Error("Trying to access the structure of a PseudoJet which has no associated structure"); 00485 return _structure(); 00486 } 00487 00488 //---------------------------------------------------------------------- 00489 // return a reference to the shared pointer to the 00490 // PseudoJetStructureBase associated with this PseudoJet 00491 const SharedPtr<PseudoJetStructureBase> & PseudoJet::structure_shared_ptr() const { 00492 return _structure; 00493 } 00494 00495 00496 //---------------------------------------------------------------------- 00497 // check if it has been recombined with another PseudoJet in which 00498 // case, return its partner through the argument. Otherwise, 00499 // 'partner' is set to 0. 00500 // 00501 // false is also returned if this PseudoJet has no associated 00502 // ClusterSequence 00503 bool PseudoJet::has_partner(PseudoJet &partner) const{ 00504 return validated_structure_ptr()->has_partner(*this, partner); 00505 } 00506 00507 //---------------------------------------------------------------------- 00508 // check if it has been recombined with another PseudoJet in which 00509 // case, return its child through the argument. Otherwise, 'child' 00510 // is set to 0. 00511 // 00512 // false is also returned if this PseudoJet has no associated 00513 // ClusterSequence, with the child set to 0 00514 bool PseudoJet::has_child(PseudoJet &child) const{ 00515 return validated_structure_ptr()->has_child(*this, child); 00516 } 00517 00518 //---------------------------------------------------------------------- 00519 // check if it is the product of a recombination, in which case 00520 // return the 2 parents through the 'parent1' and 'parent2' 00521 // arguments. Otherwise, set these to 0. 00522 // 00523 // false is also returned if this PseudoJet has no parent 00524 // ClusterSequence 00525 bool PseudoJet::has_parents(PseudoJet &parent1, PseudoJet &parent2) const{ 00526 return validated_structure_ptr()->has_parents(*this, parent1, parent2); 00527 } 00528 00529 //---------------------------------------------------------------------- 00530 // check if the current PseudoJet contains the one passed as 00531 // argument 00532 // 00533 // false is also returned if this PseudoJet has no associated 00534 // ClusterSequence. 00535 bool PseudoJet::contains(const PseudoJet &constituent) const{ 00536 return validated_structure_ptr()->object_in_jet(constituent, *this); 00537 } 00538 00539 //---------------------------------------------------------------------- 00540 // check if the current PseudoJet is contained the one passed as 00541 // argument 00542 // 00543 // false is also returned if this PseudoJet has no associated 00544 // ClusterSequence 00545 bool PseudoJet::is_inside(const PseudoJet &jet) const{ 00546 return validated_structure_ptr()->object_in_jet(*this, jet); 00547 } 00548 00549 00550 //---------------------------------------------------------------------- 00551 // returns true if the PseudoJet has constituents 00552 bool PseudoJet::has_constituents() const{ 00553 return (_structure()) && (_structure->has_constituents()); 00554 } 00555 00556 //---------------------------------------------------------------------- 00557 // retrieve the constituents. 00558 vector<PseudoJet> PseudoJet::constituents() const{ 00559 return validated_structure_ptr()->constituents(*this); 00560 } 00561 00562 00563 //---------------------------------------------------------------------- 00564 // returns true if the PseudoJet has support for exclusive subjets 00565 bool PseudoJet::has_exclusive_subjets() const{ 00566 return (_structure()) && (_structure->has_exclusive_subjets()); 00567 } 00568 00569 //---------------------------------------------------------------------- 00570 // return a vector of all subjets of the current jet (in the sense 00571 // of the exclusive algorithm) that would be obtained when running 00572 // the algorithm with the given dcut. 00573 // 00574 // Time taken is O(m ln m), where m is the number of subjets that 00575 // are found. If m gets to be of order of the total number of 00576 // constituents in the jet, this could be substantially slower than 00577 // just getting that list of constituents. 00578 // 00579 // an Error is thrown if this PseudoJet has no currently valid 00580 // associated ClusterSequence 00581 std::vector<PseudoJet> PseudoJet::exclusive_subjets (const double & dcut) const { 00582 return validated_structure_ptr()->exclusive_subjets(*this, dcut); 00583 } 00584 00585 //---------------------------------------------------------------------- 00586 // return the size of exclusive_subjets(...); still n ln n with same 00587 // coefficient, but marginally more efficient than manually taking 00588 // exclusive_subjets.size() 00589 // 00590 // an Error is thrown if this PseudoJet has no currently valid 00591 // associated ClusterSequence 00592 int PseudoJet::n_exclusive_subjets(const double & dcut) const { 00593 return validated_structure_ptr()->n_exclusive_subjets(*this, dcut); 00594 } 00595 00596 //---------------------------------------------------------------------- 00597 // return the list of subjets obtained by unclustering the supplied 00598 // jet down to n subjets (or all constituents if there are fewer 00599 // than n). 00600 // 00601 // requires n ln n time 00602 // 00603 // an Error is thrown if this PseudoJet has no currently valid 00604 // associated ClusterSequence 00605 std::vector<PseudoJet> PseudoJet::exclusive_subjets_up_to (int nsub) const { 00606 return validated_structure_ptr()->exclusive_subjets_up_to(*this, nsub); 00607 } 00608 00609 //---------------------------------------------------------------------- 00610 // Same as exclusive_subjets_up_to but throws an error if there are 00611 // fewer than nsub particles in the jet 00612 std::vector<PseudoJet> PseudoJet::exclusive_subjets (int nsub) const { 00613 vector<PseudoJet> subjets = exclusive_subjets_up_to(nsub); 00614 if (int(subjets.size()) < nsub) { 00615 ostringstream err; 00616 err << "Requested " << nsub << " exclusive subjets, but there were only " 00617 << subjets.size() << " particles in the jet"; 00618 throw Error(err.str()); 00619 } 00620 return subjets; 00621 } 00622 00623 //---------------------------------------------------------------------- 00624 // return the dij that was present in the merging nsub+1 -> nsub 00625 // subjets inside this jet. 00626 // 00627 // an Error is thrown if this PseudoJet has no currently valid 00628 // associated ClusterSequence 00629 double PseudoJet::exclusive_subdmerge(int nsub) const { 00630 return validated_structure_ptr()->exclusive_subdmerge(*this, nsub); 00631 } 00632 00633 //---------------------------------------------------------------------- 00634 // return the maximum dij that occurred in the whole event at the 00635 // stage that the nsub+1 -> nsub merge of subjets occurred inside 00636 // this jet. 00637 // 00638 // an Error is thrown if this PseudoJet has no currently valid 00639 // associated ClusterSequence 00640 double PseudoJet::exclusive_subdmerge_max(int nsub) const { 00641 return validated_structure_ptr()->exclusive_subdmerge_max(*this, nsub); 00642 } 00643 00644 00645 // returns true if a jet has pieces 00646 // 00647 // By default a single particle or a jet coming from a 00648 // ClusterSequence have no pieces and this methos will return false. 00649 bool PseudoJet::has_pieces() const{ 00650 return ((_structure()) && (_structure->has_pieces(*this))); 00651 } 00652 00653 // retrieve the pieces that make up the jet. 00654 // 00655 // By default a jet does not have pieces. 00656 // If the underlying interface supports "pieces" retrieve the 00657 // pieces from there. 00658 std::vector<PseudoJet> PseudoJet::pieces() const{ 00659 return validated_structure_ptr()->pieces(*this); 00660 // if (!has_pieces()) 00661 // throw Error("Trying to retrieve the pieces of a PseudoJet that has no support for pieces."); 00662 // 00663 // return _structure->pieces(*this); 00664 } 00665 00666 00667 //---------------------------------------------------------------------- 00668 // the following ones require a computation of the area in the 00669 // associated ClusterSequence (See ClusterSequenceAreaBase for details) 00670 //---------------------------------------------------------------------- 00671 00672 //---------------------------------------------------------------------- 00673 // if possible, return a valid ClusterSequenceAreaBase pointer; otherwise 00674 // throw an error 00675 const ClusterSequenceAreaBase * PseudoJet::validated_csab() const { 00676 const ClusterSequenceAreaBase *csab = dynamic_cast<const ClusterSequenceAreaBase*>(validated_cs()); 00677 if (csab == NULL) throw Error("you requested jet-area related information, but the PseudoJet does not have associated area information."); 00678 return csab; 00679 } 00680 00681 00682 //---------------------------------------------------------------------- 00683 // check if it has a defined area 00684 bool PseudoJet::has_area() const{ 00685 //if (! has_associated_cluster_sequence()) return false; 00686 if (! has_structure()) return false; 00687 return (validated_structure_ptr()->has_area() != 0); 00688 } 00689 00690 //---------------------------------------------------------------------- 00691 // return the jet (scalar) area. 00692 // throw an Error if there is no support for area in the associated CS 00693 double PseudoJet::area() const{ 00694 return validated_structure_ptr()->area(*this); 00695 } 00696 00697 //---------------------------------------------------------------------- 00698 // return the error (uncertainty) associated with the determination 00699 // of the area of this jet. 00700 // throws an Error if there is no support for area in the associated CS 00701 double PseudoJet::area_error() const{ 00702 return validated_structure_ptr()->area_error(*this); 00703 } 00704 00705 //---------------------------------------------------------------------- 00706 // return the jet 4-vector area 00707 // throws an Error if there is no support for area in the associated CS 00708 PseudoJet PseudoJet::area_4vector() const{ 00709 return validated_structure_ptr()->area_4vector(*this); 00710 } 00711 00712 //---------------------------------------------------------------------- 00713 // true if this jet is made exclusively of ghosts 00714 // throws an Error if there is no support for area in the associated CS 00715 bool PseudoJet::is_pure_ghost() const{ 00716 return validated_structure_ptr()->is_pure_ghost(*this); 00717 } 00718 00719 00720 //---------------------------------------------------------------------- 00721 // 00722 // end of the methods accessing the information in the associated 00723 // Cluster Sequence 00724 // 00725 //---------------------------------------------------------------------- 00726 00727 //---------------------------------------------------------------------- 00728 /// provide a meaningful error message for InexistentUserInfo 00729 PseudoJet::InexistentUserInfo::InexistentUserInfo() : Error("you attempted to perform a dynamic cast of a PseudoJet's extra info, but the extra info pointer was null") 00730 {} 00731 00732 00733 //---------------------------------------------------------------------- 00734 // sort the indices so that values[indices[0..n-1]] is sorted 00735 // into increasing order 00736 void sort_indices(vector<int> & indices, 00737 const vector<double> & values) { 00738 IndexedSortHelper index_sort_helper(&values); 00739 sort(indices.begin(), indices.end(), index_sort_helper); 00740 } 00741 00742 00743 00744 //---------------------------------------------------------------------- 00745 /// given a vector of values with a one-to-one correspondence with the 00746 /// vector of objects, sort objects into an order such that the 00747 /// associated values would be in increasing order 00748 template<class T> vector<T> objects_sorted_by_values( 00749 const vector<T> & objects, 00750 const vector<double> & values) { 00751 00752 assert(objects.size() == values.size()); 00753 00754 // get a vector of indices 00755 vector<int> indices(values.size()); 00756 for (size_t i = 0; i < indices.size(); i++) {indices[i] = i;} 00757 00758 // sort the indices 00759 sort_indices(indices, values); 00760 00761 // copy the objects 00762 vector<T> objects_sorted(objects.size()); 00763 00764 // place the objects in the correct order 00765 for (size_t i = 0; i < indices.size(); i++) { 00766 objects_sorted[i] = objects[indices[i]]; 00767 } 00768 00769 return objects_sorted; 00770 } 00771 00772 //---------------------------------------------------------------------- 00773 /// return a vector of jets sorted into decreasing kt2 00774 vector<PseudoJet> sorted_by_pt(const vector<PseudoJet> & jets) { 00775 vector<double> minus_kt2(jets.size()); 00776 for (size_t i = 0; i < jets.size(); i++) {minus_kt2[i] = -jets[i].kt2();} 00777 return objects_sorted_by_values(jets, minus_kt2); 00778 } 00779 00780 //---------------------------------------------------------------------- 00781 /// return a vector of jets sorted into increasing rapidity 00782 vector<PseudoJet> sorted_by_rapidity(const vector<PseudoJet> & jets) { 00783 vector<double> rapidities(jets.size()); 00784 for (size_t i = 0; i < jets.size(); i++) {rapidities[i] = jets[i].rap();} 00785 return objects_sorted_by_values(jets, rapidities); 00786 } 00787 00788 //---------------------------------------------------------------------- 00789 /// return a vector of jets sorted into decreasing energy 00790 vector<PseudoJet> sorted_by_E(const vector<PseudoJet> & jets) { 00791 vector<double> energies(jets.size()); 00792 for (size_t i = 0; i < jets.size(); i++) {energies[i] = -jets[i].E();} 00793 return objects_sorted_by_values(jets, energies); 00794 } 00795 00796 //---------------------------------------------------------------------- 00797 /// return a vector of jets sorted into increasing pz 00798 vector<PseudoJet> sorted_by_pz(const vector<PseudoJet> & jets) { 00799 vector<double> pz(jets.size()); 00800 for (size_t i = 0; i < jets.size(); i++) {pz[i] = jets[i].pz();} 00801 return objects_sorted_by_values(jets, pz); 00802 } 00803 00804 00805 00806 //------------------------------------------------------------------------------- 00807 // helper functions to build a jet made of pieces 00808 //------------------------------------------------------------------------------- 00809 00810 // build a "CompositeJet" from the vector of its pieces 00811 // 00812 // In this case, E-scheme recombination is assumed to compute the 00813 // total momentum 00814 PseudoJet join(const vector<PseudoJet> & pieces){ 00815 // compute the total momentum 00816 //-------------------------------------------------- 00817 PseudoJet result; // automatically initialised to 0 00818 for (unsigned int i=0; i<pieces.size(); i++) 00819 result += pieces[i]; 00820 00821 // attach a CompositeJetStructure to the result 00822 //-------------------------------------------------- 00823 CompositeJetStructure *cj_struct = new CompositeJetStructure(pieces); 00824 00825 result.set_structure_shared_ptr(SharedPtr<PseudoJetStructureBase>(cj_struct)); 00826 00827 return result; 00828 } 00829 00830 // build a "CompositeJet" from a single PseudoJet 00831 PseudoJet join(const PseudoJet & j1){ 00832 return join(vector<PseudoJet>(1,j1)); 00833 } 00834 00835 // build a "CompositeJet" from two PseudoJet 00836 PseudoJet join(const PseudoJet & j1, const PseudoJet & j2){ 00837 vector<PseudoJet> pieces; 00838 pieces.push_back(j1); 00839 pieces.push_back(j2); 00840 return join(pieces); 00841 } 00842 00843 // build a "CompositeJet" from 3 PseudoJet 00844 PseudoJet join(const PseudoJet & j1, const PseudoJet & j2, const PseudoJet & j3){ 00845 vector<PseudoJet> pieces; 00846 pieces.push_back(j1); 00847 pieces.push_back(j2); 00848 pieces.push_back(j3); 00849 return join(pieces); 00850 } 00851 00852 // build a "CompositeJet" from 4 PseudoJet 00853 PseudoJet join(const PseudoJet & j1, const PseudoJet & j2, const PseudoJet & j3, const PseudoJet & j4){ 00854 vector<PseudoJet> pieces; 00855 pieces.push_back(j1); 00856 pieces.push_back(j2); 00857 pieces.push_back(j3); 00858 pieces.push_back(j4); 00859 return join(pieces); 00860 } 00861 00862 00863 00864 00865 FASTJET_END_NAMESPACE 00866