FastJet 3.0beta1
ClusterSequence.cc
00001 //STARTHEADER
00002 // $Id: ClusterSequence.cc 2333 2011-07-01 08:59:29Z salam $
00003 //
00004 // Copyright (c) 2005-2011, 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 #include "fastjet/Error.hh"
00032 #include "fastjet/PseudoJet.hh"
00033 #include "fastjet/ClusterSequence.hh"
00034 #include "fastjet/ClusterSequenceStructure.hh"
00035 #include "fastjet/version.hh" // stores the current version number
00036 #include<iostream>
00037 #include<sstream>
00038 #include<fstream>
00039 #include<cmath>
00040 #include<cstdlib>
00041 #include<cassert>
00042 #include<string>
00043 #include<set>
00044 
00045 FASTJET_BEGIN_NAMESPACE      // defined in fastjet/internal/base.hh
00046 
00047 using namespace std;
00048 
00049 //// initialised static member has to go in the .cc code
00050 JetAlgorithm ClusterSequence::_default_jet_algorithm = kt_algorithm;
00051 //
00052 
00053 
00054 // destructor that guarantees proper bookkeeping for the CS Structure
00055 ClusterSequence::~ClusterSequence () {
00056   // set the pointer in the wrapper to this object to NULL to say that
00057   // we're going out of scope
00058   if (_structure_shared_ptr()){
00059     ClusterSequenceStructure* csi = dynamic_cast<ClusterSequenceStructure*>(_structure_shared_ptr()); 
00060     // normally the csi is purely internal so it really should not be
00061     // NULL i.e assert should be OK
00062     // (we assert rather than throw an error, since failure here is a
00063     // sign of major internal problems)
00064     assert(csi != NULL);
00065     csi->set_associated_cs(NULL);
00066 
00067     // if the user had given the CS responsibility to delete itself,
00068     // but then deletes the CS themselves, the following lines of
00069     // code will ensure that the structure_shared_ptr will have
00070     // a proper object count (so that jets associated with the CS will
00071     // throw the correct error if the user tries to access their
00072     // constituents).
00073     if (_deletes_self_when_unused) {
00074       _structure_shared_ptr.set_count(_structure_shared_ptr.use_count() 
00075                                         + _structure_use_count_after_construction);
00076     }
00077   }
00078 }
00079 
00080 //-----------
00081 void ClusterSequence::signal_imminent_self_deletion() const {
00082   // normally if the destructor is called when
00083   // _deletes_self_when_unused is true, it assumes that it's been
00084   // called by the user (and it therefore resets the shared pointer
00085   // count to the true count).
00086   //
00087   // for self deletion (called from the destructor of the CSstructure,
00088   // the shared_ptr to which has just had its pointer -> 0) you do
00089   // _not_ want to reset the pointer count (otherwise you will end up
00090   // with a double delete on the shared pointer once you start
00091   // deleting the internal structure of the CS).
00092   //
00093   // the following modification ensures that the count reset will not
00094   // take place in the destructor
00095   assert(_deletes_self_when_unused);
00096   _deletes_self_when_unused = false;
00097 }
00098 
00099 //----------------------------------------------------------------------
00100 void ClusterSequence::_initialise_and_run (
00101                                   const double & R,
00102                                   const Strategy & strategy,
00103                                   const bool & writeout_combinations) {
00104 
00105   JetDefinition jet_def(_default_jet_algorithm, R, strategy);
00106   _initialise_and_run(jet_def, writeout_combinations);
00107 }
00108 
00109 
00110 //----------------------------------------------------------------------
00111 void ClusterSequence::_initialise_and_run (
00112                                   const JetDefinition & jet_def,
00113                                   const bool & writeout_combinations) {
00114 
00115   // transfer all relevant info into internal variables
00116   _decant_options(jet_def, writeout_combinations);
00117 
00118   // set up the history entries for the initial particles (those
00119   // currently in _jets)
00120   _fill_initial_history();
00121 
00122   // don't run anything if the event is empty
00123   if (n_particles() == 0) return;
00124 
00125   // ----- deal with special cases: plugins & e+e- ------
00126   if (_jet_algorithm == plugin_algorithm) {
00127     // allows plugin_xyz() functions to modify cluster sequence
00128     _plugin_activated = true;
00129     // let the plugin do its work here
00130     _jet_def.plugin()->run_clustering( (*this) );
00131     _plugin_activated = false;
00132     _update_structure_use_count();
00133     return;
00134   } else if (_jet_algorithm == ee_kt_algorithm ||
00135              _jet_algorithm == ee_genkt_algorithm) {
00136     // ignore requested strategy
00137     _strategy = N2Plain;
00138     if (_jet_algorithm == ee_kt_algorithm) {
00139       // make sure that R is large enough so that "beam" recomb only
00140       // occurs when a single particle is left
00141       // Normally, this should be automatically set to 4 from JetDefinition
00142       assert(_Rparam > 2.0); 
00143       // this is used to renormalise the dij to get a "standard" form
00144       // and our convention in e+e- will be different from that
00145       // in long.inv case; NB: _invR2 name should be changed -> _renorm_dij?
00146       _invR2 = 1.0;
00147     } else {
00148       // as of 2009-01-09, choose R to be an angular distance, in
00149       // radians.  Since the algorithm uses 2(1-cos(theta)) as its
00150       // squared angular measure, make sure that the _R2 is defined
00151       // in a similar way.
00152       if (_Rparam > pi) {
00153         // choose a value that ensures that back-to-back particles will
00154         // always recombine 
00155         //_R2 = 4.0000000000001;
00156         _R2 = 2 * ( 3.0 + cos(_Rparam) );
00157       } else {
00158         _R2    = 2 * ( 1.0 - cos(_Rparam) );
00159       }
00160       _invR2 = 1.0/_R2;
00161     }
00162     _simple_N2_cluster_EEBriefJet();
00163     return;
00164   } else if (_jet_algorithm == undefined_jet_algorithm) {
00165     throw Error("A ClusterSequence cannot be created with an uninitialised JetDefinition");
00166   }
00167 
00168 
00169   // automatically redefine the strategy according to N if that is
00170   // what the user requested -- transition points (and especially
00171   // their R-dependence) are based on empirical observations for a
00172   // R=0.4, 0.7 and 1.0, running on toth (3.4GHz, Pentium IV D [dual
00173   // core] with 2MB of cache).
00174   if (_strategy == Best) {
00175     int N = _jets.size();
00176     if (N <= 55*max(0.5,min(1.0,_Rparam))) {// empirical scaling with R
00177       _strategy = N2Plain;
00178     } else if (N > 6200/pow(_Rparam,2.0) && jet_def.jet_algorithm() == cambridge_algorithm) {
00179       _strategy = NlnNCam;
00180 #ifndef DROP_CGAL
00181     } else if ((N > 16000/pow(_Rparam,1.15) && jet_def.jet_algorithm() != antikt_algorithm)
00182                || N > 35000/pow(_Rparam,1.15)) {
00183       _strategy = NlnN;
00184 #endif  // DROP_CGAL
00185     } else if (N <= 450) {
00186       _strategy = N2Tiled;
00187     } else {                   
00188       _strategy = N2MinHeapTiled;
00189     }
00190   }
00191 
00192   // R >= 2pi is not supported by all clustering strategies owing to
00193   // periodicity issues (a particle might cluster with itself). When
00194   // R>=2pi, we therefore automatically switch to a strategy that is
00195   // known to work.
00196   if (_Rparam >= twopi) {
00197     if (   _strategy == NlnN
00198         || _strategy == NlnN3pi
00199         || _strategy == NlnNCam
00200         || _strategy == NlnNCam2pi2R
00201         || _strategy == NlnNCam4pi) {
00202 #ifdef DROP_CGAL
00203       _strategy = N2MinHeapTiled;
00204 #else
00205       _strategy = NlnN4pi;
00206 #endif    
00207     }
00208     if (jet_def.strategy() != Best && _strategy != jet_def.strategy()) {
00209       ostringstream oss;
00210       oss << "Cluster strategy " << strategy_string(jet_def.strategy())
00211           << " automatically changed to " << strategy_string()
00212           << " because the former is not supported for R = " << _Rparam
00213           << " >= 2pi";
00214       _changed_strategy_warning.warn(oss.str());
00215     }
00216   }
00217 
00218 
00219   // run the code containing the selected strategy
00220   // 
00221   // We order the strategies stqrting from the ones used by the Best
00222   // strategy in the order of increasing N, then the remaining ones
00223   // again in the order of increasing N.
00224   if (_strategy == N2Plain) {
00225     // BriefJet provides standard long.invariant kt alg.
00226     this->_simple_N2_cluster_BriefJet();
00227   } else if (_strategy == N2Tiled) {
00228     this->_faster_tiled_N2_cluster();
00229   } else if (_strategy == N2MinHeapTiled) {
00230     this->_minheap_faster_tiled_N2_cluster();
00231   } else if (_strategy == NlnN) {
00232     this->_delaunay_cluster();
00233   } else if (_strategy == NlnNCam) {
00234     this->_CP2DChan_cluster_2piMultD();
00235   } else if (_strategy == NlnN3pi || _strategy == NlnN4pi ) {
00236     this->_delaunay_cluster();
00237   } else if (_strategy ==  N3Dumb ) {
00238     this->_really_dumb_cluster();
00239   } else if (_strategy == N2PoorTiled) {
00240     this->_tiled_N2_cluster();
00241   } else if (_strategy == NlnNCam4pi) {
00242     this->_CP2DChan_cluster();
00243   } else if (_strategy == NlnNCam2pi2R) {
00244     this->_CP2DChan_cluster_2pi2R();
00245   } else {
00246     ostringstream err;
00247     err << "Unrecognised value for strategy: "<<_strategy;
00248     throw Error(err.str());
00249   }
00250 
00251 }
00252 
00253 
00254 // these needs to be defined outside the class definition.
00255 bool ClusterSequence::_first_time = true;
00256 int ClusterSequence::_n_exclusive_warnings = 0;
00257 
00258 
00259 //----------------------------------------------------------------------
00260 // the version string
00261 string fastjet_version_string() {
00262   return "FastJet version "+string(fastjet_version);
00263 }
00264 
00265 
00266 //----------------------------------------------------------------------
00267 // prints a banner on the first call
00268 void ClusterSequence::_print_banner() {
00269 
00270   if (!_first_time) {return;}
00271   _first_time = false;
00272   
00273   
00274   //Symp. Discr. Alg, p.472 (2002) and  CGAL (http://www.cgal.org);
00275 
00276   cout << "#--------------------------------------------------------------------------\n";
00277   cout << "#                         FastJet release " << fastjet_version << endl;
00278   cout << "#            Written by M. Cacciari, G.P. Salam and G. Soyez            \n"; 
00279   cout << "#                         http://www.fastjet.fr                         \n"; 
00280   cout << "#                                                                       \n";
00281   cout << "# Longitudinally invariant Kt, anti-Kt, and inclusive Cambridge/Aachen  \n";
00282   cout << "# clustering using fast geometric algorithms, with jet areas and optional\n";
00283   cout << "# external jet-finder plugins. If you use this code towards a scientific \n";
00284   cout << "# publication please cite Phys. Lett. B641 (2006) [hep-ph/0512210] and   \n";
00285   cout << "# M. Cacciari, G.P. Salam and G. Soyez, http://fastjet.fr/              \n";
00286   cout << "#                                                                       \n";
00287   cout << "# This package uses T.Chan's closest pair algorithm, Proc.13th ACM-SIAM \n";
00288   cout << "# Symp. Discr. Alg, p.472 (2002), S.Fortune's Voronoi algorithm and code" ;
00289 #ifndef DROP_CGAL
00290   cout << endl << "# and CGAL: http://www.cgal.org/";
00291 #endif  // DROP_CGAL
00292   cout << ".\n";
00293   cout << "#-------------------------------------------------------------------------\n";
00294   // make sure we really have the output done.
00295   cout.flush();
00296 }
00297 
00298 //----------------------------------------------------------------------
00299 // transfer all relevant info into internal variables
00300 void ClusterSequence::_decant_options(const JetDefinition & jet_def,
00301                                       const bool & writeout_combinations) {
00302   // let the user know what's going on
00303   _print_banner();
00304 
00305   // make a local copy of the jet definition (for future use?)
00306   _jet_def = jet_def;
00307   
00308   _writeout_combinations = writeout_combinations;
00309   _jet_algorithm = jet_def.jet_algorithm();
00310   _Rparam = jet_def.R();  _R2 = _Rparam*_Rparam; _invR2 = 1.0/_R2;
00311   _strategy = jet_def.strategy();
00312 
00313   // disallow interference from the plugin
00314   _plugin_activated = false;
00315 
00316   // initialised the wrapper to the current CS
00317   _structure_shared_ptr.reset(new ClusterSequenceStructure(this));
00318   _update_structure_use_count(); // make sure it's correct already here
00319 }
00320 
00321 
00322 //----------------------------------------------------------------------
00323 // initialise the history in a standard way
00324 void ClusterSequence::_fill_initial_history () {
00325 
00326   //if (_jets.size() == 0) {throw Error("Cannot run jet-finder on empty event");}
00327 
00328   // reserve sufficient space for everything
00329   _jets.reserve(_jets.size()*2);
00330   _history.reserve(_jets.size()*2);
00331 
00332   _Qtot = 0;
00333 
00334   for (int i = 0; i < static_cast<int>(_jets.size()) ; i++) {
00335     history_element element;
00336     element.parent1 = InexistentParent;
00337     element.parent2 = InexistentParent;
00338     element.child   = Invalid;
00339     element.jetp_index = i;
00340     element.dij     = 0.0;
00341     element.max_dij_so_far = 0.0;
00342 
00343     _history.push_back(element);
00344     
00345     // do any momentum preprocessing needed by the recombination scheme
00346     _jet_def.recombiner()->preprocess(_jets[i]);
00347 
00348     // get cross-referencing right from PseudoJets
00349     _jets[i].set_cluster_hist_index(i);
00350     _set_structure_shared_ptr(_jets[i]);
00351 
00352     // determine the total energy in the event
00353     _Qtot += _jets[i].E();
00354   }
00355   _initial_n = _jets.size();
00356   _deletes_self_when_unused = false;
00357 }
00358 
00359 
00360 //----------------------------------------------------------------------
00361 // Return the component corresponding to the specified index.
00362 // taken from CLHEP
00363 string ClusterSequence::strategy_string (Strategy strategy_in)  const {
00364   string strategy;
00365   switch(strategy_in) {
00366   case NlnN:
00367     strategy = "NlnN"; break;
00368   case NlnN3pi:
00369     strategy = "NlnN3pi"; break;
00370   case NlnN4pi:
00371     strategy = "NlnN4pi"; break;
00372   case N2Plain:
00373     strategy = "N2Plain"; break;
00374   case N2Tiled:
00375     strategy = "N2Tiled"; break;
00376   case N2MinHeapTiled:
00377     strategy = "N2MinHeapTiled"; break;
00378   case N2PoorTiled:
00379     strategy = "N2PoorTiled"; break;
00380   case N3Dumb:
00381     strategy = "N3Dumb"; break;
00382   case NlnNCam4pi:
00383     strategy = "NlnNCam4pi"; break;
00384   case NlnNCam2pi2R:
00385     strategy = "NlnNCam2pi2R"; break;
00386   case NlnNCam:
00387     strategy = "NlnNCam"; break; // 2piMultD
00388   case plugin_strategy:
00389     strategy = "plugin strategy"; break;
00390   default:
00391     strategy = "Unrecognized";
00392   }
00393   return strategy;
00394 }  
00395 
00396 
00397 double ClusterSequence::jet_scale_for_algorithm(
00398                                   const PseudoJet & jet) const {
00399   if (_jet_algorithm == kt_algorithm)             {return jet.kt2();}
00400   else if (_jet_algorithm == cambridge_algorithm) {return 1.0;}
00401   else if (_jet_algorithm == antikt_algorithm) {
00402     double kt2=jet.kt2();
00403     return kt2 > 1e-300 ? 1.0/kt2 : 1e300;
00404   } else if (_jet_algorithm == genkt_algorithm) {
00405     double kt2 = jet.kt2();
00406     double p   = jet_def().extra_param();
00407     if (p <= 0 && kt2 < 1e-300) kt2 = 1e-300; // dodgy safety check
00408     return pow(kt2, p);
00409   } else if (_jet_algorithm == cambridge_for_passive_algorithm) {
00410     double kt2 = jet.kt2();
00411     double lim = _jet_def.extra_param();
00412     if (kt2 < lim*lim && kt2 != 0.0) {
00413       return 1.0/kt2;
00414     } else {return 1.0;}
00415   } else {throw Error("Unrecognised jet algorithm");}
00416 }
00417 
00418 
00419 // //----------------------------------------------------------------------
00420 // /// transfer the sequence contained in other_seq into our own;
00421 // /// any plugin "extras" contained in the from_seq will be lost
00422 // /// from there.
00423 // void ClusterSequence::transfer_from_sequence(ClusterSequence & from_seq) {
00424 // 
00425 //   if (will_delete_self_when_unused()) 
00426 //     throw(Error("cannot use CS::transfer_from_sequence after a call to delete_self_when_unused()"));
00427 // 
00428 //   // the metadata
00429 //   _jet_def                 = from_seq._jet_def                ;
00430 //   _writeout_combinations   = from_seq._writeout_combinations  ;
00431 //   _initial_n               = from_seq._initial_n              ;
00432 //   _Rparam                  = from_seq._Rparam                 ;
00433 //   _R2                      = from_seq._R2                     ;
00434 //   _invR2                   = from_seq._invR2                  ;
00435 //   _strategy                = from_seq._strategy               ;
00436 //   _jet_algorithm           = from_seq._jet_algorithm          ;
00437 //   _plugin_activated        = from_seq._plugin_activated       ;
00438 // 
00439 //   // the data
00440 //   _jets     = from_seq._jets;
00441 //   _history  = from_seq._history;
00442 //   // the following transfers ownership of the extras from the from_seq
00443 //   _extras   = from_seq._extras;
00444 // 
00445 //   // transfer of ownership
00446 //   if (_structure_shared_ptr()) {
00447 //     // anything that is currently associated with the cluster sequence
00448 //     // should be told that its cluster sequence no longer exists
00449 //     ClusterSequenceStructure* csi = dynamic_cast<ClusterSequenceStructure*>(_structure_shared_ptr()); 
00450 //     assert(csi != NULL);
00451 //     csi->set_associated_cs(NULL);
00452 //   }
00453 //   // create a new _structure_shared_ptr to reflect the fact that
00454 //   // this CS is essentially a new one
00455 //   _structure_shared_ptr.reset(new ClusterSequenceStructure(this));
00456 //   _update_structure_use_count();
00457 //   
00458 //   for (vector<PseudoJet>::iterator jit = _jets.begin(); jit != _jets.end(); jit++)
00459 //     _set_structure_shared_ptr(*jit);
00460 // }
00461 
00462 
00463 //----------------------------------------------------------------------
00464 // transfer the sequence contained in other_seq into our own;
00465 // any plugin "extras" contained in the from_seq will be lost
00466 // from there.
00467 //
00468 // It also sets the ClusterSequence pointers of the PseudoJets in
00469 // the history to point to this ClusterSequence
00470 //
00471 // The second argument is an action that will be applied on every
00472 // jets in the resulting ClusterSequence
00473 void ClusterSequence::transfer_from_sequence(const ClusterSequence & from_seq,
00474                                              const FunctionOfPseudoJet<PseudoJet> * action_on_jets){
00475 
00476   if (will_delete_self_when_unused()) 
00477     throw(Error("cannot use CS::transfer_from_sequence after a call to delete_self_when_unused()"));
00478 
00479   // the metadata
00480   _jet_def                 = from_seq._jet_def                ;
00481   _writeout_combinations   = from_seq._writeout_combinations  ;
00482   _initial_n               = from_seq._initial_n              ;
00483   _Rparam                  = from_seq._Rparam                 ;
00484   _R2                      = from_seq._R2                     ;
00485   _invR2                   = from_seq._invR2                  ;
00486   _strategy                = from_seq._strategy               ;
00487   _jet_algorithm           = from_seq._jet_algorithm          ;
00488   _plugin_activated        = from_seq._plugin_activated       ;
00489 
00490   // the data
00491 
00492   // apply the transformation on the jets if needed
00493   if (action_on_jets)
00494     _jets     = (*action_on_jets)(from_seq._jets);
00495   else
00496     _jets     = from_seq._jets;
00497   _history  = from_seq._history;
00498   // the following shares ownership of the extras with the from_seq;
00499   // no transformations will be applied to the extras
00500   _extras   = from_seq._extras;
00501 
00502   // clean up existing structure
00503   if (_structure_shared_ptr()) {
00504     // If there are jets associated with an old version of the CS and
00505     // a new one, keeping track of when to delete the CS becomes more
00506     // complex; so we don't allow this situation to occur.
00507     if (_deletes_self_when_unused) throw Error("transfer_from_sequence cannot be used for a cluster sequence that deletes self when unused");
00508     
00509     // anything that is currently associated with the cluster sequence
00510     // should be told that its cluster sequence no longer exists
00511     ClusterSequenceStructure* csi = dynamic_cast<ClusterSequenceStructure*>(_structure_shared_ptr()); 
00512     assert(csi != NULL);
00513     csi->set_associated_cs(NULL);
00514   }
00515   // create a new _structure_shared_ptr to reflect the fact that
00516   // this CS is essentially a new one
00517   _structure_shared_ptr.reset(new ClusterSequenceStructure(this));
00518   _update_structure_use_count();
00519   
00520   for (unsigned int i=0; i<_jets.size(); i++){
00521     // we reset the cluster history index in case action_on_jets
00522     // messed up with it
00523     _jets[i].set_cluster_hist_index(from_seq._jets[i].cluster_hist_index());
00524 
00525     // reset the structure pointer
00526     _set_structure_shared_ptr(_jets[i]);
00527   }
00528 }
00529 
00530 
00531 //----------------------------------------------------------------------
00532 // record an ij recombination and reset the _jets[newjet_k] momentum and
00533 // user index to be those of newjet
00534 void ClusterSequence::plugin_record_ij_recombination(
00535            int jet_i, int jet_j, double dij, 
00536            const PseudoJet & newjet, int & newjet_k) {
00537 
00538   plugin_record_ij_recombination(jet_i, jet_j, dij, newjet_k);
00539 
00540   // now transfer newjet into place
00541   int tmp_index = _jets[newjet_k].cluster_hist_index();
00542   _jets[newjet_k] = newjet;
00543   _jets[newjet_k].set_cluster_hist_index(tmp_index);
00544   _set_structure_shared_ptr(_jets[newjet_k]);
00545 }
00546 
00547 
00548 //----------------------------------------------------------------------
00549 // return all inclusive jets with pt > ptmin
00550 vector<PseudoJet> ClusterSequence::inclusive_jets (const double & ptmin) const{
00551   double dcut = ptmin*ptmin;
00552   int i = _history.size() - 1; // last jet
00553   vector<PseudoJet> jets;
00554   if (_jet_algorithm == kt_algorithm) {
00555     while (i >= 0) {
00556       // with our specific definition of dij and diB (i.e. R appears only in 
00557       // dij), then dij==diB is the same as the jet.perp2() and we can exploit
00558       // this in selecting the jets...
00559       if (_history[i].max_dij_so_far < dcut) {break;}
00560       if (_history[i].parent2 == BeamJet && _history[i].dij >= dcut) {
00561         // for beam jets
00562         int parent1 = _history[i].parent1;
00563         jets.push_back(_jets[_history[parent1].jetp_index]);}
00564       i--;
00565     }
00566   } else if (_jet_algorithm == cambridge_algorithm) {
00567     while (i >= 0) {
00568       // inclusive jets are all at end of clustering sequence in the
00569       // Cambridge algorithm -- so if we find a non-exclusive jet, then
00570       // we can exit
00571       if (_history[i].parent2 != BeamJet) {break;}
00572       int parent1 = _history[i].parent1;
00573       const PseudoJet & jet = _jets[_history[parent1].jetp_index];
00574       if (jet.perp2() >= dcut) {jets.push_back(jet);}
00575       i--;
00576     }
00577   } else if (_jet_algorithm == plugin_algorithm 
00578              || _jet_algorithm == ee_kt_algorithm
00579              || _jet_algorithm == antikt_algorithm
00580              || _jet_algorithm == genkt_algorithm
00581              || _jet_algorithm == ee_genkt_algorithm
00582              || _jet_algorithm == cambridge_for_passive_algorithm) {
00583     // for inclusive jets with a plugin algorithm, we make no
00584     // assumptions about anything (relation of dij to momenta,
00585     // ordering of the dij, etc.)
00586     while (i >= 0) {
00587       if (_history[i].parent2 == BeamJet) {
00588         int parent1 = _history[i].parent1;
00589         const PseudoJet & jet = _jets[_history[parent1].jetp_index];
00590         if (jet.perp2() >= dcut) {jets.push_back(jet);}
00591       }
00592       i--;
00593     }
00594   } else {throw Error("cs::inclusive_jets(...): Unrecognized jet algorithm");}
00595   return jets;
00596 }
00597 
00598 
00599 //----------------------------------------------------------------------
00600 // return the number of exclusive jets that would have been obtained
00601 // running the algorithm in exclusive mode with the given dcut
00602 int ClusterSequence::n_exclusive_jets (const double & dcut) const {
00603 
00604   // first locate the point where clustering would have stopped (i.e. the
00605   // first time max_dij_so_far > dcut)
00606   int i = _history.size() - 1; // last jet
00607   while (i >= 0) {
00608     if (_history[i].max_dij_so_far <= dcut) {break;}
00609     i--;
00610   }
00611   int stop_point = i + 1;
00612   // relation between stop_point, njets assumes one extra jet disappears
00613   // at each clustering.
00614   int njets = 2*_initial_n - stop_point;
00615   return njets;
00616 }
00617 
00618 //----------------------------------------------------------------------
00619 // return all exclusive jets that would have been obtained running
00620 // the algorithm in exclusive mode with the given dcut
00621 vector<PseudoJet> ClusterSequence::exclusive_jets (const double & dcut) const {
00622   int njets = n_exclusive_jets(dcut);
00623   return exclusive_jets(njets);
00624 }
00625 
00626 
00627 //----------------------------------------------------------------------
00628 // return the jets obtained by clustering the event to n jets.
00629 // Throw an error if there are fewer than n particles.
00630 vector<PseudoJet> ClusterSequence::exclusive_jets (const int & njets) const {
00631 
00632   // make sure the user does not ask for more than jets than there
00633   // were particles in the first place.
00634   if (njets > _initial_n) {
00635     ostringstream err;
00636     err << "Requested " << njets << " exclusive jets, but there were only " 
00637         << _initial_n << " particles in the event";
00638     throw Error(err.str());
00639   }
00640 
00641   return exclusive_jets_up_to(njets);
00642 }
00643 
00644 //----------------------------------------------------------------------
00645 // return the jets obtained by clustering the event to n jets.
00646 // If there are fewer than n particles, simply return all particles
00647 vector<PseudoJet> ClusterSequence::exclusive_jets_up_to (const int & njets) const {
00648 
00649   // provide a warning when extracting exclusive jets for algorithms 
00650   // that does not support it explicitly.
00651   // Native algorithm that support it are: kt, ee_kt, cambridge, 
00652   //   genkt and ee_genkt (both with p>=0)
00653   // For plugins, we check Plugin::exclusive_sequence_meaningful()
00654   if (( _jet_def.jet_algorithm() != kt_algorithm) &&
00655       ( _jet_def.jet_algorithm() != cambridge_algorithm) &&
00656       ( _jet_def.jet_algorithm() != ee_kt_algorithm) &&
00657       (((_jet_def.jet_algorithm() != genkt_algorithm) && 
00658         (_jet_def.jet_algorithm() != ee_genkt_algorithm)) || 
00659        (_jet_def.extra_param() <0)) &&
00660       ((_jet_def.jet_algorithm() != plugin_algorithm) ||
00661        (!_jet_def.plugin()->exclusive_sequence_meaningful())) &&
00662       (_n_exclusive_warnings < 5)) {
00663     _n_exclusive_warnings++;
00664     cerr << "FastJet WARNING: dcut and exclusive jets for jet-finders other than kt should be interpreted with care." << endl;
00665   }
00666 
00667 
00668   // calculate the point where we have to stop the clustering.
00669   // relation between stop_point, njets assumes one extra jet disappears
00670   // at each clustering.
00671   int stop_point = 2*_initial_n - njets;
00672   // make sure it's safe when more jets are requested than there are particles
00673   if (stop_point < _initial_n) stop_point = _initial_n;
00674 
00675   // some sanity checking to make sure that e+e- does not give us
00676   // surprises (should we ever implement e+e-)...
00677   if (2*_initial_n != static_cast<int>(_history.size())) {
00678     ostringstream err;
00679     err << "2*_initial_n != _history.size() -- this endangers internal assumptions!\n";
00680     throw Error(err.str());
00681     //assert(false);
00682   }
00683 
00684   // now go forwards and reconstitute the jets that we have --
00685   // basically for any history element, see if the parent jets to
00686   // which it refers were created before the stopping point -- if they
00687   // were then add them to the list, otherwise they are subsequent
00688   // recombinations of the jets that we are looking for.
00689   vector<PseudoJet> jets;
00690   for (unsigned int i = stop_point; i < _history.size(); i++) {
00691     int parent1 = _history[i].parent1;
00692     if (parent1 < stop_point) {
00693       jets.push_back(_jets[_history[parent1].jetp_index]);
00694     }
00695     int parent2 = _history[i].parent2;
00696     if (parent2 < stop_point && parent2 > 0) {
00697       jets.push_back(_jets[_history[parent2].jetp_index]);
00698     }
00699     
00700   }
00701 
00702   // sanity check...
00703   if (int(jets.size()) != min(_initial_n, njets)) {
00704     ostringstream err;
00705     err << "ClusterSequence::exclusive_jets: size of returned vector ("
00706          <<jets.size()<<") does not coincide with requested number of jets ("
00707          <<njets<<")";
00708     throw Error(err.str());
00709   }
00710 
00711   return jets;
00712 }
00713 
00714 //----------------------------------------------------------------------
00715 /// return the dmin corresponding to the recombination that went from
00716 /// n+1 to n jets
00717 double ClusterSequence::exclusive_dmerge (const int & njets) const {
00718   assert(njets >= 0);
00719   if (njets >= _initial_n) {return 0.0;}
00720   return _history[2*_initial_n-njets-1].dij;
00721 }
00722 
00723 
00724 //----------------------------------------------------------------------
00725 /// return the maximum of the dmin encountered during all recombinations 
00726 /// up to the one that led to an n-jet final state; identical to
00727 /// exclusive_dmerge, except in cases where the dmin do not increase
00728 /// monotonically.
00729 double ClusterSequence::exclusive_dmerge_max (const int & njets) const {
00730   assert(njets >= 0);
00731   if (njets >= _initial_n) {return 0.0;}
00732   return _history[2*_initial_n-njets-1].max_dij_so_far;
00733 }
00734 
00735 
00736 //----------------------------------------------------------------------
00737 /// return a vector of all subjets of the current jet (in the sense
00738 /// of the exclusive algorithm) that would be obtained when running
00739 /// the algorithm with the given dcut.
00740 std::vector<PseudoJet> ClusterSequence::exclusive_subjets 
00741    (const PseudoJet & jet, const double & dcut) const {
00742 
00743   set<const history_element*> subhist;
00744 
00745   // get the set of history elements that correspond to subjets at
00746   // scale dcut
00747   get_subhist_set(subhist, jet, dcut, 0);
00748 
00749   // now transfer this into a sequence of jets
00750   vector<PseudoJet> subjets;
00751   subjets.reserve(subhist.size());
00752   for (set<const history_element*>::iterator elem = subhist.begin(); 
00753        elem != subhist.end(); elem++) {
00754     subjets.push_back(_jets[(*elem)->jetp_index]);
00755   }
00756   return subjets;
00757 }
00758 
00759 //----------------------------------------------------------------------
00760 /// return the size of exclusive_subjets(...); still n ln n with same
00761 /// coefficient, but marginally more efficient than manually taking
00762 /// exclusive_subjets.size()
00763 int ClusterSequence::n_exclusive_subjets(const PseudoJet & jet, 
00764                         const double & dcut) const {
00765   set<const history_element*> subhist;
00766   // get the set of history elements that correspond to subjets at
00767   // scale dcut
00768   get_subhist_set(subhist, jet, dcut, 0);
00769   return subhist.size();
00770 }
00771 
00772 //----------------------------------------------------------------------
00773 /// return the list of subjets obtained by unclustering the supplied
00774 /// jet down to nsub subjets. Throws an error if there are fewer than
00775 /// nsub particles in the jet.
00776 std::vector<PseudoJet> ClusterSequence::exclusive_subjets
00777    (const PseudoJet & jet, int nsub) const {
00778   vector<PseudoJet> subjets = exclusive_subjets_up_to(jet, nsub);
00779   if (int(subjets.size()) < nsub) {
00780     ostringstream err;
00781     err << "Requested " << nsub << " exclusive subjets, but there were only " 
00782         << subjets.size() << " particles in the jet";
00783     throw Error(err.str());
00784   }
00785   return subjets;
00786 
00787 }
00788 
00789 //----------------------------------------------------------------------
00790 /// return the list of subjets obtained by unclustering the supplied
00791 /// jet down to nsub subjets (or all constituents if there are fewer
00792 /// than nsub).
00793 std::vector<PseudoJet> ClusterSequence::exclusive_subjets_up_to
00794    (const PseudoJet & jet, int nsub) const {
00795 
00796   set<const history_element*> subhist;
00797 
00798   // prepare the vector into which we'll put the result
00799   vector<PseudoJet> subjets;
00800   if (nsub <  0) throw Error("Requested a negative number of subjets. This is nonsensical.");
00801   if (nsub == 0) return subjets;
00802 
00803   // get the set of history elements that correspond to subjets at
00804   // scale dcut
00805   get_subhist_set(subhist, jet, -1.0, nsub);
00806 
00807   // now transfer this into a sequence of jets
00808   subjets.reserve(subhist.size());
00809   for (set<const history_element*>::iterator elem = subhist.begin(); 
00810        elem != subhist.end(); elem++) {
00811     subjets.push_back(_jets[(*elem)->jetp_index]);
00812   }
00813   return subjets;
00814 }
00815 
00816 
00817 //----------------------------------------------------------------------
00818 /// return the dij that was present in the merging nsub+1 -> nsub 
00819 /// subjets inside this jet.
00820 /// 
00821 /// If the jet has nsub or fewer constituents, it will return 0.
00822 double ClusterSequence::exclusive_subdmerge(const PseudoJet & jet, int nsub) const {
00823   set<const history_element*> subhist;
00824 
00825   // get the set of history elements that correspond to subjets at
00826   // scale dcut
00827   get_subhist_set(subhist, jet, -1.0, nsub);
00828   
00829   set<const history_element*>::iterator highest = subhist.end();
00830   highest--;
00831   /// will be zero if nconst <= nsub, since highest will be an original 
00832   /// particle have zero dij
00833   return (*highest)->dij;
00834 }
00835 
00836 
00837 //----------------------------------------------------------------------
00838 /// return the maximum dij that occurred in the whole event at the
00839 /// stage that the nsub+1 -> nsub merge of subjets occurred inside 
00840 /// this jet.
00841 ///
00842 /// If the jet has nsub or fewer constituents, it will return 0.
00843 double ClusterSequence::exclusive_subdmerge_max(const PseudoJet & jet, int nsub) const {
00844 
00845   set<const history_element*> subhist;
00846 
00847   // get the set of history elements that correspond to subjets at
00848   // scale dcut
00849   get_subhist_set(subhist, jet, -1.0, nsub);
00850   
00851   set<const history_element*>::iterator highest = subhist.end();
00852   highest--;
00853   /// will be zero if nconst <= nsub, since highest will be an original 
00854   /// particle have zero dij
00855   return (*highest)->max_dij_so_far;
00856 }
00857 
00858 
00859 
00860 //----------------------------------------------------------------------
00861 /// return a set of pointers to history entries corresponding to the
00862 /// subjets of this jet; one stops going working down through the
00863 /// subjets either when 
00864 ///   - there is no further to go
00865 ///   - one has found maxjet entries
00866 ///   - max_dij_so_far <= dcut
00867 void ClusterSequence::get_subhist_set(set<const history_element*> & subhist,
00868                                      const  PseudoJet & jet, 
00869                                      double dcut, int maxjet) const {
00870   assert(contains(jet));
00871   
00872   subhist.clear();
00873   subhist.insert(&(_history[jet.cluster_hist_index()]));
00874 
00875   // establish the set of jets that are relevant
00876   int njet = 1;
00877   while (true) {
00878     // first find out if we need to probe deeper into jet.
00879     // Get history element closest to end of sequence
00880     set<const history_element*>::iterator highest = subhist.end();
00881     assert (highest != subhist.begin()); 
00882     highest--;
00883     const history_element* elem = *highest;
00884     // make sure we haven't got too many jets
00885     if (njet == maxjet) break;
00886     // make sure it has parents
00887     if (elem->parent1 < 0)            break;
00888     // make sure that we still resolve it at scale dcut
00889     if (elem->max_dij_so_far <= dcut) break;
00890 
00891     // then do so: replace "highest" with its two parents
00892     subhist.erase(highest);
00893     subhist.insert(&(_history[elem->parent1]));
00894     subhist.insert(&(_history[elem->parent2]));
00895     njet++;
00896   }
00897 }
00898 
00899 //----------------------------------------------------------------------
00900 // work through the object's history until
00901 bool ClusterSequence::object_in_jet(const PseudoJet & object, 
00902                                     const PseudoJet & jet) const {
00903 
00904   // make sure the object conceivably belongs to this clustering
00905   // sequence
00906   assert(contains(object) && contains(jet));
00907 
00908   const PseudoJet * this_object = &object;
00909   const PseudoJet * childp;
00910   while(true) {
00911     if (this_object->cluster_hist_index() == jet.cluster_hist_index()) {
00912       return true;
00913     } else if (has_child(*this_object, childp)) {
00914       this_object = childp;
00915     } else {
00916       return false;
00917     }
00918   }
00919 }
00920 
00921 //----------------------------------------------------------------------
00922 /// if the jet has parents in the clustering, it returns true
00923 /// and sets parent1 and parent2 equal to them.
00924 ///
00925 /// if it has no parents it returns false and sets parent1 and
00926 /// parent2 to zero
00927 bool ClusterSequence::has_parents(const PseudoJet & jet, PseudoJet & parent1, 
00928                               PseudoJet & parent2) const {
00929 
00930   const history_element & hist = _history[jet.cluster_hist_index()];
00931 
00932   // make sure we do not run into any unexpected situations --
00933   // i.e. both parents valid, or neither
00934   assert ((hist.parent1 >= 0 && hist.parent2 >= 0) || 
00935           (hist.parent1 < 0 && hist.parent2 < 0));
00936 
00937   if (hist.parent1 < 0) {
00938     parent1 = PseudoJet(0.0,0.0,0.0,0.0);
00939     parent2 = parent1;
00940     return false;
00941   } else {
00942     parent1 = _jets[_history[hist.parent1].jetp_index];
00943     parent2 = _jets[_history[hist.parent2].jetp_index];
00944     // order the parents in decreasing pt
00945     if (parent1.perp2() < parent2.perp2()) std::swap(parent1,parent2);
00946     return true;
00947   }
00948 }
00949 
00950 //----------------------------------------------------------------------
00951 /// if the jet has a child then return true and give the child jet
00952 /// otherwise return false and set the child to zero
00953 bool ClusterSequence::has_child(const PseudoJet & jet, PseudoJet & child) const {
00954 
00955   //const history_element & hist = _history[jet.cluster_hist_index()];
00956   //
00957   //if (hist.child >= 0) {
00958   //  child = _jets[_history[hist.child].jetp_index];
00959   //  return true;
00960   //} else {
00961   //  child = PseudoJet(0.0,0.0,0.0,0.0);
00962   //  return false;
00963   //}
00964   const PseudoJet * childp;
00965   bool res = has_child(jet, childp);
00966   if (res) {
00967     child = *childp;
00968     return true;
00969   } else {
00970     child = PseudoJet(0.0,0.0,0.0,0.0);
00971     return false;
00972   }
00973 }
00974 
00975 bool ClusterSequence::has_child(const PseudoJet & jet, const PseudoJet * & childp) const {
00976 
00977   const history_element & hist = _history[jet.cluster_hist_index()];
00978 
00979   // check that this jet has a child and that the child corresponds to
00980   // a true jet [RETHINK-IF-CHANGE-NUMBERING: what is the right
00981   // behaviour if the child is the same jet but made inclusive...?]
00982   if (hist.child >= 0 && _history[hist.child].jetp_index >= 0) {
00983     childp = &(_jets[_history[hist.child].jetp_index]);
00984     return true;
00985   } else {
00986     childp = NULL;
00987     return false;
00988   }
00989 }
00990 
00991 
00992 //----------------------------------------------------------------------
00993 /// if this jet has a child (and so a partner) return true
00994 /// and give the partner, otherwise return false and set the
00995 /// partner to zero
00996 bool ClusterSequence::has_partner(const PseudoJet & jet, 
00997                               PseudoJet & partner) const {
00998 
00999   const history_element & hist = _history[jet.cluster_hist_index()];
01000 
01001   // make sure we have a child and that the child does not correspond
01002   // to a clustering with the beam (or some other invalid quantity)
01003   if (hist.child >= 0 && _history[hist.child].parent2 >= 0) {
01004     const history_element & child_hist = _history[hist.child];
01005     if (child_hist.parent1 == jet.cluster_hist_index()) {
01006       // partner will be child's parent2 -- for iB clustering
01007       // parent2 will not be valid
01008       partner = _jets[_history[child_hist.parent2].jetp_index];
01009     } else {
01010       // partner will be child's parent1
01011       partner = _jets[_history[child_hist.parent1].jetp_index];
01012     }
01013     return true;
01014   } else {
01015     partner = PseudoJet(0.0,0.0,0.0,0.0);
01016     return false;
01017   }
01018 }
01019 
01020 
01021 //----------------------------------------------------------------------
01022 // return a vector of the particles that make up a jet
01023 vector<PseudoJet> ClusterSequence::constituents (const PseudoJet & jet) const {
01024   vector<PseudoJet> subjets;
01025   add_constituents(jet, subjets);
01026   return subjets;
01027 }
01028 
01029 //----------------------------------------------------------------------
01030 /// output the supplied vector of jets in a format that can be read
01031 /// by an appropriate root script; the format is:
01032 /// jet-n jet-px jet-py jet-pz jet-E 
01033 ///   particle-n particle-rap particle-phi particle-pt
01034 ///   particle-n particle-rap particle-phi particle-pt
01035 ///   ...
01036 /// #END
01037 /// ... [i.e. above repeated]
01038 void ClusterSequence::print_jets_for_root(const std::vector<PseudoJet> & jets, 
01039                                           ostream & ostr) const {
01040   for (unsigned i = 0; i < jets.size(); i++) {
01041     ostr << i  << " "
01042          << jets[i].px() << " "
01043          << jets[i].py() << " "
01044          << jets[i].pz() << " "
01045          << jets[i].E() << endl;
01046     vector<PseudoJet> cst = constituents(jets[i]);
01047     for (unsigned j = 0; j < cst.size() ; j++) {
01048       ostr << " " << j << " "
01049            << cst[j].rap() << " "
01050            << cst[j].phi() << " "
01051            << cst[j].perp() << endl;
01052     }
01053     ostr << "#END" << endl;
01054   }
01055 }
01056 
01057 void ClusterSequence::print_jets_for_root(const std::vector<PseudoJet> & jets, 
01058                                           const std::string & filename,
01059                                           const std::string & comment ) const {
01060   std::ofstream ostr(filename.c_str());
01061   if (comment != "") ostr << "# " << comment << endl;
01062   print_jets_for_root(jets, ostr);
01063 }
01064 
01065 
01066 // Not yet. Perhaps in a future release
01067 // //----------------------------------------------------------------------
01068 // // print out all inclusive jets with pt > ptmin
01069 // void ClusterSequence::print_jets (const double & ptmin) const{
01070 //     vector<PseudoJet> jets = sorted_by_pt(inclusive_jets(ptmin));
01071 // 
01072 //     for (size_t j = 0; j < jets.size(); j++) {
01073 //        printf("%5u %7.3f %7.3f %9.3f\n",
01074 //        j,jets[j].rap(),jets[j].phi(),jets[j].perp());
01075 //     }
01076 // }
01077 
01078 //----------------------------------------------------------------------
01079 /// returns a vector of size n_particles() which indicates, for 
01080 /// each of the initial particles (in the order in which they were
01081 /// supplied), which of the supplied jets it belongs to; if it does
01082 /// not belong to any of the supplied jets, the index is set to -1;
01083 vector<int> ClusterSequence::particle_jet_indices(
01084                         const vector<PseudoJet> & jets) const {
01085 
01086   vector<int> indices(n_particles());
01087 
01088   // first label all particles as not belonging to any jets
01089   for (unsigned ipart = 0; ipart < n_particles(); ipart++) 
01090     indices[ipart] = -1;
01091 
01092   // then for each of the jets relabel its consituents as belonging to
01093   // that jet
01094   for (unsigned ijet = 0; ijet < jets.size(); ijet++) {
01095 
01096     vector<PseudoJet> jet_constituents(constituents(jets[ijet]));
01097 
01098     for (unsigned ip = 0; ip < jet_constituents.size(); ip++) {
01099       // a safe (if slightly redundant) way of getting the particle
01100       // index (for initial particles it is actually safe to assume
01101       // ipart=iclust).
01102       unsigned iclust = jet_constituents[ip].cluster_hist_index();
01103       unsigned ipart = history()[iclust].jetp_index;
01104       indices[ipart] = ijet;
01105     }
01106   }
01107 
01108   return indices;
01109 }
01110 
01111 
01112 //----------------------------------------------------------------------
01113 // recursive routine that adds on constituents of jet to the subjet_vector
01114 void ClusterSequence::add_constituents (
01115            const PseudoJet & jet, vector<PseudoJet> & subjet_vector) const {
01116   // find out position in cluster history
01117   int i = jet.cluster_hist_index();
01118   int parent1 = _history[i].parent1;
01119   int parent2 = _history[i].parent2;
01120 
01121   if (parent1 == InexistentParent) {
01122     // It is an original particle (labelled by its parent having value
01123     // InexistentParent), therefore add it on to the subjet vector
01124     // Note: we add the initial particle and not simply 'jet' so that
01125     //       calling add_constituents with a subtracted jet containing
01126     //       only one particle will work.
01127     subjet_vector.push_back(_jets[i]);
01128     return;
01129   } 
01130 
01131   // add parent 1
01132   add_constituents(_jets[_history[parent1].jetp_index], subjet_vector);
01133 
01134   // see if parent2 is a real jet; if it is then add its constituents
01135   if (parent2 != BeamJet) {
01136     add_constituents(_jets[_history[parent2].jetp_index], subjet_vector);
01137   }
01138 }
01139 
01140 
01141 
01142 //----------------------------------------------------------------------
01143 // initialise the history in a standard way
01144 void ClusterSequence::_add_step_to_history (
01145                const int & step_number, const int & parent1, 
01146                const int & parent2, const int & jetp_index,
01147                const double & dij) {
01148 
01149   history_element element;
01150   element.parent1 = parent1;
01151   element.parent2 = parent2;
01152   element.jetp_index = jetp_index;
01153   element.child = Invalid;
01154   element.dij   = dij;
01155   element.max_dij_so_far = max(dij,_history[_history.size()-1].max_dij_so_far);
01156   _history.push_back(element);
01157 
01158   int local_step = _history.size()-1;
01159   assert(local_step == step_number);
01160 
01161   assert(parent1 >= 0);
01162   _history[parent1].child = local_step;
01163   if (parent2 >= 0) {_history[parent2].child = local_step;}
01164 
01165   // get cross-referencing right from PseudoJets
01166   if (jetp_index != Invalid) {
01167     assert(jetp_index >= 0);
01168     //cout << _jets.size() <<" "<<jetp_index<<"\n";
01169     _jets[jetp_index].set_cluster_hist_index(local_step);
01170     _set_structure_shared_ptr(_jets[jetp_index]);
01171   }
01172 
01173   if (_writeout_combinations) {
01174     cout << local_step << ": " 
01175          << parent1 << " with " << parent2
01176          << "; y = "<< dij<<endl;
01177   }
01178 
01179 }
01180 
01181 
01182 
01183 
01184 //======================================================================
01185 // Return an order in which to read the history such that _history[order[i]] 
01186 // will always correspond to the same set of consituent particles if 
01187 // two branching histories are equivalent in terms of the particles
01188 // contained in any given pseudojet.
01189 vector<int> ClusterSequence::unique_history_order() const {
01190 
01191   // first construct an array that will tell us the lowest constituent
01192   // of a given jet -- this will always be one of the original
01193   // particles, whose order is well defined and so will help us to
01194   // follow the tree in a unique manner.
01195   valarray<int> lowest_constituent(_history.size());
01196   int hist_n = _history.size();
01197   lowest_constituent = hist_n; // give it a large number
01198   for (int i = 0; i < hist_n; i++) {
01199     // sets things up for the initial partons
01200     lowest_constituent[i] = min(lowest_constituent[i],i); 
01201     // propagates them through to the children of this parton
01202     if (_history[i].child > 0) lowest_constituent[_history[i].child] 
01203       = min(lowest_constituent[_history[i].child],lowest_constituent[i]);
01204   }
01205 
01206   // establish an array for what we have and have not extracted so far
01207   valarray<bool> extracted(_history.size()); extracted = false;
01208   vector<int> unique_tree;
01209   unique_tree.reserve(_history.size());
01210 
01211   // now work our way through the tree
01212   for (unsigned i = 0; i < n_particles(); i++) {
01213     if (!extracted[i]) {
01214       unique_tree.push_back(i);
01215       extracted[i] = true;
01216       _extract_tree_children(i, extracted, lowest_constituent, unique_tree);
01217     }
01218   }
01219 
01220   return unique_tree;
01221 }
01222 
01223 //======================================================================
01224 // helper for unique_history_order
01225 void ClusterSequence::_extract_tree_children(
01226        int position, 
01227        valarray<bool> & extracted, 
01228        const valarray<int> & lowest_constituent,
01229        vector<int> & unique_tree) const {
01230   if (!extracted[position]) {
01231     // that means we may have unidentified parents around, so go and
01232     // collect them (extracted[position]) will then be made true)
01233     _extract_tree_parents(position,extracted,lowest_constituent,unique_tree);
01234   } 
01235   
01236   // now look after the children...
01237   int child = _history[position].child;
01238   if (child  >= 0) _extract_tree_children(child,extracted,lowest_constituent,unique_tree);
01239 }
01240 
01241 
01242 //======================================================================
01243 // return the list of unclustered particles
01244 vector<PseudoJet> ClusterSequence::unclustered_particles() const {
01245   vector<PseudoJet> unclustered;
01246   for (unsigned i = 0; i < n_particles() ; i++) {
01247     if (_history[i].child == Invalid) 
01248       unclustered.push_back(_jets[_history[i].jetp_index]);
01249   }
01250   return unclustered;
01251 }
01252 
01253 
01254 
01255 //----------------------------------------------------------------------
01256 // returns true if the cluster sequence contains this jet (i.e. jet's
01257 // structure is this cluster sequence's and the cluster history index
01258 // is in a consistent range)
01259 bool ClusterSequence::contains(const PseudoJet & jet) const {
01260   return jet.cluster_hist_index() >= 0 
01261     &&   jet.cluster_hist_index() < int(_history.size())
01262     &&   jet.has_valid_cluster_sequence()
01263     &&   jet.associated_cluster_sequence() == this;
01264 }
01265 
01266 
01267 
01268 //======================================================================
01269 // helper for unique_history_order
01270 void ClusterSequence::_extract_tree_parents(
01271        int position, 
01272        valarray<bool> & extracted, 
01273        const valarray<int> & lowest_constituent,
01274        vector<int> & unique_tree) const {
01275 
01276   if (!extracted[position]) {
01277     int parent1 = _history[position].parent1;
01278     int parent2 = _history[position].parent2;
01279     // where relevant order parents so that we will first treat the
01280     // one containing the smaller "lowest_constituent"
01281     if (parent1 >= 0 && parent2 >= 0) {
01282       if (lowest_constituent[parent1] > lowest_constituent[parent2]) 
01283         std::swap(parent1, parent2);
01284     }
01285     // then actually run through the parents to extract the constituents...
01286     if (parent1 >= 0 && !extracted[parent1]) 
01287       _extract_tree_parents(parent1,extracted,lowest_constituent,unique_tree);
01288     if (parent2 >= 0 && !extracted[parent2]) 
01289       _extract_tree_parents(parent2,extracted,lowest_constituent,unique_tree);
01290     // finally declare this position to be accounted for and push it
01291     // onto our list.
01292     unique_tree.push_back(position);
01293     extracted[position] = true;
01294   }
01295 }
01296 
01297 
01298 //======================================================================
01299 /// carries out the bookkeeping associated with the step of recombining
01300 /// jet_i and jet_j (assuming a distance dij) and returns the index
01301 /// of the recombined jet, newjet_k.
01302 void ClusterSequence::_do_ij_recombination_step(
01303                                const int & jet_i, const int & jet_j, 
01304                                const double & dij, 
01305                                int & newjet_k) {
01306 
01307   // Create the new jet by recombining the first two.
01308   //
01309   // For efficiency reasons, use a ctr that initialises only the
01310   // shared pointers, since the rest of the info will anyway be dealt
01311   // with by the recombiner.
01312   PseudoJet newjet(false); 
01313   _jet_def.recombiner()->recombine(_jets[jet_i], _jets[jet_j], newjet);
01314   _jets.push_back(newjet);
01315   // original version...
01316   //_jets.push_back(_jets[jet_i] + _jets[jet_j]);
01317 
01318   // get its index
01319   newjet_k = _jets.size()-1;
01320 
01321   // get history index
01322   int newstep_k = _history.size();
01323   // and provide jet with the info
01324   _jets[newjet_k].set_cluster_hist_index(newstep_k);
01325 
01326   // finally sort out the history 
01327   int hist_i = _jets[jet_i].cluster_hist_index();
01328   int hist_j = _jets[jet_j].cluster_hist_index();
01329 
01330   _add_step_to_history(newstep_k, min(hist_i, hist_j), max(hist_i,hist_j),
01331                        newjet_k, dij);
01332 
01333 }
01334 
01335 
01336 //======================================================================
01337 /// carries out the bookkeeping associated with the step of recombining
01338 /// jet_i with the beam
01339 void ClusterSequence::_do_iB_recombination_step(
01340                                   const int & jet_i, const double & diB) {
01341   // get history index
01342   int newstep_k = _history.size();
01343 
01344   // recombine the jet with the beam
01345   _add_step_to_history(newstep_k,_jets[jet_i].cluster_hist_index(),BeamJet,
01346                        Invalid, diB);
01347 
01348 }
01349 
01350 
01351 // make sure the static member _changed_strategy_warning is defined. 
01352 LimitedWarning ClusterSequence::_changed_strategy_warning;
01353 
01354 
01355 //----------------------------------------------------------------------
01356 void ClusterSequence::_set_structure_shared_ptr(PseudoJet & j) {
01357   j.set_structure_shared_ptr(_structure_shared_ptr);
01358   // record the use count of the structure shared point to help
01359   // in case we want to ask the CS to handle its own memory
01360   _update_structure_use_count();
01361 }
01362 
01363 
01364 //----------------------------------------------------------------------
01365 void ClusterSequence::_update_structure_use_count() {
01366   // record the use count of the structure shared point to help
01367   // in case we want to ask the CS to handle its own memory
01368   _structure_use_count_after_construction = _structure_shared_ptr.use_count();
01369 }
01370 
01371 //----------------------------------------------------------------------
01372 /// by calling this routine you tell the ClusterSequence to delete
01373 /// itself when all the Pseudojets associated with it have gone out
01374 /// of scope. 
01375 void ClusterSequence::delete_self_when_unused() {
01376   // the trick we use to handle this is to modify the use count; 
01377   // that way the structure will be deleted when there are no external
01378   // objects left associated the CS and the structure's destructor will then
01379   // look after deleting the cluster sequence
01380   
01381   // first make sure that there is at least one other object
01382   // associated with the CS
01383   int new_count = _structure_shared_ptr.use_count() - _structure_use_count_after_construction;
01384   if (new_count <= 0) {
01385     throw Error("delete_self_when_unused may only be called if at least one object outside the CS (e.g. a jet) is already associated with the CS");
01386   }
01387 
01388   _structure_shared_ptr.set_count(new_count);
01389   _deletes_self_when_unused = true;
01390 }
01391 
01392 
01393 FASTJET_END_NAMESPACE
01394 
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends