00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031 #include "fastjet/Error.hh"
00032 #include "fastjet/PseudoJet.hh"
00033 #include "fastjet/ClusterSequence.hh"
00034 #include "fastjet/version.hh"
00035 #include<iostream>
00036 #include<sstream>
00037 #include<fstream>
00038 #include<cmath>
00039 #include<cstdlib>
00040 #include<cassert>
00041 #include<string>
00042 #include<set>
00043
00044 FASTJET_BEGIN_NAMESPACE
00045
00046 using namespace std;
00047
00049 JetAlgorithm ClusterSequence::_default_jet_algorithm = kt_algorithm;
00050
00051
00052
00053
00054 ClusterSequence::~ClusterSequence () {
00055
00056
00057 if (_wrapper_to_this()) _wrapper_to_this->_cs = NULL;
00058
00059 }
00060
00061
00062 void ClusterSequence::_initialise_and_run (
00063 const double & R,
00064 const Strategy & strategy,
00065 const bool & writeout_combinations) {
00066
00067 JetDefinition jet_def(_default_jet_algorithm, R, strategy);
00068 _initialise_and_run(jet_def, writeout_combinations);
00069 }
00070
00071
00072
00073 void ClusterSequence::_initialise_and_run (
00074 const JetDefinition & jet_def,
00075 const bool & writeout_combinations) {
00076
00077
00078 _decant_options(jet_def, writeout_combinations);
00079
00080
00081
00082 _fill_initial_history();
00083
00084
00085 if (n_particles() == 0) return;
00086
00087
00088 if (_jet_algorithm == plugin_algorithm) {
00089
00090 _plugin_activated = true;
00091
00092 _jet_def.plugin()->run_clustering( (*this) );
00093 _plugin_activated = false;
00094 return;
00095 } else if (_jet_algorithm == ee_kt_algorithm ||
00096 _jet_algorithm == ee_genkt_algorithm) {
00097
00098 _strategy = N2Plain;
00099 if (_jet_algorithm == ee_kt_algorithm) {
00100
00101
00102
00103 assert(_Rparam > 2.0);
00104
00105
00106
00107 _invR2 = 1.0;
00108 } else {
00109
00110
00111
00112
00113 if (_Rparam > pi) {
00114
00115
00116
00117 _R2 = 2 * ( 3.0 + cos(_Rparam) );
00118 } else {
00119 _R2 = 2 * ( 1.0 - cos(_Rparam) );
00120 }
00121 _invR2 = 1.0/_R2;
00122 }
00123
00124 _simple_N2_cluster_EEBriefJet();
00125 return;
00126 }
00127
00128
00129
00130
00131
00132
00133
00134 if (_strategy == Best) {
00135 int N = _jets.size();
00136 if (N > 6200/pow(_Rparam,2.0)
00137 && jet_def.jet_algorithm() == cambridge_algorithm) {
00138 _strategy = NlnNCam;}
00139 else
00140 #ifndef DROP_CGAL
00141 if ((N > 16000/pow(_Rparam,1.15) && jet_def.jet_algorithm() != antikt_algorithm)
00142 || N > 35000/pow(_Rparam,1.15)) {
00143 _strategy = NlnN; }
00144 else
00145 #endif // DROP_CGAL
00146 if (N > 450) {
00147 _strategy = N2MinHeapTiled;
00148 }
00149 else if (N > 55*max(0.5,min(1.0,_Rparam))) {
00150 _strategy = N2Tiled;
00151 } else {
00152 _strategy = N2Plain;
00153 }
00154 }
00155
00156
00157
00158
00159
00160 if (_Rparam >= twopi) {
00161 if ( _strategy == NlnN
00162 || _strategy == NlnN3pi
00163 || _strategy == NlnNCam
00164 || _strategy == NlnNCam2pi2R
00165 || _strategy == NlnNCam4pi) {
00166 #ifdef DROP_CGAL
00167 _strategy = N2MinHeapTiled;
00168 #else
00169 _strategy = NlnN4pi;
00170 #endif
00171 }
00172 if (jet_def.strategy() != Best && _strategy != jet_def.strategy()) {
00173 ostringstream oss;
00174 oss << "Cluster strategy " << strategy_string(jet_def.strategy())
00175 << " automatically changed to " << strategy_string()
00176 << " because the former is not supported for R = " << _Rparam
00177 << " >= 2pi";
00178 _changed_strategy_warning.warn(oss.str());
00179 }
00180 }
00181
00182
00183
00184 if (_strategy == NlnN || _strategy == NlnN3pi
00185 || _strategy == NlnN4pi ) {
00186 this->_delaunay_cluster();
00187 } else if (_strategy == N3Dumb ) {
00188 this->_really_dumb_cluster();
00189 } else if (_strategy == N2Tiled) {
00190 this->_faster_tiled_N2_cluster();
00191 } else if (_strategy == N2PoorTiled) {
00192 this->_tiled_N2_cluster();
00193 } else if (_strategy == N2Plain) {
00194
00195
00196 this->_simple_N2_cluster_BriefJet();
00197 } else if (_strategy == N2MinHeapTiled) {
00198 this->_minheap_faster_tiled_N2_cluster();
00199 } else if (_strategy == NlnNCam4pi) {
00200 this->_CP2DChan_cluster();
00201 } else if (_strategy == NlnNCam2pi2R) {
00202 this->_CP2DChan_cluster_2pi2R();
00203 } else if (_strategy == NlnNCam) {
00204 this->_CP2DChan_cluster_2piMultD();
00205 } else {
00206 ostringstream err;
00207 err << "Unrecognised value for strategy: "<<_strategy;
00208 throw Error(err.str());
00209
00210 }
00211 }
00212
00213
00214
00215 bool ClusterSequence::_first_time = true;
00216 int ClusterSequence::_n_exclusive_warnings = 0;
00217
00218
00219
00220
00221 string fastjet_version_string() {
00222 return "FastJet version "+string(fastjet_version);
00223 }
00224
00225
00226
00227
00228 void ClusterSequence::_print_banner() {
00229
00230 if (!_first_time) {return;}
00231 _first_time = false;
00232
00233
00234
00235
00236 cout << "#--------------------------------------------------------------------------\n";
00237 cout << "# FastJet release " << fastjet_version << endl;
00238 cout << "# Written by M. Cacciari, G.P. Salam and G. Soyez \n";
00239 cout << "# http://www.fastjet.fr \n";
00240 cout << "# \n";
00241 cout << "# Longitudinally invariant Kt, anti-Kt, and inclusive Cambridge/Aachen \n";
00242 cout << "# clustering using fast geometric algorithms, with area measures and optional\n";
00243 cout << "# external jet-finder plugins. \n";
00244 cout << "# Please cite Phys. Lett. B641 (2006) [hep-ph/0512210] if you use this code.\n";
00245 cout << "# \n";
00246 cout << "# This package uses T.Chan's closest pair algorithm, Proc.13th ACM-SIAM \n";
00247 cout << "# Symp. Discr. Alg, p.472 (2002), S.Fortune's Voronoi algorithm and code " ;
00248 #ifndef DROP_CGAL
00249 cout << endl << "# and CGAL: http://www.cgal.org/";
00250 #endif // DROP_CGAL
00251 cout << ".\n";
00252 cout << "#-------------------------------------------------------------------------\n";
00253
00254 cout.flush();
00255 }
00256
00257
00258
00259 void ClusterSequence::_decant_options(const JetDefinition & jet_def,
00260 const bool & writeout_combinations) {
00261
00262 _print_banner();
00263
00264
00265 _jet_def = jet_def;
00266
00267 _writeout_combinations = writeout_combinations;
00268 _jet_algorithm = jet_def.jet_algorithm();
00269 _Rparam = jet_def.R(); _R2 = _Rparam*_Rparam; _invR2 = 1.0/_R2;
00270 _strategy = jet_def.strategy();
00271
00272
00273 _plugin_activated = false;
00274
00275
00276 _wrapper_to_this.reset(new ClusterSequenceWrapper(this));
00277 }
00278
00279
00280
00281
00282 void ClusterSequence::_fill_initial_history () {
00283
00284
00285
00286
00287 _jets.reserve(_jets.size()*2);
00288 _history.reserve(_jets.size()*2);
00289
00290 _Qtot = 0;
00291
00292 for (int i = 0; i < static_cast<int>(_jets.size()) ; i++) {
00293 history_element element;
00294 element.parent1 = InexistentParent;
00295 element.parent2 = InexistentParent;
00296 element.child = Invalid;
00297 element.jetp_index = i;
00298 element.dij = 0.0;
00299 element.max_dij_so_far = 0.0;
00300
00301 _history.push_back(element);
00302
00303
00304 _jet_def.recombiner()->preprocess(_jets[i]);
00305
00306
00307 _jets[i].set_cluster_hist_index(i);
00308 _jets[i].set_associated_csw(_wrapper_to_this);
00309
00310
00311 _Qtot += _jets[i].E();
00312 }
00313 _initial_n = _jets.size();
00314 }
00315
00316
00317
00318
00319
00320 string ClusterSequence::strategy_string (Strategy strategy_in) const {
00321 string strategy;
00322 switch(strategy_in) {
00323 case NlnN:
00324 strategy = "NlnN"; break;
00325 case NlnN3pi:
00326 strategy = "NlnN3pi"; break;
00327 case NlnN4pi:
00328 strategy = "NlnN4pi"; break;
00329 case N2Plain:
00330 strategy = "N2Plain"; break;
00331 case N2Tiled:
00332 strategy = "N2Tiled"; break;
00333 case N2MinHeapTiled:
00334 strategy = "N2MinHeapTiled"; break;
00335 case N2PoorTiled:
00336 strategy = "N2PoorTiled"; break;
00337 case N3Dumb:
00338 strategy = "N3Dumb"; break;
00339 case NlnNCam4pi:
00340 strategy = "NlnNCam4pi"; break;
00341 case NlnNCam2pi2R:
00342 strategy = "NlnNCam2pi2R"; break;
00343 case NlnNCam:
00344 strategy = "NlnNCam"; break;
00345 case plugin_strategy:
00346 strategy = "plugin strategy"; break;
00347 default:
00348 strategy = "Unrecognized";
00349 }
00350 return strategy;
00351 }
00352
00353
00354 double ClusterSequence::jet_scale_for_algorithm(
00355 const PseudoJet & jet) const {
00356 if (_jet_algorithm == kt_algorithm) {return jet.kt2();}
00357 else if (_jet_algorithm == cambridge_algorithm) {return 1.0;}
00358 else if (_jet_algorithm == antikt_algorithm) {
00359 double kt2=jet.kt2();
00360 return kt2 > 1e-300 ? 1.0/kt2 : 1e300;
00361 } else if (_jet_algorithm == genkt_algorithm) {
00362 double kt2 = jet.kt2();
00363 double p = jet_def().extra_param();
00364 if (p <= 0 && kt2 < 1e-300) kt2 = 1e-300;
00365 return pow(kt2, p);
00366 } else if (_jet_algorithm == cambridge_for_passive_algorithm) {
00367 double kt2 = jet.kt2();
00368 double lim = _jet_def.extra_param();
00369 if (kt2 < lim*lim && kt2 != 0.0) {
00370 return 1.0/kt2;
00371 } else {return 1.0;}
00372 } else {throw Error("Unrecognised jet algorithm");}
00373 }
00374
00375
00376
00380 void ClusterSequence::transfer_from_sequence(ClusterSequence & from_seq, bool transfer_ownership) {
00381
00382
00383 _jet_def = from_seq._jet_def ;
00384 _writeout_combinations = from_seq._writeout_combinations ;
00385 _initial_n = from_seq._initial_n ;
00386 _Rparam = from_seq._Rparam ;
00387 _R2 = from_seq._R2 ;
00388 _invR2 = from_seq._invR2 ;
00389 _strategy = from_seq._strategy ;
00390 _jet_algorithm = from_seq._jet_algorithm ;
00391 _plugin_activated = from_seq._plugin_activated ;
00392
00393
00394 _jets = from_seq._jets;
00395 _history = from_seq._history;
00396
00397 _extras = from_seq._extras;
00398
00399
00400 if (transfer_ownership){
00401
00402 if (! _wrapper_to_this()) _wrapper_to_this.reset(new ClusterSequenceWrapper(this));
00403
00404 for (vector<PseudoJet>::iterator jit = _jets.begin(); jit != _jets.end(); jit++)
00405 jit->set_associated_csw(_wrapper_to_this);
00406 }
00407 }
00408
00409
00410
00411
00412 void ClusterSequence::plugin_record_ij_recombination(
00413 int jet_i, int jet_j, double dij,
00414 const PseudoJet & newjet, int & newjet_k) {
00415
00416 plugin_record_ij_recombination(jet_i, jet_j, dij, newjet_k);
00417
00418
00419 int tmp_index = _jets[newjet_k].cluster_hist_index();
00420 _jets[newjet_k] = newjet;
00421 _jets[newjet_k].set_cluster_hist_index(tmp_index);
00422 _jets[newjet_k].set_associated_csw(_wrapper_to_this);
00423
00424 }
00425
00426
00427
00428
00429 vector<PseudoJet> ClusterSequence::inclusive_jets (const double & ptmin) const{
00430 double dcut = ptmin*ptmin;
00431 int i = _history.size() - 1;
00432 vector<PseudoJet> jets;
00433 if (_jet_algorithm == kt_algorithm) {
00434 while (i >= 0) {
00435
00436
00437
00438 if (_history[i].max_dij_so_far < dcut) {break;}
00439 if (_history[i].parent2 == BeamJet && _history[i].dij >= dcut) {
00440
00441 int parent1 = _history[i].parent1;
00442 jets.push_back(_jets[_history[parent1].jetp_index]);}
00443 i--;
00444 }
00445 } else if (_jet_algorithm == cambridge_algorithm) {
00446 while (i >= 0) {
00447
00448
00449
00450 if (_history[i].parent2 != BeamJet) {break;}
00451 int parent1 = _history[i].parent1;
00452 const PseudoJet & jet = _jets[_history[parent1].jetp_index];
00453 if (jet.perp2() >= dcut) {jets.push_back(jet);}
00454 i--;
00455 }
00456 } else if (_jet_algorithm == plugin_algorithm
00457 || _jet_algorithm == ee_kt_algorithm
00458 || _jet_algorithm == antikt_algorithm
00459 || _jet_algorithm == genkt_algorithm
00460 || _jet_algorithm == ee_genkt_algorithm
00461 || _jet_algorithm == cambridge_for_passive_algorithm) {
00462
00463
00464
00465 while (i >= 0) {
00466 if (_history[i].parent2 == BeamJet) {
00467 int parent1 = _history[i].parent1;
00468 const PseudoJet & jet = _jets[_history[parent1].jetp_index];
00469 if (jet.perp2() >= dcut) {jets.push_back(jet);}
00470 }
00471 i--;
00472 }
00473 } else {throw Error("cs::inclusive_jets(...): Unrecognized jet algorithm");}
00474 return jets;
00475 }
00476
00477
00478
00479
00480
00481 int ClusterSequence::n_exclusive_jets (const double & dcut) const {
00482
00483
00484
00485 int i = _history.size() - 1;
00486 while (i >= 0) {
00487 if (_history[i].max_dij_so_far <= dcut) {break;}
00488 i--;
00489 }
00490 int stop_point = i + 1;
00491
00492
00493 int njets = 2*_initial_n - stop_point;
00494 return njets;
00495 }
00496
00497
00498
00499
00500 vector<PseudoJet> ClusterSequence::exclusive_jets (const double & dcut) const {
00501 int njets = n_exclusive_jets(dcut);
00502 return exclusive_jets(njets);
00503 }
00504
00505
00506
00507
00508 vector<PseudoJet> ClusterSequence::exclusive_jets (const int & njets) const {
00509
00510
00511
00512 assert (njets <= _initial_n);
00513
00514
00515
00516
00517
00518
00519 if (( _jet_def.jet_algorithm() != kt_algorithm) &&
00520 ( _jet_def.jet_algorithm() != cambridge_algorithm) &&
00521 ( _jet_def.jet_algorithm() != ee_kt_algorithm) &&
00522 (((_jet_def.jet_algorithm() != genkt_algorithm) &&
00523 (_jet_def.jet_algorithm() != ee_genkt_algorithm)) ||
00524 (_jet_def.extra_param() <0)) &&
00525 ((_jet_def.jet_algorithm() != plugin_algorithm) ||
00526 (!_jet_def.plugin()->exclusive_sequence_meaningful())) &&
00527 (_n_exclusive_warnings < 5)) {
00528 _n_exclusive_warnings++;
00529 cerr << "FastJet WARNING: dcut and exclusive jets for jet-finders other than kt should be interpreted with care." << endl;
00530 }
00531
00532
00533
00534
00535
00536 int stop_point = 2*_initial_n - njets;
00537
00538
00539
00540 if (2*_initial_n != static_cast<int>(_history.size())) {
00541 ostringstream err;
00542 err << "2*_initial_n != _history.size() -- this endangers internal assumptions!\n";
00543 throw Error(err.str());
00544
00545 }
00546
00547
00548
00549
00550
00551
00552 vector<PseudoJet> jets;
00553 for (unsigned int i = stop_point; i < _history.size(); i++) {
00554 int parent1 = _history[i].parent1;
00555 if (parent1 < stop_point) {
00556 jets.push_back(_jets[_history[parent1].jetp_index]);
00557 }
00558 int parent2 = _history[i].parent2;
00559 if (parent2 < stop_point && parent2 > 0) {
00560 jets.push_back(_jets[_history[parent2].jetp_index]);
00561 }
00562
00563 }
00564
00565
00566 if (static_cast<int>(jets.size()) != njets) {
00567 ostringstream err;
00568 err << "ClusterSequence::exclusive_jets: size of returned vector ("
00569 <<jets.size()<<") does not coincide with requested number of jets ("
00570 <<njets<<")";
00571 throw Error(err.str());
00572 }
00573
00574 return jets;
00575 }
00576
00577
00580 double ClusterSequence::exclusive_dmerge (const int & njets) const {
00581 assert(njets >= 0);
00582 if (njets >= _initial_n) {return 0.0;}
00583 return _history[2*_initial_n-njets-1].dij;
00584 }
00585
00586
00587
00592 double ClusterSequence::exclusive_dmerge_max (const int & njets) const {
00593 assert(njets >= 0);
00594 if (njets >= _initial_n) {return 0.0;}
00595 return _history[2*_initial_n-njets-1].max_dij_so_far;
00596 }
00597
00598
00599
00603 std::vector<PseudoJet> ClusterSequence::exclusive_subjets
00604 (const PseudoJet & jet, const double & dcut) const {
00605
00606 set<const history_element*> subhist;
00607
00608
00609
00610 get_subhist_set(subhist, jet, dcut, 0);
00611
00612
00613 vector<PseudoJet> subjets;
00614 subjets.reserve(subhist.size());
00615 for (set<const history_element*>::iterator elem = subhist.begin();
00616 elem != subhist.end(); elem++) {
00617 subjets.push_back(_jets[(*elem)->jetp_index]);
00618 }
00619 return subjets;
00620 }
00621
00622
00626 int ClusterSequence::n_exclusive_subjets(const PseudoJet & jet,
00627 const double & dcut) const {
00628 set<const history_element*> subhist;
00629
00630
00631 get_subhist_set(subhist, jet, dcut, 0);
00632 return subhist.size();
00633 }
00634
00635
00639 std::vector<PseudoJet> ClusterSequence::exclusive_subjets
00640 (const PseudoJet & jet, int n) const {
00641
00642 set<const history_element*> subhist;
00643
00644
00645
00646 get_subhist_set(subhist, jet, -1.0, n);
00647
00648
00649 vector<PseudoJet> subjets;
00650 subjets.reserve(subhist.size());
00651 for (set<const history_element*>::iterator elem = subhist.begin();
00652 elem != subhist.end(); elem++) {
00653 subjets.push_back(_jets[(*elem)->jetp_index]);
00654 }
00655 return subjets;
00656 }
00657
00658
00659
00664 double ClusterSequence::exclusive_subdmerge(const PseudoJet & jet, int nsub) const {
00665 set<const history_element*> subhist;
00666
00667
00668
00669 get_subhist_set(subhist, jet, -1.0, nsub);
00670
00671 set<const history_element*>::iterator highest = subhist.end();
00672 highest--;
00675 return (*highest)->dij;
00676 }
00677
00678
00679
00685 double ClusterSequence::exclusive_subdmerge_max(const PseudoJet & jet, int nsub) const {
00686
00687 set<const history_element*> subhist;
00688
00689
00690
00691 get_subhist_set(subhist, jet, -1.0, nsub);
00692
00693 set<const history_element*>::iterator highest = subhist.end();
00694 highest--;
00697 return (*highest)->max_dij_so_far;
00698 }
00699
00700
00701
00702
00709 void ClusterSequence::get_subhist_set(set<const history_element*> & subhist,
00710 const PseudoJet & jet,
00711 double dcut, int maxjet) const {
00712 subhist.clear();
00713 subhist.insert(&(_history[jet.cluster_hist_index()]));
00714
00715
00716 int njet = 1;
00717 while (true) {
00718
00719
00720 set<const history_element*>::iterator highest = subhist.end();
00721 assert (highest != subhist.begin());
00722 highest--;
00723 const history_element* elem = *highest;
00724
00725 if (njet == maxjet) break;
00726
00727 if (elem->parent1 < 0) break;
00728
00729 if (elem->max_dij_so_far <= dcut) break;
00730
00731
00732 subhist.erase(highest);
00733 subhist.insert(&(_history[elem->parent1]));
00734 subhist.insert(&(_history[elem->parent2]));
00735 njet++;
00736 }
00737 }
00738
00739
00740
00741 bool ClusterSequence::object_in_jet(const PseudoJet & object,
00742 const PseudoJet & jet) const {
00743
00744
00745
00746 assert(_potentially_valid(object) && _potentially_valid(jet));
00747
00748 const PseudoJet * this_object = &object;
00749 const PseudoJet * childp;
00750 while(true) {
00751 if (this_object->cluster_hist_index() == jet.cluster_hist_index()) {
00752 return true;
00753 } else if (has_child(*this_object, childp)) {this_object = childp;}
00754 else {return false;}
00755 }
00756 }
00757
00758
00764 bool ClusterSequence::has_parents(const PseudoJet & jet, PseudoJet & parent1,
00765 PseudoJet & parent2) const {
00766
00767 const history_element & hist = _history[jet.cluster_hist_index()];
00768
00769
00770
00771 assert ((hist.parent1 >= 0 && hist.parent2 >= 0) ||
00772 (hist.parent1 < 0 && hist.parent2 < 0));
00773
00774 if (hist.parent1 < 0) {
00775 parent1 = PseudoJet(0.0,0.0,0.0,0.0);
00776 parent2 = parent1;
00777 return false;
00778 } else {
00779 parent1 = _jets[_history[hist.parent1].jetp_index];
00780 parent2 = _jets[_history[hist.parent2].jetp_index];
00781
00782 if (parent1.perp2() < parent2.perp2()) std::swap(parent1,parent2);
00783 return true;
00784 }
00785 }
00786
00787
00790 bool ClusterSequence::has_child(const PseudoJet & jet, PseudoJet & child) const {
00791
00792
00793
00794
00795
00796
00797
00798
00799
00800
00801 const PseudoJet * childp;
00802 bool res = has_child(jet, childp);
00803 if (res) {
00804 child = *childp;
00805 return true;
00806 } else {
00807 child = PseudoJet(0.0,0.0,0.0,0.0);
00808 return false;
00809 }
00810 }
00811
00812 bool ClusterSequence::has_child(const PseudoJet & jet, const PseudoJet * & childp) const {
00813
00814 const history_element & hist = _history[jet.cluster_hist_index()];
00815
00816
00817
00818
00819 if (hist.child >= 0 && _history[hist.child].jetp_index >= 0) {
00820 childp = &(_jets[_history[hist.child].jetp_index]);
00821 return true;
00822 } else {
00823 childp = NULL;
00824 return false;
00825 }
00826 }
00827
00828
00829
00833 bool ClusterSequence::has_partner(const PseudoJet & jet,
00834 PseudoJet & partner) const {
00835
00836 const history_element & hist = _history[jet.cluster_hist_index()];
00837
00838
00839
00840 if (hist.child >= 0 && _history[hist.child].parent2 >= 0) {
00841 const history_element & child_hist = _history[hist.child];
00842 if (child_hist.parent1 == jet.cluster_hist_index()) {
00843
00844
00845 partner = _jets[_history[child_hist.parent2].jetp_index];
00846 } else {
00847
00848 partner = _jets[_history[child_hist.parent1].jetp_index];
00849 }
00850 return true;
00851 } else {
00852 partner = PseudoJet(0.0,0.0,0.0,0.0);
00853 return false;
00854 }
00855 }
00856
00857
00858
00859
00860 vector<PseudoJet> ClusterSequence::constituents (const PseudoJet & jet) const {
00861 vector<PseudoJet> subjets;
00862 add_constituents(jet, subjets);
00863 return subjets;
00864 }
00865
00866
00875 void ClusterSequence::print_jets_for_root(const std::vector<PseudoJet> & jets,
00876 ostream & ostr) const {
00877 for (unsigned i = 0; i < jets.size(); i++) {
00878 ostr << i << " "
00879 << jets[i].px() << " "
00880 << jets[i].py() << " "
00881 << jets[i].pz() << " "
00882 << jets[i].E() << endl;
00883 vector<PseudoJet> cst = constituents(jets[i]);
00884 for (unsigned j = 0; j < cst.size() ; j++) {
00885 ostr << " " << j << " "
00886 << cst[j].rap() << " "
00887 << cst[j].phi() << " "
00888 << cst[j].perp() << endl;
00889 }
00890 ostr << "#END" << endl;
00891 }
00892 }
00893
00894 void ClusterSequence::print_jets_for_root(const std::vector<PseudoJet> & jets,
00895 const std::string & filename,
00896 const std::string & comment ) const {
00897 std::ofstream ostr(filename.c_str());
00898 if (comment != "") ostr << "# " << comment << endl;
00899 print_jets_for_root(jets, ostr);
00900 }
00901
00902
00903
00904
00905
00906
00907
00908
00909
00910
00911
00912
00913
00914
00915
00920 vector<int> ClusterSequence::particle_jet_indices(
00921 const vector<PseudoJet> & jets) const {
00922
00923 vector<int> indices(n_particles());
00924
00925
00926 for (unsigned ipart = 0; ipart < n_particles(); ipart++)
00927 indices[ipart] = -1;
00928
00929
00930
00931 for (unsigned ijet = 0; ijet < jets.size(); ijet++) {
00932
00933 vector<PseudoJet> jet_constituents(constituents(jets[ijet]));
00934
00935 for (unsigned ip = 0; ip < jet_constituents.size(); ip++) {
00936
00937
00938
00939 unsigned iclust = jet_constituents[ip].cluster_hist_index();
00940 unsigned ipart = history()[iclust].jetp_index;
00941 indices[ipart] = ijet;
00942 }
00943 }
00944
00945 return indices;
00946 }
00947
00948
00949
00950
00951 void ClusterSequence::add_constituents (
00952 const PseudoJet & jet, vector<PseudoJet> & subjet_vector) const {
00953
00954 int i = jet.cluster_hist_index();
00955 int parent1 = _history[i].parent1;
00956 int parent2 = _history[i].parent2;
00957
00958 if (parent1 == InexistentParent) {
00959
00960
00961
00962
00963
00964 subjet_vector.push_back(_jets[i]);
00965 return;
00966 }
00967
00968
00969 add_constituents(_jets[_history[parent1].jetp_index], subjet_vector);
00970
00971
00972 if (parent2 != BeamJet) {
00973 add_constituents(_jets[_history[parent2].jetp_index], subjet_vector);
00974 }
00975 }
00976
00977
00978
00979
00980
00981 void ClusterSequence::_add_step_to_history (
00982 const int & step_number, const int & parent1,
00983 const int & parent2, const int & jetp_index,
00984 const double & dij) {
00985
00986 history_element element;
00987 element.parent1 = parent1;
00988 element.parent2 = parent2;
00989 element.jetp_index = jetp_index;
00990 element.child = Invalid;
00991 element.dij = dij;
00992 element.max_dij_so_far = max(dij,_history[_history.size()-1].max_dij_so_far);
00993 _history.push_back(element);
00994
00995 int local_step = _history.size()-1;
00996 assert(local_step == step_number);
00997
00998 assert(parent1 >= 0);
00999 _history[parent1].child = local_step;
01000 if (parent2 >= 0) {_history[parent2].child = local_step;}
01001
01002
01003 if (jetp_index != Invalid) {
01004 assert(jetp_index >= 0);
01005
01006 _jets[jetp_index].set_cluster_hist_index(local_step);
01007 _jets[jetp_index].set_associated_csw(_wrapper_to_this);
01008 }
01009
01010 if (_writeout_combinations) {
01011 cout << local_step << ": "
01012 << parent1 << " with " << parent2
01013 << "; y = "<< dij<<endl;
01014 }
01015
01016 }
01017
01018
01019
01020
01021
01022
01023
01024
01025
01026 vector<int> ClusterSequence::unique_history_order() const {
01027
01028
01029
01030
01031
01032 valarray<int> lowest_constituent(_history.size());
01033 int hist_n = _history.size();
01034 lowest_constituent = hist_n;
01035 for (int i = 0; i < hist_n; i++) {
01036
01037 lowest_constituent[i] = min(lowest_constituent[i],i);
01038
01039 if (_history[i].child > 0) lowest_constituent[_history[i].child]
01040 = min(lowest_constituent[_history[i].child],lowest_constituent[i]);
01041 }
01042
01043
01044 valarray<bool> extracted(_history.size()); extracted = false;
01045 vector<int> unique_tree;
01046 unique_tree.reserve(_history.size());
01047
01048
01049 for (unsigned i = 0; i < n_particles(); i++) {
01050 if (!extracted[i]) {
01051 unique_tree.push_back(i);
01052 extracted[i] = true;
01053 _extract_tree_children(i, extracted, lowest_constituent, unique_tree);
01054 }
01055 }
01056
01057 return unique_tree;
01058 }
01059
01060
01061
01062 void ClusterSequence::_extract_tree_children(
01063 int position,
01064 valarray<bool> & extracted,
01065 const valarray<int> & lowest_constituent,
01066 vector<int> & unique_tree) const {
01067 if (!extracted[position]) {
01068
01069
01070 _extract_tree_parents(position,extracted,lowest_constituent,unique_tree);
01071 }
01072
01073
01074 int child = _history[position].child;
01075 if (child >= 0) _extract_tree_children(child,extracted,lowest_constituent,unique_tree);
01076 }
01077
01078
01079
01080
01081 vector<PseudoJet> ClusterSequence::unclustered_particles() const {
01082 vector<PseudoJet> unclustered;
01083 for (unsigned i = 0; i < n_particles() ; i++) {
01084 if (_history[i].child == Invalid)
01085 unclustered.push_back(_jets[_history[i].jetp_index]);
01086 }
01087 return unclustered;
01088 }
01089
01090
01091
01092
01093
01094 void ClusterSequence::_extract_tree_parents(
01095 int position,
01096 valarray<bool> & extracted,
01097 const valarray<int> & lowest_constituent,
01098 vector<int> & unique_tree) const {
01099
01100 if (!extracted[position]) {
01101 int parent1 = _history[position].parent1;
01102 int parent2 = _history[position].parent2;
01103
01104
01105 if (parent1 >= 0 && parent2 >= 0) {
01106 if (lowest_constituent[parent1] > lowest_constituent[parent2])
01107 std::swap(parent1, parent2);
01108 }
01109
01110 if (parent1 >= 0 && !extracted[parent1])
01111 _extract_tree_parents(parent1,extracted,lowest_constituent,unique_tree);
01112 if (parent2 >= 0 && !extracted[parent2])
01113 _extract_tree_parents(parent2,extracted,lowest_constituent,unique_tree);
01114
01115
01116 unique_tree.push_back(position);
01117 extracted[position] = true;
01118 }
01119 }
01120
01121
01122
01126 void ClusterSequence::_do_ij_recombination_step(
01127 const int & jet_i, const int & jet_j,
01128 const double & dij,
01129 int & newjet_k) {
01130
01131
01132 PseudoJet newjet;
01133 _jet_def.recombiner()->recombine(_jets[jet_i], _jets[jet_j], newjet);
01134 _jets.push_back(newjet);
01135
01136
01137
01138
01139 newjet_k = _jets.size()-1;
01140
01141
01142 int newstep_k = _history.size();
01143
01144 _jets[newjet_k].set_cluster_hist_index(newstep_k);
01145
01146
01147 int hist_i = _jets[jet_i].cluster_hist_index();
01148 int hist_j = _jets[jet_j].cluster_hist_index();
01149
01150 _add_step_to_history(newstep_k, min(hist_i, hist_j), max(hist_i,hist_j),
01151 newjet_k, dij);
01152
01153 }
01154
01155
01156
01159 void ClusterSequence::_do_iB_recombination_step(
01160 const int & jet_i, const double & diB) {
01161
01162 int newstep_k = _history.size();
01163
01164
01165 _add_step_to_history(newstep_k,_jets[jet_i].cluster_hist_index(),BeamJet,
01166 Invalid, diB);
01167
01168 }
01169
01170
01171
01172 LimitedWarning ClusterSequence::_changed_strategy_warning;
01173
01174
01175 FASTJET_END_NAMESPACE
01176