31#include "fastjet/Error.hh"
32#include "fastjet/PseudoJet.hh"
33#include "fastjet/ClusterSequence.hh"
34#include "fastjet/ClusterSequenceStructure.hh"
35#include "fastjet/version.hh"
36#include "fastjet/internal/LazyTiling9Alt.hh"
37#include "fastjet/internal/LazyTiling9.hh"
38#include "fastjet/internal/LazyTiling25.hh"
40#include "fastjet/internal/LazyTiling9SeparateGhosts.hh"
51FASTJET_BEGIN_NAMESPACE
147#ifdef FASTJET_HAVE_LIMITED_THREAD_SAFETY
148atomic<ostream *> ClusterSequence::_fastjet_banner_ostr{&cout};
150ostream * ClusterSequence::_fastjet_banner_ostr = &cout;
155ClusterSequence::~ClusterSequence () {
158 if (_structure_shared_ptr){
159 ClusterSequenceStructure* csi =
dynamic_cast<ClusterSequenceStructure*
>(_structure_shared_ptr.get());
165 csi->set_associated_cs(NULL);
177 if (_deletes_self_when_unused) {
178 _structure_shared_ptr.set_count(_structure_shared_ptr.use_count()
179 + _structure_use_count_after_construction);
216void ClusterSequence::signal_imminent_self_deletion()
const {
230 assert(_deletes_self_when_unused);
231 _deletes_self_when_unused =
false;
248void ClusterSequence::_initialise_and_run (
250 const bool & writeout_combinations) {
253 _decant_options(jet_def_in, writeout_combinations);
256 _initialise_and_run_no_decant();
260void ClusterSequence::_initialise_and_run_no_decant () {
264 _fill_initial_history();
267 if (n_particles() == 0)
return;
270 if (_jet_algorithm == plugin_algorithm) {
272 _plugin_activated =
true;
274 _jet_def.plugin()->run_clustering( (*
this) );
275 _plugin_activated =
false;
276 _update_structure_use_count();
278 }
else if (_jet_algorithm == ee_kt_algorithm ||
279 _jet_algorithm == ee_genkt_algorithm) {
280 if (_jet_algorithm == ee_kt_algorithm) {
284 assert(_Rparam > 2.0);
298 _R2 = 2 * ( 3.0 + cos(_Rparam) );
300 _R2 = 2 * ( 1.0 - cos(_Rparam) );
305 if (_strategy == N2PlainEEAccurate) {
306 _simple_N2_cluster_EEAccurateBriefJet();
310 _simple_N2_cluster_EEBriefJet();
313 }
else if (_jet_algorithm == undefined_jet_algorithm) {
314 throw Error(
"A ClusterSequence cannot be created with an uninitialised JetDefinition");
328 if (_strategy == Best) {
329 _strategy = _best_strategy();
334 }
else if (_strategy == BestFJ30) {
335 int N = _jets.size();
341 if (min(1.0,max(0.1,_Rparam)*3.3)*N <= 30) {
343 }
else if (N > 6200/pow(_Rparam,2.0) && _jet_def.jet_algorithm() == cambridge_algorithm) {
346 }
else if ((N > 16000/pow(_Rparam,1.15) && _jet_def.jet_algorithm() != antikt_algorithm)
347 || N > 35000/pow(_Rparam,1.15)) {
350 }
else if (N <= 450) {
361 if (_Rparam >= twopi) {
362 if ( _strategy == NlnN
363 || _strategy == NlnN3pi
364 || _strategy == NlnNCam
365 || _strategy == NlnNCam2pi2R
366 || _strategy == NlnNCam4pi) {
373 if (_jet_def.strategy() != Best && _strategy != _jet_def.strategy()) {
375 oss <<
"Cluster strategy " << strategy_string(_jet_def.strategy())
376 <<
" automatically changed to " << strategy_string()
377 <<
" because the former is not supported for R = " << _Rparam
379 _changed_strategy_warning.warn(oss.str());
389 if (_strategy == N2Plain) {
391 this->_simple_N2_cluster_BriefJet();
392 }
else if (_strategy == N2Tiled) {
393 this->_faster_tiled_N2_cluster();
394 }
else if (_strategy == N2MinHeapTiled) {
395 this->_minheap_faster_tiled_N2_cluster();
396 }
else if (_strategy == N2MHTLazy9Alt) {
399 _plugin_activated =
true;
400 LazyTiling9Alt tiling(*
this);
402 _plugin_activated =
false;
404 }
else if (_strategy == N2MHTLazy25) {
407 _plugin_activated =
true;
408 LazyTiling25 tiling(*
this);
410 _plugin_activated =
false;
412 }
else if (_strategy == N2MHTLazy9) {
415 _plugin_activated =
true;
416 LazyTiling9 tiling(*
this);
418 _plugin_activated =
false;
420 }
else if (_strategy == N2MHTLazy9AntiKtSeparateGhosts) {
424 _plugin_activated =
true;
425 LazyTiling9SeparateGhosts tiling(*
this);
427 _plugin_activated =
false;
429 throw Error(
"N2MHTLazy9AntiKtSeparateGhosts strategy not supported with FJCORE");
432 }
else if (_strategy == NlnN) {
433 this->_delaunay_cluster();
434 }
else if (_strategy == NlnNCam) {
435 this->_CP2DChan_cluster_2piMultD();
436 }
else if (_strategy == NlnN3pi || _strategy == NlnN4pi ) {
437 this->_delaunay_cluster();
438 }
else if (_strategy == N3Dumb ) {
439 this->_really_dumb_cluster();
440 }
else if (_strategy == N2PoorTiled) {
441 this->_tiled_N2_cluster();
442 }
else if (_strategy == NlnNCam4pi) {
443 this->_CP2DChan_cluster();
444 }
else if (_strategy == NlnNCam2pi2R) {
445 this->_CP2DChan_cluster_2pi2R();
448 err <<
"Unrecognised value for strategy: "<<_strategy;
449 throw Error(err.str());
456thread_safety_helpers::FirstTimeTrue ClusterSequence::_first_time;
457LimitedWarning ClusterSequence::_exclusive_warnings;
463 return "FastJet version "+string(fastjet_version);
469void ClusterSequence::print_banner() {
471 if (!_first_time())
return;
474 ostream * ostr = _fastjet_banner_ostr;
477 (*ostr) <<
"#--------------------------------------------------------------------------\n";
478 (*ostr) <<
"# FastJet release " << fastjet_version << endl;
479 (*ostr) <<
"# M. Cacciari, G.P. Salam and G. Soyez \n";
480 (*ostr) <<
"# A software package for jet finding and analysis at colliders \n";
481 (*ostr) <<
"# http://fastjet.fr \n";
483 (*ostr) <<
"# Please cite EPJC72(2012)1896 [arXiv:1111.6097] if you use this package\n";
484 (*ostr) <<
"# for scientific work and optionally PLB641(2006)57 [hep-ph/0512210]. \n";
486 (*ostr) <<
"# FastJet is provided without warranty under the GNU GPL v2 or higher. \n";
487 (*ostr) <<
"# It uses T. Chan's closest pair algorithm, S. Fortune's Voronoi code";
489 (*ostr) <<
",\n# CGAL ";
493 (*ostr) <<
"and 3rd party plugin jet algorithms. See COPYING file for details.\n";
494 (*ostr) <<
"#--------------------------------------------------------------------------\n";
502 const bool & writeout_combinations) {
504 _jet_def = jet_def_in;
505 _writeout_combinations = writeout_combinations;
509 _decant_options_partial();
514void ClusterSequence::_decant_options_partial() {
518 _jet_algorithm = _jet_def.jet_algorithm();
519 _Rparam = _jet_def.R(); _R2 = _Rparam*_Rparam; _invR2 = 1.0/_R2;
520 _strategy = _jet_def.strategy();
523 _plugin_activated =
false;
527 _update_structure_use_count();
533void ClusterSequence::_fill_initial_history () {
538 _jets.reserve(_jets.size()*2);
539 _history.reserve(_jets.size()*2);
543 for (
int i = 0; i < static_cast<int>(_jets.size()) ; i++) {
545 element.
parent1 = InexistentParent;
546 element.
parent2 = InexistentParent;
547 element.
child = Invalid;
552 _history.push_back(element);
555 _jet_def.recombiner()->preprocess(_jets[i]);
558 _jets[i].set_cluster_hist_index(i);
559 _set_structure_shared_ptr(_jets[i]);
562 _Qtot += _jets[i].E();
564 _initial_n = _jets.size();
565 _deletes_self_when_unused =
false;
570string ClusterSequence::strategy_string (
Strategy strategy_in)
const {
572 switch(strategy_in) {
574 strategy =
"NlnN";
break;
576 strategy =
"NlnN3pi";
break;
578 strategy =
"NlnN4pi";
break;
580 strategy =
"N2Plain";
break;
582 strategy =
"N2PlainEEAccurate";
break;
584 strategy =
"N2Tiled";
break;
586 strategy =
"N2MinHeapTiled";
break;
588 strategy =
"N2PoorTiled";
break;
590 strategy =
"N2MHTLazy9";
break;
592 strategy =
"N2MHTLazy9Alt";
break;
594 strategy =
"N2MHTLazy25";
break;
596 strategy =
"N2MHTLazy9AntiKtSeparateGhosts";
break;
598 strategy =
"N3Dumb";
break;
600 strategy =
"NlnNCam4pi";
break;
602 strategy =
"NlnNCam2pi2R";
break;
604 strategy =
"NlnNCam";
break;
606 strategy =
"plugin strategy";
break;
608 strategy =
"Unrecognized";
614double ClusterSequence::jet_scale_for_algorithm(
619 double kt2=jet.
kt2();
620 return kt2 > 1e-300 ? 1.0/kt2 : 1e300;
622 double kt2 = jet.
kt2();
623 double p = jet_def().extra_param();
624 if (p <= 0 && kt2 < 1e-300) kt2 = 1e-300;
627 double kt2 = jet.
kt2();
628 double lim = _jet_def.extra_param();
629 if (kt2 < lim*lim && kt2 != 0.0) {
632 }
else {
throw Error(
"Unrecognised jet algorithm");}
651 int N = _jets.size();
654 double bounded_R = max(_Rparam, 0.1);
658 if (N <= 30 || N <= 39.0/(bounded_R + 0.6)) {
669 const static _Parabola N_Tiled_to_MHT_lowR (-45.4947,54.3528,44.6283);
670 const static _Parabola L_MHT_to_MHTLazy9_lowR (0.677807,-1.05006,10.6994);
671 const static _Parabola L_MHTLazy9_to_MHTLazy25_akt_lowR(0.169967,-0.512589,12.1572);
672 const static _Parabola L_MHTLazy9_to_MHTLazy25_kt_lowR (0.16237,-0.484612,12.3373);
673 const static _Parabola L_MHTLazy9_to_MHTLazy25_cam_lowR = L_MHTLazy9_to_MHTLazy25_kt_lowR;
674 const static _Parabola L_MHTLazy25_to_NlnN_akt_lowR (0.0472051,-0.22043,15.9196);
675 const static _Parabola L_MHTLazy25_to_NlnN_kt_lowR (0.118609,-0.326811,14.8287);
676 const static _Parabola L_MHTLazy25_to_NlnN_cam_lowR (0.10119,-0.295748,14.3924);
678 const static _Line L_Tiled_to_MHTLazy9_medR (-1.31304,7.29621);
679 const static _Parabola L_MHTLazy9_to_MHTLazy25_akt_medR = L_MHTLazy9_to_MHTLazy25_akt_lowR;
680 const static _Parabola L_MHTLazy9_to_MHTLazy25_kt_medR = L_MHTLazy9_to_MHTLazy25_kt_lowR;
681 const static _Parabola L_MHTLazy9_to_MHTLazy25_cam_medR = L_MHTLazy9_to_MHTLazy25_cam_lowR;
682 const static _Parabola L_MHTLazy25_to_NlnN_akt_medR = L_MHTLazy25_to_NlnN_akt_lowR;
683 const static _Parabola L_MHTLazy25_to_NlnN_kt_medR = L_MHTLazy25_to_NlnN_kt_lowR;
684 const static _Parabola L_MHTLazy25_to_NlnN_cam_medR = L_MHTLazy25_to_NlnN_cam_lowR;
686 const static double N_Plain_to_MHTLazy9_largeR = 75;
687 const static double N_MHTLazy9_to_MHTLazy25_akt_largeR = 700;
688 const static double N_MHTLazy9_to_MHTLazy25_kt_largeR = 1000;
689 const static double N_MHTLazy9_to_MHTLazy25_cam_largeR = 1000;
690 const static double N_MHTLazy25_to_NlnN_akt_largeR = 100000;
691 const static double N_MHTLazy25_to_NlnN_kt_largeR = 40000;
692 const static double N_MHTLazy25_to_NlnN_cam_largeR = 15000;
703 double p = jet_def().extra_param();
711 jet_algorithm = _jet_algorithm;
714 if (bounded_R < 0.65) {
716 if (N < N_Tiled_to_MHT_lowR(bounded_R))
return N2Tiled;
717 double logN = log(
double(N));
718 if (logN < L_MHT_to_MHTLazy9_lowR(bounded_R))
return N2MinHeapTiled;
721 if (logN < L_MHTLazy9_to_MHTLazy25_akt_lowR(bounded_R))
return N2MHTLazy9;
722 else if (logN < L_MHTLazy25_to_NlnN_akt_lowR(bounded_R))
return N2MHTLazy25;
725 if (logN < L_MHTLazy9_to_MHTLazy25_kt_lowR(bounded_R))
return N2MHTLazy9;
726 else if (logN < L_MHTLazy25_to_NlnN_kt_lowR(bounded_R))
return N2MHTLazy25;
729 if (logN < L_MHTLazy9_to_MHTLazy25_cam_lowR(bounded_R))
return N2MHTLazy9;
730 else if (logN < L_MHTLazy25_to_NlnN_cam_lowR(bounded_R))
return N2MHTLazy25;
734 }
else if (bounded_R < 0.5*pi) {
736 double logN = log(
double(N));
737 if (logN < L_Tiled_to_MHTLazy9_medR(bounded_R))
return N2Tiled;
740 if (logN < L_MHTLazy9_to_MHTLazy25_akt_medR(bounded_R))
return N2MHTLazy9;
741 else if (logN < L_MHTLazy25_to_NlnN_akt_medR(bounded_R))
return N2MHTLazy25;
744 if (logN < L_MHTLazy9_to_MHTLazy25_kt_medR(bounded_R))
return N2MHTLazy9;
745 else if (logN < L_MHTLazy25_to_NlnN_kt_medR(bounded_R))
return N2MHTLazy25;
748 if (logN < L_MHTLazy9_to_MHTLazy25_cam_medR(bounded_R))
return N2MHTLazy9;
749 else if (logN < L_MHTLazy25_to_NlnN_cam_medR(bounded_R))
return N2MHTLazy25;
755 if (N < N_Plain_to_MHTLazy9_largeR)
return N2Plain;
758 if (N < N_MHTLazy9_to_MHTLazy25_akt_largeR)
return N2MHTLazy9;
759 else if (N < N_MHTLazy25_to_NlnN_akt_largeR)
return N2MHTLazy25;
762 if (N < N_MHTLazy9_to_MHTLazy25_kt_largeR)
return N2MHTLazy9;
763 else if (N < N_MHTLazy25_to_NlnN_kt_largeR)
return N2MHTLazy25;
766 if (N < N_MHTLazy9_to_MHTLazy25_cam_largeR)
return N2MHTLazy9;
767 else if (N < N_MHTLazy25_to_NlnN_cam_largeR)
return N2MHTLazy25;
776 assert(0 &&
"Code should never reach here");
830 _deletes_self_when_unused =
false;
831 transfer_from_sequence(cs);
849 if (will_delete_self_when_unused())
850 throw(
Error(
"cannot use CS::transfer_from_sequence after a call to delete_self_when_unused()"));
853 _jet_def = from_seq._jet_def ;
854 _writeout_combinations = from_seq._writeout_combinations ;
855 _initial_n = from_seq._initial_n ;
856 _Rparam = from_seq._Rparam ;
858 _invR2 = from_seq._invR2 ;
859 _strategy = from_seq._strategy ;
860 _jet_algorithm = from_seq._jet_algorithm ;
861 _plugin_activated = from_seq._plugin_activated ;
867 _jets = (*action_on_jets)(from_seq.
_jets);
869 _jets = from_seq.
_jets;
873 _extras = from_seq._extras;
876 if (_structure_shared_ptr) {
880 if (_deletes_self_when_unused)
throw Error(
"transfer_from_sequence cannot be used for a cluster sequence that deletes self when unused");
891 _update_structure_use_count();
893 for (
unsigned int i=0; i<_jets.size(); i++){
896 _jets[i].set_cluster_hist_index(from_seq.
_jets[i].cluster_hist_index());
899 _set_structure_shared_ptr(_jets[i]);
907void ClusterSequence::plugin_record_ij_recombination(
908 int jet_i,
int jet_j,
double dij,
909 const PseudoJet & newjet,
int & newjet_k) {
911 plugin_record_ij_recombination(jet_i, jet_j, dij, newjet_k);
914 int tmp_index = _jets[newjet_k].cluster_hist_index();
915 _jets[newjet_k] = newjet;
917 _set_structure_shared_ptr(_jets[newjet_k]);
923vector<PseudoJet> ClusterSequence::inclusive_jets (
const double ptmin)
const{
924 double dcut = ptmin*ptmin;
925 int i = _history.size() - 1;
926 vector<PseudoJet> jets_local;
932 if (_history[i].max_dij_so_far < dcut) {
break;}
933 if (_history[i].parent2 == BeamJet && _history[i].dij >= dcut) {
935 int parent1 = _history[i].parent1;
936 jets_local.push_back(_jets[_history[parent1].jetp_index]);}
944 if (_history[i].parent2 != BeamJet) {
break;}
945 int parent1 = _history[i].parent1;
946 const PseudoJet & jet = _jets[_history[parent1].jetp_index];
947 if (jet.
perp2() >= dcut) {jets_local.push_back(jet);}
960 if (_history[i].parent2 == BeamJet) {
961 int parent1 = _history[i].parent1;
962 const PseudoJet & jet = _jets[_history[parent1].jetp_index];
963 if (jet.
perp2() >= dcut) {jets_local.push_back(jet);}
967 }
else {
throw Error(
"cs::inclusive_jets(...): Unrecognized jet algorithm");}
975int ClusterSequence::n_exclusive_jets (
const double dcut)
const {
979 int i = _history.size() - 1;
981 if (_history[i].max_dij_so_far <= dcut) {
break;}
984 int stop_point = i + 1;
987 int njets = 2*_initial_n - stop_point;
994vector<PseudoJet> ClusterSequence::exclusive_jets (
const double dcut)
const {
995 int njets = n_exclusive_jets(dcut);
996 return exclusive_jets(njets);
1003vector<PseudoJet> ClusterSequence::exclusive_jets (
const int njets)
const {
1007 if (njets > _initial_n) {
1009 err <<
"Requested " << njets <<
" exclusive jets, but there were only "
1010 << _initial_n <<
" particles in the event";
1011 throw Error(err.str());
1014 return exclusive_jets_up_to(njets);
1020vector<PseudoJet> ClusterSequence::exclusive_jets_up_to (
const int njets)
const {
1032 (_jet_def.extra_param() <0)) &&
1034 (!_jet_def.plugin()->exclusive_sequence_meaningful()))) {
1035 _exclusive_warnings.warn(
"dcut and exclusive jets for jet-finders other than kt, C/A or genkt with p>=0 should be interpreted with care.");
1042 int stop_point = 2*_initial_n - njets;
1044 if (stop_point < _initial_n) stop_point = _initial_n;
1048 if (2*_initial_n !=
static_cast<int>(_history.size())) {
1050 err <<
"2*_initial_n != _history.size() -- this endangers internal assumptions!\n";
1051 throw Error(err.str());
1060 vector<PseudoJet> jets_local;
1061 for (
unsigned int i = stop_point; i < _history.size(); i++) {
1062 int parent1 = _history[i].parent1;
1063 if (parent1 < stop_point) {
1064 jets_local.push_back(_jets[_history[parent1].jetp_index]);
1066 int parent2 = _history[i].parent2;
1067 if (parent2 < stop_point && parent2 > 0) {
1068 jets_local.push_back(_jets[_history[parent2].jetp_index]);
1074 if (
int(jets_local.size()) != min(_initial_n, njets)) {
1076 err <<
"ClusterSequence::exclusive_jets: size of returned vector ("
1077 <<jets_local.size()<<
") does not coincide with requested number of jets ("
1079 throw Error(err.str());
1088double ClusterSequence::exclusive_dmerge (
const int njets)
const {
1090 if (njets >= _initial_n) {
return 0.0;}
1091 return _history[2*_initial_n-njets-1].dij;
1100double ClusterSequence::exclusive_dmerge_max (
const int njets)
const {
1102 if (njets >= _initial_n) {
return 0.0;}
1103 return _history[2*_initial_n-njets-1].max_dij_so_far;
1111std::vector<PseudoJet> ClusterSequence::exclusive_subjets
1112 (
const PseudoJet & jet,
const double dcut)
const {
1114 set<const history_element*> subhist;
1118 get_subhist_set(subhist, jet, dcut, 0);
1121 vector<PseudoJet> subjets;
1122 subjets.reserve(subhist.size());
1123 for (set<const history_element*>::iterator elem = subhist.begin();
1124 elem != subhist.end(); elem++) {
1125 subjets.push_back(_jets[(*elem)->jetp_index]);
1134int ClusterSequence::n_exclusive_subjets(
const PseudoJet & jet,
1135 const double dcut)
const {
1136 set<const history_element*> subhist;
1139 get_subhist_set(subhist, jet, dcut, 0);
1140 return subhist.size();
1147std::vector<PseudoJet> ClusterSequence::exclusive_subjets
1148 (
const PseudoJet & jet,
int nsub)
const {
1149 vector<PseudoJet> subjets = exclusive_subjets_up_to(jet, nsub);
1150 if (
int(subjets.size()) < nsub) {
1152 err <<
"Requested " << nsub <<
" exclusive subjets, but there were only "
1153 << subjets.size() <<
" particles in the jet";
1154 throw Error(err.str());
1164std::vector<PseudoJet> ClusterSequence::exclusive_subjets_up_to
1165 (
const PseudoJet & jet,
int nsub)
const {
1167 set<const history_element*> subhist;
1170 vector<PseudoJet> subjets;
1171 if (nsub < 0)
throw Error(
"Requested a negative number of subjets. This is nonsensical.");
1172 if (nsub == 0)
return subjets;
1176 get_subhist_set(subhist, jet, -1.0, nsub);
1179 subjets.reserve(subhist.size());
1180 for (set<const history_element*>::iterator elem = subhist.begin();
1181 elem != subhist.end(); elem++) {
1182 subjets.push_back(_jets[(*elem)->jetp_index]);
1193double ClusterSequence::exclusive_subdmerge(
const PseudoJet & jet,
int nsub)
const {
1194 set<const history_element*> subhist;
1198 get_subhist_set(subhist, jet, -1.0, nsub);
1200 set<const history_element*>::iterator highest = subhist.end();
1204 return (*highest)->dij;
1214double ClusterSequence::exclusive_subdmerge_max(
const PseudoJet & jet,
int nsub)
const {
1216 set<const history_element*> subhist;
1220 get_subhist_set(subhist, jet, -1.0, nsub);
1222 set<const history_element*>::iterator highest = subhist.end();
1226 return (*highest)->max_dij_so_far;
1238void ClusterSequence::get_subhist_set(set<const history_element*> & subhist,
1240 double dcut,
int maxjet)
const {
1241 assert(contains(jet));
1251 set<const history_element*>::iterator highest = subhist.end();
1252 assert (highest != subhist.begin());
1256 if (njet == maxjet)
break;
1263 subhist.erase(highest);
1264 subhist.insert(&(_history[elem->
parent1]));
1265 subhist.insert(&(_history[elem->
parent2]));
1272bool ClusterSequence::object_in_jet(
const PseudoJet &
object,
1277 assert(contains(
object) && contains(jet));
1279 const PseudoJet * this_object = &object;
1284 }
else if (has_child(*this_object, childp)) {
1285 this_object = childp;
1313 parent1 = _jets[_history[hist.
parent1].jetp_index];
1314 parent2 = _jets[_history[hist.
parent2].jetp_index];
1316 if (parent1.
perp2() < parent2.
perp2()) std::swap(parent1,parent2);
1336 bool res = has_child(jet, childp);
1353 if (hist.
child >= 0 && _history[hist.
child].jetp_index >= 0) {
1354 childp = &(_jets[_history[hist.
child].jetp_index]);
1374 if (hist.
child >= 0 && _history[hist.
child].parent2 >= 0) {
1379 partner = _jets[_history[child_hist.
parent2].jetp_index];
1382 partner = _jets[_history[child_hist.
parent1].jetp_index];
1394vector<PseudoJet> ClusterSequence::constituents (
const PseudoJet & jet)
const {
1395 vector<PseudoJet> subjets;
1396 add_constituents(jet, subjets);
1409void ClusterSequence::print_jets_for_root(
const std::vector<PseudoJet> & jets_in,
1410 ostream & ostr)
const {
1411 for (
unsigned i = 0; i < jets_in.size(); i++) {
1413 << jets_in[i].px() <<
" "
1414 << jets_in[i].py() <<
" "
1415 << jets_in[i].pz() <<
" "
1416 << jets_in[i].E() << endl;
1417 vector<PseudoJet> cst = constituents(jets_in[i]);
1418 for (
unsigned j = 0; j < cst.size() ; j++) {
1419 ostr <<
" " << j <<
" "
1420 << cst[j].rap() <<
" "
1421 << cst[j].phi() <<
" "
1422 << cst[j].perp() << endl;
1424 ostr <<
"#END" << endl;
1428void ClusterSequence::print_jets_for_root(
const std::vector<PseudoJet> & jets_in,
1429 const std::string & filename,
1430 const std::string & comment )
const {
1431 std::ofstream ostr(filename.c_str());
1432 if (comment !=
"") ostr <<
"# " << comment << endl;
1433 print_jets_for_root(jets_in, ostr);
1454vector<int> ClusterSequence::particle_jet_indices(
1455 const vector<PseudoJet> & jets_in)
const {
1457 vector<int> indices(n_particles());
1460 for (
unsigned ipart = 0; ipart < n_particles(); ipart++)
1461 indices[ipart] = -1;
1465 for (
unsigned ijet = 0; ijet < jets_in.size(); ijet++) {
1467 vector<PseudoJet> jet_constituents(constituents(jets_in[ijet]));
1469 for (
unsigned ip = 0; ip < jet_constituents.size(); ip++) {
1473 unsigned iclust = jet_constituents[ip].cluster_hist_index();
1474 unsigned ipart = history()[iclust].jetp_index;
1475 indices[ipart] = ijet;
1485void ClusterSequence::add_constituents (
1486 const PseudoJet & jet, vector<PseudoJet> & subjet_vector)
const {
1489 int parent1 = _history[i].parent1;
1490 int parent2 = _history[i].parent2;
1492 if (parent1 == InexistentParent) {
1498 subjet_vector.push_back(_jets[i]);
1503 add_constituents(_jets[_history[parent1].jetp_index], subjet_vector);
1506 if (parent2 != BeamJet) {
1507 add_constituents(_jets[_history[parent2].jetp_index], subjet_vector);
1515void ClusterSequence::_add_step_to_history (
1518 const int parent2,
const int jetp_index,
1521 history_element element;
1522 element.parent1 = parent1;
1523 element.parent2 = parent2;
1524 element.jetp_index = jetp_index;
1525 element.child = Invalid;
1527 element.max_dij_so_far = max(dij,_history[_history.size()-1].max_dij_so_far);
1528 _history.push_back(element);
1530 int local_step = _history.size()-1;
1542 assert(parent1 >= 0);
1543 if (_history[parent1].child != Invalid){
1544 throw InternalError(
"trying to recombine an object that has previously been recombined");
1546 _history[parent1].child = local_step;
1548 if (_history[parent2].child != Invalid){
1549 throw InternalError(
"trying to recombine an object that has previously been recombined");
1551 _history[parent2].child = local_step;
1555 if (jetp_index != Invalid) {
1556 assert(jetp_index >= 0);
1558 _jets[jetp_index].set_cluster_hist_index(local_step);
1559 _set_structure_shared_ptr(_jets[jetp_index]);
1562 if (_writeout_combinations) {
1563 cout << local_step <<
": "
1564 << parent1 <<
" with " << parent2
1565 <<
"; y = "<< dij<<endl;
1578vector<int> ClusterSequence::unique_history_order()
const {
1584 valarray<int> lowest_constituent(_history.size());
1585 int hist_n = _history.size();
1586 lowest_constituent = hist_n;
1587 for (
int i = 0; i < hist_n; i++) {
1589 lowest_constituent[i] = min(lowest_constituent[i],i);
1591 if (_history[i].child > 0) lowest_constituent[_history[i].child]
1592 = min(lowest_constituent[_history[i].child],lowest_constituent[i]);
1596 valarray<bool> extracted(_history.size()); extracted =
false;
1597 vector<int> unique_tree;
1598 unique_tree.reserve(_history.size());
1601 for (
unsigned i = 0; i < n_particles(); i++) {
1602 if (!extracted[i]) {
1603 unique_tree.push_back(i);
1604 extracted[i] =
true;
1605 _extract_tree_children(i, extracted, lowest_constituent, unique_tree);
1614void ClusterSequence::_extract_tree_children(
1616 valarray<bool> & extracted,
1617 const valarray<int> & lowest_constituent,
1618 vector<int> & unique_tree)
const {
1619 if (!extracted[position]) {
1622 _extract_tree_parents(position,extracted,lowest_constituent,unique_tree);
1626 int child = _history[position].child;
1627 if (child >= 0) _extract_tree_children(child,extracted,lowest_constituent,unique_tree);
1633vector<PseudoJet> ClusterSequence::unclustered_particles()
const {
1634 vector<PseudoJet> unclustered;
1635 for (
unsigned i = 0; i < n_particles() ; i++) {
1636 if (_history[i].child == Invalid)
1637 unclustered.push_back(_jets[_history[i].jetp_index]);
1647vector<PseudoJet> ClusterSequence::childless_pseudojets()
const {
1648 vector<PseudoJet> unclustered;
1649 for (
unsigned i = 0; i < _history.size() ; i++) {
1650 if ((_history[i].child == Invalid) && (_history[i].parent2 != BeamJet))
1651 unclustered.push_back(_jets[_history[i].jetp_index]);
1662bool ClusterSequence::contains(
const PseudoJet & jet)
const {
1673void ClusterSequence::_extract_tree_parents(
1675 valarray<bool> & extracted,
1676 const valarray<int> & lowest_constituent,
1677 vector<int> & unique_tree)
const {
1679 if (!extracted[position]) {
1680 int parent1 = _history[position].parent1;
1681 int parent2 = _history[position].parent2;
1684 if (parent1 >= 0 && parent2 >= 0) {
1685 if (lowest_constituent[parent1] > lowest_constituent[parent2])
1686 std::swap(parent1, parent2);
1689 if (parent1 >= 0 && !extracted[parent1])
1690 _extract_tree_parents(parent1,extracted,lowest_constituent,unique_tree);
1691 if (parent2 >= 0 && !extracted[parent2])
1692 _extract_tree_parents(parent2,extracted,lowest_constituent,unique_tree);
1695 unique_tree.push_back(position);
1696 extracted[position] =
true;
1705void ClusterSequence::_do_ij_recombination_step(
1706 const int jet_i,
const int jet_j,
1716 _jet_def.recombiner()->recombine(_jets[jet_i], _jets[jet_j], newjet);
1717 _jets.push_back(newjet);
1722 newjet_k = _jets.size()-1;
1725 int newstep_k = _history.size();
1727 _jets[newjet_k].set_cluster_hist_index(newstep_k);
1730 int hist_i = _jets[jet_i].cluster_hist_index();
1731 int hist_j = _jets[jet_j].cluster_hist_index();
1733 _add_step_to_history(min(hist_i, hist_j), max(hist_i,hist_j),
1746void ClusterSequence::_do_iB_recombination_step(
1747 const int jet_i,
const double diB) {
1749 _add_step_to_history(_jets[jet_i].cluster_hist_index(),BeamJet,
1766void ClusterSequence::_set_structure_shared_ptr(
PseudoJet & j) {
1770 _update_structure_use_count();
1775void ClusterSequence::_update_structure_use_count() {
1778 _structure_use_count_after_construction = _structure_shared_ptr.use_count();
1785void ClusterSequence::delete_self_when_unused() {
1793 int new_count = _structure_shared_ptr.use_count() - _structure_use_count_after_construction;
1794 if (new_count <= 0) {
1795 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");
1802 _structure_shared_ptr.set_count(new_count);
1804 _deletes_self_when_unused =
true;
1808FASTJET_END_NAMESPACE
Contains any information related to the clustering that should be directly accessible to PseudoJet.
virtual void set_associated_cs(const ClusterSequence *new_cs)
set the associated csw
std::vector< history_element > _history
this vector will contain the branching history; for each stage, _history[i].jetp_index indicates wher...
std::vector< PseudoJet > _jets
This contains the physical PseudoJets; for each PseudoJet one can find the corresponding position in ...
base class corresponding to errors that can be thrown by FastJet
base class providing interface for a generic function of a PseudoJet
class corresponding to critical internal errors
class that is intended to hold a full definition of the jet clusterer
class to provide facilities for giving warnings up to some maximum number of times and to provide glo...
Class to contain pseudojets, including minimal information of use to jet-clustering routines.
double perp2() const
returns the squared transverse momentum
void set_cluster_hist_index(const int index)
set the cluster_hist_index, intended to be used by clustering routines.
double kt2() const
returns the squared transverse momentum
bool has_valid_cluster_sequence() const
returns true if this PseudoJet has an associated and still valid(ated) ClusterSequence.
int cluster_hist_index() const
return the cluster_hist_index, intended to be used by clustering routines.
void set_structure_shared_ptr(const SharedPtr< PseudoJetStructureBase > &structure_in)
set the associated structure
const ClusterSequence * associated_cluster_sequence() const
get a (const) pointer to the parent ClusterSequence (NULL if inexistent)
Strategy
the various options for the algorithmic strategy to adopt in clustering events with kt and cambridge ...
@ N2Plain
fastest below 50
@ N3Dumb
worse even than the usual N^3 algorithms
@ N2MHTLazy9AntiKtSeparateGhosts
Like N2MHTLazy9 in a number of respects, but does not calculate ghost-ghost distances and so does not...
@ N2Tiled
fastest from about 50..500
@ NlnN3pi
legacy N ln N using 3pi coverage of cylinder.
@ N2MHTLazy9Alt
Like to N2MHTLazy9 but uses slightly different optimizations, e.g.
@ N2PlainEEAccurate
a variant of N2Plain strategy for native algorithms (eekt, ee_genkt) that uses a more accurate calcu...
@ plugin_strategy
the plugin has been used...
@ N2MHTLazy9
only looks into a neighbouring tile for a particle's nearest neighbour (NN) if that particle's in-til...
@ NlnN
best of the NlnN variants – best overall for N>10^4.
@ NlnN4pi
legacy N ln N using 4pi coverage of cylinder
@ NlnNCam
Chan's closest pair method (in a variant with 2pi+minimal extra variant), for use exclusively with th...
@ NlnNCam2pi2R
Chan's closest pair method (in a variant with 2pi+2R coverage), for use exclusively with the Cambridg...
@ N2MHTLazy25
Similar to N2MHTLazy9, but uses tiles of size R/2 and a 5x5 tile grid around the particle.
@ NlnNCam4pi
Chan's closest pair method (in a variant with 4pi coverage), for use exclusively with the Cambridge a...
@ N2MinHeapTiled
faster that N2Tiled above about 500 particles; differs from it by retainig the di(closest j) distance...
JetAlgorithm
the various families of jet-clustering algorithm
@ ee_genkt_algorithm
the e+e- genkt algorithm (R > 2 and p=1 gives ee_kt)
@ genkt_algorithm
like the k_t but with distance measures dij = min(kti^{2p},ktj^{2p}) Delta R_{ij}^2 / R^2 diB = 1/kti...
@ plugin_algorithm
any plugin algorithm supplied by the user
@ ee_kt_algorithm
the e+e- kt algorithm
@ cambridge_algorithm
the longitudinally invariant variant of the cambridge algorithm (aka Aachen algoithm).
@ cambridge_for_passive_algorithm
a version of cambridge with a special distance measure for particles whose pt is < extra_param(); thi...
@ antikt_algorithm
like the k_t but with distance measures dij = min(1/kti^2,1/ktj^2) Delta R_{ij}^2 / R^2 diB = 1/kti^2
@ kt_algorithm
the longitudinally invariant kt algorithm
string fastjet_version_string()
return a string containing information about the release
a single element in the clustering history
int parent1
index in _history where first parent of this jet was created (InexistentParent if this jet is an orig...
int jetp_index
index in the _jets vector where we will find the PseudoJet object corresponding to this jet (i....
int parent2
index in _history where second parent of this jet was created (InexistentParent if this jet is an ori...
int child
index in _history where the current jet is
double max_dij_so_far
the largest recombination distance seen so far in the clustering history.
double dij
the distance corresponding to the recombination at this stage of the clustering.