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