31#include "fastjet/config.h"
32#include "fastjet/Error.hh"
33#include "fastjet/PseudoJet.hh"
34#include "fastjet/ClusterSequence.hh"
35#include "fastjet/ClusterSequenceStructure.hh"
36#include "fastjet/version.hh"
37#include "fastjet/internal/LazyTiling9Alt.hh"
38#include "fastjet/internal/LazyTiling9.hh"
39#include "fastjet/internal/LazyTiling25.hh"
41#include "fastjet/internal/LazyTiling9SeparateGhosts.hh"
52FASTJET_BEGIN_NAMESPACE
148#ifdef FASTJET_HAVE_LIMITED_THREAD_SAFETY
149atomic<ostream *> ClusterSequence::_fastjet_banner_ostr{&cout};
151ostream * ClusterSequence::_fastjet_banner_ostr = &cout;
156ClusterSequence::~ClusterSequence () {
159 if (_structure_shared_ptr){
160 ClusterSequenceStructure* csi =
dynamic_cast<ClusterSequenceStructure*
>(_structure_shared_ptr.get());
166 csi->set_associated_cs(NULL);
178 if (_deletes_self_when_unused) {
179 _structure_shared_ptr.set_count(_structure_shared_ptr.use_count()
180 + _structure_use_count_after_construction);
217void ClusterSequence::signal_imminent_self_deletion()
const {
231 assert(_deletes_self_when_unused);
232 _deletes_self_when_unused =
false;
249void ClusterSequence::_initialise_and_run (
251 const bool & writeout_combinations) {
254 _decant_options(jet_def_in, writeout_combinations);
257 _initialise_and_run_no_decant();
261void ClusterSequence::_initialise_and_run_no_decant () {
265 _fill_initial_history();
268 if (n_particles() == 0)
return;
271 if (_jet_algorithm == plugin_algorithm) {
273 _plugin_activated =
true;
275 _jet_def.plugin()->run_clustering( (*
this) );
276 _plugin_activated =
false;
277 _update_structure_use_count();
279 }
else if (_jet_algorithm == ee_kt_algorithm ||
280 _jet_algorithm == ee_genkt_algorithm) {
281 if (_jet_algorithm == ee_kt_algorithm) {
285 assert(_Rparam > 2.0);
299 _R2 = 2 * ( 3.0 + cos(_Rparam) );
301 _R2 = 2 * ( 1.0 - cos(_Rparam) );
306 if (_strategy == N2PlainEEAccurate) {
307 _simple_N2_cluster_EEAccurateBriefJet();
311 _simple_N2_cluster_EEBriefJet();
314 }
else if (_jet_algorithm == undefined_jet_algorithm) {
315 throw Error(
"A ClusterSequence cannot be created with an uninitialised JetDefinition");
329 if (_strategy == Best) {
330 _strategy = _best_strategy();
335 }
else if (_strategy == BestFJ30) {
336 int N = _jets.size();
342 if (min(1.0,max(0.1,_Rparam)*3.3)*N <= 30) {
344 }
else if (N > 6200/pow(_Rparam,2.0) && _jet_def.jet_algorithm() == cambridge_algorithm) {
347 }
else if ((N > 16000/pow(_Rparam,1.15) && _jet_def.jet_algorithm() != antikt_algorithm)
348 || N > 35000/pow(_Rparam,1.15)) {
351 }
else if (N <= 450) {
362 if (_Rparam >= twopi) {
363 if ( _strategy == NlnN
364 || _strategy == NlnN3pi
365 || _strategy == NlnNCam
366 || _strategy == NlnNCam2pi2R
367 || _strategy == NlnNCam4pi) {
374 if (_jet_def.strategy() != Best && _strategy != _jet_def.strategy()) {
376 oss <<
"Cluster strategy " << strategy_string(_jet_def.strategy())
377 <<
" automatically changed to " << strategy_string()
378 <<
" because the former is not supported for R = " << _Rparam
380 _changed_strategy_warning.warn(oss.str());
390 if (_strategy == N2Plain) {
392 this->_simple_N2_cluster_BriefJet();
393 }
else if (_strategy == N2Tiled) {
394 this->_faster_tiled_N2_cluster();
395 }
else if (_strategy == N2MinHeapTiled) {
396 this->_minheap_faster_tiled_N2_cluster();
397 }
else if (_strategy == N2MHTLazy9Alt) {
400 _plugin_activated =
true;
401 LazyTiling9Alt tiling(*
this);
403 _plugin_activated =
false;
405 }
else if (_strategy == N2MHTLazy25) {
408 _plugin_activated =
true;
409 LazyTiling25 tiling(*
this);
411 _plugin_activated =
false;
413 }
else if (_strategy == N2MHTLazy9) {
416 _plugin_activated =
true;
417 LazyTiling9 tiling(*
this);
419 _plugin_activated =
false;
421 }
else if (_strategy == N2MHTLazy9AntiKtSeparateGhosts) {
425 _plugin_activated =
true;
426 LazyTiling9SeparateGhosts tiling(*
this);
428 _plugin_activated =
false;
430 throw Error(
"N2MHTLazy9AntiKtSeparateGhosts strategy not supported with FJCORE");
433 }
else if (_strategy == NlnN) {
434 this->_delaunay_cluster();
435 }
else if (_strategy == NlnNCam) {
436 this->_CP2DChan_cluster_2piMultD();
437 }
else if (_strategy == NlnN3pi || _strategy == NlnN4pi ) {
438 this->_delaunay_cluster();
439 }
else if (_strategy == N3Dumb ) {
440 this->_really_dumb_cluster();
441 }
else if (_strategy == N2PoorTiled) {
442 this->_tiled_N2_cluster();
443 }
else if (_strategy == NlnNCam4pi) {
444 this->_CP2DChan_cluster();
445 }
else if (_strategy == NlnNCam2pi2R) {
446 this->_CP2DChan_cluster_2pi2R();
449 err <<
"Unrecognised value for strategy: "<<_strategy;
450 throw Error(err.str());
457thread_safety_helpers::FirstTimeTrue ClusterSequence::_first_time;
458LimitedWarning ClusterSequence::_exclusive_warnings;
464 return "FastJet version "+string(fastjet_version);
470void ClusterSequence::print_banner() {
472 if (!_first_time())
return;
475 ostream * ostr = _fastjet_banner_ostr;
478 (*ostr) <<
"#--------------------------------------------------------------------------\n";
479 (*ostr) <<
"# FastJet release " << fastjet_version << endl;
480 (*ostr) <<
"# M. Cacciari, G.P. Salam and G. Soyez \n";
481 (*ostr) <<
"# A software package for jet finding and analysis at colliders \n";
482 (*ostr) <<
"# https://fastjet.fr \n";
484 (*ostr) <<
"# Please cite EPJC72(2012)1896 [arXiv:1111.6097] if you use this package\n";
485 (*ostr) <<
"# for scientific work and optionally PLB641(2006)57 [hep-ph/0512210]. \n";
487 (*ostr) <<
"# FastJet is provided without warranty under the GNU GPL v2 or higher. \n";
488 (*ostr) <<
"# It uses T. Chan's closest pair algorithm, S. Fortune's Voronoi code";
490 (*ostr) <<
",\n# CGAL ";
494 (*ostr) <<
"and 3rd party plugin jet algorithms. See COPYING file for details.\n";
495 (*ostr) <<
"#--------------------------------------------------------------------------\n";
503 const bool & writeout_combinations) {
505 _jet_def = jet_def_in;
506 _writeout_combinations = writeout_combinations;
510 _decant_options_partial();
515void ClusterSequence::_decant_options_partial() {
519 _jet_algorithm = _jet_def.jet_algorithm();
520 _Rparam = _jet_def.R(); _R2 = _Rparam*_Rparam; _invR2 = 1.0/_R2;
521 _strategy = _jet_def.strategy();
524 _plugin_activated =
false;
528 _update_structure_use_count();
534void ClusterSequence::_fill_initial_history () {
539 _jets.reserve(_jets.size()*2);
540 _history.reserve(_jets.size()*2);
544 for (
int i = 0; i < static_cast<int>(_jets.size()) ; i++) {
546 element.
parent1 = InexistentParent;
547 element.
parent2 = InexistentParent;
548 element.
child = Invalid;
553 _history.push_back(element);
556 _jet_def.recombiner()->preprocess(_jets[i]);
559 _jets[i].set_cluster_hist_index(i);
560 _set_structure_shared_ptr(_jets[i]);
563 _Qtot += _jets[i].E();
565 _initial_n = _jets.size();
566 _deletes_self_when_unused =
false;
571string ClusterSequence::strategy_string (
Strategy strategy_in)
const {
573 switch(strategy_in) {
575 strategy =
"NlnN";
break;
577 strategy =
"NlnN3pi";
break;
579 strategy =
"NlnN4pi";
break;
581 strategy =
"N2Plain";
break;
583 strategy =
"N2PlainEEAccurate";
break;
585 strategy =
"N2Tiled";
break;
587 strategy =
"N2MinHeapTiled";
break;
589 strategy =
"N2PoorTiled";
break;
591 strategy =
"N2MHTLazy9";
break;
593 strategy =
"N2MHTLazy9Alt";
break;
595 strategy =
"N2MHTLazy25";
break;
597 strategy =
"N2MHTLazy9AntiKtSeparateGhosts";
break;
599 strategy =
"N3Dumb";
break;
601 strategy =
"NlnNCam4pi";
break;
603 strategy =
"NlnNCam2pi2R";
break;
605 strategy =
"NlnNCam";
break;
607 strategy =
"plugin strategy";
break;
609 strategy =
"Unrecognized";
615double ClusterSequence::jet_scale_for_algorithm(
620 double kt2=jet.
kt2();
621 return kt2 > 1e-300 ? 1.0/kt2 : 1e300;
623 double kt2 = jet.
kt2();
624 double p = jet_def().extra_param();
625 if (p <= 0 && kt2 < 1e-300) kt2 = 1e-300;
628 double kt2 = jet.
kt2();
629 double lim = _jet_def.extra_param();
630 if (kt2 < lim*lim && kt2 != 0.0) {
633 }
else {
throw Error(
"Unrecognised jet algorithm");}
652 int N = _jets.size();
655 double bounded_R = max(_Rparam, 0.1);
659 if (N <= 30 || N <= 39.0/(bounded_R + 0.6)) {
670 const static _Parabola N_Tiled_to_MHT_lowR (-45.4947,54.3528,44.6283);
671 const static _Parabola L_MHT_to_MHTLazy9_lowR (0.677807,-1.05006,10.6994);
672 const static _Parabola L_MHTLazy9_to_MHTLazy25_akt_lowR(0.169967,-0.512589,12.1572);
673 const static _Parabola L_MHTLazy9_to_MHTLazy25_kt_lowR (0.16237,-0.484612,12.3373);
674 const static _Parabola L_MHTLazy9_to_MHTLazy25_cam_lowR = L_MHTLazy9_to_MHTLazy25_kt_lowR;
675 const static _Parabola L_MHTLazy25_to_NlnN_akt_lowR (0.0472051,-0.22043,15.9196);
676 const static _Parabola L_MHTLazy25_to_NlnN_kt_lowR (0.118609,-0.326811,14.8287);
677 const static _Parabola L_MHTLazy25_to_NlnN_cam_lowR (0.10119,-0.295748,14.3924);
679 const static _Line L_Tiled_to_MHTLazy9_medR (-1.31304,7.29621);
680 const static _Parabola L_MHTLazy9_to_MHTLazy25_akt_medR = L_MHTLazy9_to_MHTLazy25_akt_lowR;
681 const static _Parabola L_MHTLazy9_to_MHTLazy25_kt_medR = L_MHTLazy9_to_MHTLazy25_kt_lowR;
682 const static _Parabola L_MHTLazy9_to_MHTLazy25_cam_medR = L_MHTLazy9_to_MHTLazy25_cam_lowR;
683 const static _Parabola L_MHTLazy25_to_NlnN_akt_medR = L_MHTLazy25_to_NlnN_akt_lowR;
684 const static _Parabola L_MHTLazy25_to_NlnN_kt_medR = L_MHTLazy25_to_NlnN_kt_lowR;
685 const static _Parabola L_MHTLazy25_to_NlnN_cam_medR = L_MHTLazy25_to_NlnN_cam_lowR;
687 const static double N_Plain_to_MHTLazy9_largeR = 75;
688 const static double N_MHTLazy9_to_MHTLazy25_akt_largeR = 700;
689 const static double N_MHTLazy9_to_MHTLazy25_kt_largeR = 1000;
690 const static double N_MHTLazy9_to_MHTLazy25_cam_largeR = 1000;
691 const static double N_MHTLazy25_to_NlnN_akt_largeR = 100000;
692 const static double N_MHTLazy25_to_NlnN_kt_largeR = 40000;
693 const static double N_MHTLazy25_to_NlnN_cam_largeR = 15000;
704 double p = jet_def().extra_param();
712 jet_algorithm = _jet_algorithm;
715 if (bounded_R < 0.65) {
717 if (N < N_Tiled_to_MHT_lowR(bounded_R))
return N2Tiled;
718 double logN = log(
double(N));
719 if (logN < L_MHT_to_MHTLazy9_lowR(bounded_R))
return N2MinHeapTiled;
722 if (logN < L_MHTLazy9_to_MHTLazy25_akt_lowR(bounded_R))
return N2MHTLazy9;
723 else if (logN < L_MHTLazy25_to_NlnN_akt_lowR(bounded_R))
return N2MHTLazy25;
726 if (logN < L_MHTLazy9_to_MHTLazy25_kt_lowR(bounded_R))
return N2MHTLazy9;
727 else if (logN < L_MHTLazy25_to_NlnN_kt_lowR(bounded_R))
return N2MHTLazy25;
730 if (logN < L_MHTLazy9_to_MHTLazy25_cam_lowR(bounded_R))
return N2MHTLazy9;
731 else if (logN < L_MHTLazy25_to_NlnN_cam_lowR(bounded_R))
return N2MHTLazy25;
735 }
else if (bounded_R < 0.5*pi) {
737 double logN = log(
double(N));
738 if (logN < L_Tiled_to_MHTLazy9_medR(bounded_R))
return N2Tiled;
741 if (logN < L_MHTLazy9_to_MHTLazy25_akt_medR(bounded_R))
return N2MHTLazy9;
742 else if (logN < L_MHTLazy25_to_NlnN_akt_medR(bounded_R))
return N2MHTLazy25;
745 if (logN < L_MHTLazy9_to_MHTLazy25_kt_medR(bounded_R))
return N2MHTLazy9;
746 else if (logN < L_MHTLazy25_to_NlnN_kt_medR(bounded_R))
return N2MHTLazy25;
749 if (logN < L_MHTLazy9_to_MHTLazy25_cam_medR(bounded_R))
return N2MHTLazy9;
750 else if (logN < L_MHTLazy25_to_NlnN_cam_medR(bounded_R))
return N2MHTLazy25;
756 if (N < N_Plain_to_MHTLazy9_largeR)
return N2Plain;
759 if (N < N_MHTLazy9_to_MHTLazy25_akt_largeR)
return N2MHTLazy9;
760 else if (N < N_MHTLazy25_to_NlnN_akt_largeR)
return N2MHTLazy25;
763 if (N < N_MHTLazy9_to_MHTLazy25_kt_largeR)
return N2MHTLazy9;
764 else if (N < N_MHTLazy25_to_NlnN_kt_largeR)
return N2MHTLazy25;
767 if (N < N_MHTLazy9_to_MHTLazy25_cam_largeR)
return N2MHTLazy9;
768 else if (N < N_MHTLazy25_to_NlnN_cam_largeR)
return N2MHTLazy25;
777 assert(0 &&
"Code should never reach here");
831 _deletes_self_when_unused =
false;
832 transfer_from_sequence(cs);
850 if (will_delete_self_when_unused())
851 throw(
Error(
"cannot use CS::transfer_from_sequence after a call to delete_self_when_unused()"));
854 _jet_def = from_seq._jet_def ;
855 _writeout_combinations = from_seq._writeout_combinations ;
856 _initial_n = from_seq._initial_n ;
857 _Rparam = from_seq._Rparam ;
859 _invR2 = from_seq._invR2 ;
860 _strategy = from_seq._strategy ;
861 _jet_algorithm = from_seq._jet_algorithm ;
862 _plugin_activated = from_seq._plugin_activated ;
868 _jets = (*action_on_jets)(from_seq.
_jets);
870 _jets = from_seq.
_jets;
874 _extras = from_seq._extras;
877 if (_structure_shared_ptr) {
881 if (_deletes_self_when_unused)
throw Error(
"transfer_from_sequence cannot be used for a cluster sequence that deletes self when unused");
892 _update_structure_use_count();
894 for (
unsigned int i=0; i<_jets.size(); i++){
897 _jets[i].set_cluster_hist_index(from_seq.
_jets[i].cluster_hist_index());
900 _set_structure_shared_ptr(_jets[i]);
908void ClusterSequence::plugin_record_ij_recombination(
909 int jet_i,
int jet_j,
double dij,
910 const PseudoJet & newjet,
int & newjet_k) {
912 plugin_record_ij_recombination(jet_i, jet_j, dij, newjet_k);
915 int tmp_index = _jets[newjet_k].cluster_hist_index();
916 _jets[newjet_k] = newjet;
918 _set_structure_shared_ptr(_jets[newjet_k]);
924vector<PseudoJet> ClusterSequence::inclusive_jets (
const double ptmin)
const{
925 double dcut = ptmin*ptmin;
926 int i = _history.size() - 1;
927 vector<PseudoJet> jets_local;
933 if (_history[i].max_dij_so_far < dcut) {
break;}
934 if (_history[i].parent2 == BeamJet && _history[i].dij >= dcut) {
936 int parent1 = _history[i].parent1;
937 jets_local.push_back(_jets[_history[parent1].jetp_index]);}
945 if (_history[i].parent2 != BeamJet) {
break;}
946 int parent1 = _history[i].parent1;
947 const PseudoJet & jet = _jets[_history[parent1].jetp_index];
948 if (jet.
perp2() >= dcut) {jets_local.push_back(jet);}
961 if (_history[i].parent2 == BeamJet) {
962 int parent1 = _history[i].parent1;
963 const PseudoJet & jet = _jets[_history[parent1].jetp_index];
964 if (jet.
perp2() >= dcut) {jets_local.push_back(jet);}
968 }
else {
throw Error(
"cs::inclusive_jets(...): Unrecognized jet algorithm");}
976int ClusterSequence::n_exclusive_jets (
const double dcut)
const {
980 int i = _history.size() - 1;
982 if (_history[i].max_dij_so_far <= dcut) {
break;}
985 int stop_point = i + 1;
988 int njets = 2*_initial_n - stop_point;
995vector<PseudoJet> ClusterSequence::exclusive_jets (
const double dcut)
const {
996 int njets = n_exclusive_jets(dcut);
997 return exclusive_jets(njets);
1004vector<PseudoJet> ClusterSequence::exclusive_jets (
const int njets)
const {
1008 if (njets > _initial_n) {
1010 err <<
"Requested " << njets <<
" exclusive jets, but there were only "
1011 << _initial_n <<
" particles in the event";
1012 throw Error(err.str());
1015 return exclusive_jets_up_to(njets);
1021vector<PseudoJet> ClusterSequence::exclusive_jets_up_to (
const int njets)
const {
1033 (_jet_def.extra_param() <0)) &&
1035 (!_jet_def.plugin()->exclusive_sequence_meaningful()))) {
1036 _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.");
1043 int stop_point = 2*_initial_n - njets;
1045 if (stop_point < _initial_n) stop_point = _initial_n;
1049 if (2*_initial_n !=
static_cast<int>(_history.size())) {
1051 err <<
"2*_initial_n != _history.size() -- this endangers internal assumptions!\n";
1052 throw Error(err.str());
1061 vector<PseudoJet> jets_local;
1062 for (
unsigned int i = stop_point; i < _history.size(); i++) {
1063 int parent1 = _history[i].parent1;
1064 if (parent1 < stop_point) {
1065 jets_local.push_back(_jets[_history[parent1].jetp_index]);
1067 int parent2 = _history[i].parent2;
1068 if (parent2 < stop_point && parent2 > 0) {
1069 jets_local.push_back(_jets[_history[parent2].jetp_index]);
1075 if (
int(jets_local.size()) != min(_initial_n, njets)) {
1077 err <<
"ClusterSequence::exclusive_jets: size of returned vector ("
1078 <<jets_local.size()<<
") does not coincide with requested number of jets ("
1080 throw Error(err.str());
1089double ClusterSequence::exclusive_dmerge (
const int njets)
const {
1091 if (njets >= _initial_n) {
return 0.0;}
1092 return _history[2*_initial_n-njets-1].dij;
1101double ClusterSequence::exclusive_dmerge_max (
const int njets)
const {
1103 if (njets >= _initial_n) {
return 0.0;}
1104 return _history[2*_initial_n-njets-1].max_dij_so_far;
1112std::vector<PseudoJet> ClusterSequence::exclusive_subjets
1113 (
const PseudoJet & jet,
const double dcut)
const {
1115 set<const history_element*> subhist;
1119 get_subhist_set(subhist, jet, dcut, 0);
1122 vector<PseudoJet> subjets;
1123 subjets.reserve(subhist.size());
1124 for (set<const history_element*>::iterator elem = subhist.begin();
1125 elem != subhist.end(); elem++) {
1126 subjets.push_back(_jets[(*elem)->jetp_index]);
1135int ClusterSequence::n_exclusive_subjets(
const PseudoJet & jet,
1136 const double dcut)
const {
1137 set<const history_element*> subhist;
1140 get_subhist_set(subhist, jet, dcut, 0);
1141 return subhist.size();
1148std::vector<PseudoJet> ClusterSequence::exclusive_subjets
1149 (
const PseudoJet & jet,
int nsub)
const {
1150 vector<PseudoJet> subjets = exclusive_subjets_up_to(jet, nsub);
1151 if (
int(subjets.size()) < nsub) {
1153 err <<
"Requested " << nsub <<
" exclusive subjets, but there were only "
1154 << subjets.size() <<
" particles in the jet";
1155 throw Error(err.str());
1165std::vector<PseudoJet> ClusterSequence::exclusive_subjets_up_to
1166 (
const PseudoJet & jet,
int nsub)
const {
1168 set<const history_element*> subhist;
1171 vector<PseudoJet> subjets;
1172 if (nsub < 0)
throw Error(
"Requested a negative number of subjets. This is nonsensical.");
1173 if (nsub == 0)
return subjets;
1177 get_subhist_set(subhist, jet, -1.0, nsub);
1180 subjets.reserve(subhist.size());
1181 for (set<const history_element*>::iterator elem = subhist.begin();
1182 elem != subhist.end(); elem++) {
1183 subjets.push_back(_jets[(*elem)->jetp_index]);
1194double ClusterSequence::exclusive_subdmerge(
const PseudoJet & jet,
int nsub)
const {
1195 set<const history_element*> subhist;
1199 get_subhist_set(subhist, jet, -1.0, nsub);
1201 set<const history_element*>::iterator highest = subhist.end();
1205 return (*highest)->dij;
1215double ClusterSequence::exclusive_subdmerge_max(
const PseudoJet & jet,
int nsub)
const {
1217 set<const history_element*> subhist;
1221 get_subhist_set(subhist, jet, -1.0, nsub);
1223 set<const history_element*>::iterator highest = subhist.end();
1227 return (*highest)->max_dij_so_far;
1239void ClusterSequence::get_subhist_set(set<const history_element*> & subhist,
1241 double dcut,
int maxjet)
const {
1242 assert(contains(jet));
1252 set<const history_element*>::iterator highest = subhist.end();
1253 assert (highest != subhist.begin());
1257 if (njet == maxjet)
break;
1264 subhist.erase(highest);
1265 subhist.insert(&(_history[elem->
parent1]));
1266 subhist.insert(&(_history[elem->
parent2]));
1273bool ClusterSequence::object_in_jet(
const PseudoJet &
object,
1278 assert(contains(
object) && contains(jet));
1280 const PseudoJet * this_object = &object;
1285 }
else if (has_child(*this_object, childp)) {
1286 this_object = childp;
1314 parent1 = _jets[_history[hist.
parent1].jetp_index];
1315 parent2 = _jets[_history[hist.
parent2].jetp_index];
1317 if (parent1.
perp2() < parent2.
perp2()) std::swap(parent1,parent2);
1337 bool res = has_child(jet, childp);
1354 if (hist.
child >= 0 && _history[hist.
child].jetp_index >= 0) {
1355 childp = &(_jets[_history[hist.
child].jetp_index]);
1375 if (hist.
child >= 0 && _history[hist.
child].parent2 >= 0) {
1380 partner = _jets[_history[child_hist.
parent2].jetp_index];
1383 partner = _jets[_history[child_hist.
parent1].jetp_index];
1395vector<PseudoJet> ClusterSequence::constituents (
const PseudoJet & jet)
const {
1396 vector<PseudoJet> subjets;
1397 add_constituents(jet, subjets);
1410void ClusterSequence::print_jets_for_root(
const std::vector<PseudoJet> & jets_in,
1411 ostream & ostr)
const {
1412 for (
unsigned i = 0; i < jets_in.size(); i++) {
1414 << jets_in[i].px() <<
" "
1415 << jets_in[i].py() <<
" "
1416 << jets_in[i].pz() <<
" "
1417 << jets_in[i].E() << endl;
1418 vector<PseudoJet> cst = constituents(jets_in[i]);
1419 for (
unsigned j = 0; j < cst.size() ; j++) {
1420 ostr <<
" " << j <<
" "
1421 << cst[j].rap() <<
" "
1422 << cst[j].phi() <<
" "
1423 << cst[j].perp() << endl;
1425 ostr <<
"#END" << endl;
1429void ClusterSequence::print_jets_for_root(
const std::vector<PseudoJet> & jets_in,
1430 const std::string & filename,
1431 const std::string & comment )
const {
1432 std::ofstream ostr(filename.c_str());
1433 if (comment !=
"") ostr <<
"# " << comment << endl;
1434 print_jets_for_root(jets_in, ostr);
1455vector<int> ClusterSequence::particle_jet_indices(
1456 const vector<PseudoJet> & jets_in)
const {
1458 vector<int> indices(n_particles());
1461 for (
unsigned ipart = 0; ipart < n_particles(); ipart++)
1462 indices[ipart] = -1;
1466 for (
unsigned ijet = 0; ijet < jets_in.size(); ijet++) {
1468 vector<PseudoJet> jet_constituents(constituents(jets_in[ijet]));
1470 for (
unsigned ip = 0; ip < jet_constituents.size(); ip++) {
1474 unsigned iclust = jet_constituents[ip].cluster_hist_index();
1475 unsigned ipart = history()[iclust].jetp_index;
1476 indices[ipart] = ijet;
1486void ClusterSequence::add_constituents (
1487 const PseudoJet & jet, vector<PseudoJet> & subjet_vector)
const {
1490 int parent1 = _history[i].parent1;
1491 int parent2 = _history[i].parent2;
1493 if (parent1 == InexistentParent) {
1499 subjet_vector.push_back(_jets[i]);
1504 add_constituents(_jets[_history[parent1].jetp_index], subjet_vector);
1507 if (parent2 != BeamJet) {
1508 add_constituents(_jets[_history[parent2].jetp_index], subjet_vector);
1516void ClusterSequence::_add_step_to_history (
1519 const int parent2,
const int jetp_index,
1522 history_element element;
1523 element.parent1 = parent1;
1524 element.parent2 = parent2;
1525 element.jetp_index = jetp_index;
1526 element.child = Invalid;
1528 element.max_dij_so_far = max(dij,_history[_history.size()-1].max_dij_so_far);
1529 _history.push_back(element);
1531 int local_step = _history.size()-1;
1543 assert(parent1 >= 0);
1544 if (_history[parent1].child != Invalid){
1545 throw InternalError(
"trying to recombine an object that has previously been recombined");
1547 _history[parent1].child = local_step;
1549 if (_history[parent2].child != Invalid){
1550 throw InternalError(
"trying to recombine an object that has previously been recombined");
1552 _history[parent2].child = local_step;
1556 if (jetp_index != Invalid) {
1557 assert(jetp_index >= 0);
1559 _jets[jetp_index].set_cluster_hist_index(local_step);
1560 _set_structure_shared_ptr(_jets[jetp_index]);
1563 if (_writeout_combinations) {
1564 cout << local_step <<
": "
1565 << parent1 <<
" with " << parent2
1566 <<
"; y = "<< dij<<endl;
1579vector<int> ClusterSequence::unique_history_order()
const {
1585 valarray<int> lowest_constituent(_history.size());
1586 int hist_n = _history.size();
1587 lowest_constituent = hist_n;
1588 for (
int i = 0; i < hist_n; i++) {
1590 lowest_constituent[i] = min(lowest_constituent[i],i);
1592 if (_history[i].child > 0) lowest_constituent[_history[i].child]
1593 = min(lowest_constituent[_history[i].child],lowest_constituent[i]);
1597 valarray<bool> extracted(_history.size()); extracted =
false;
1598 vector<int> unique_tree;
1599 unique_tree.reserve(_history.size());
1602 for (
unsigned i = 0; i < n_particles(); i++) {
1603 if (!extracted[i]) {
1604 unique_tree.push_back(i);
1605 extracted[i] =
true;
1606 _extract_tree_children(i, extracted, lowest_constituent, unique_tree);
1615void ClusterSequence::_extract_tree_children(
1617 valarray<bool> & extracted,
1618 const valarray<int> & lowest_constituent,
1619 vector<int> & unique_tree)
const {
1620 if (!extracted[position]) {
1623 _extract_tree_parents(position,extracted,lowest_constituent,unique_tree);
1627 int child = _history[position].child;
1628 if (child >= 0) _extract_tree_children(child,extracted,lowest_constituent,unique_tree);
1634vector<PseudoJet> ClusterSequence::unclustered_particles()
const {
1635 vector<PseudoJet> unclustered;
1636 for (
unsigned i = 0; i < n_particles() ; i++) {
1637 if (_history[i].child == Invalid)
1638 unclustered.push_back(_jets[_history[i].jetp_index]);
1648vector<PseudoJet> ClusterSequence::childless_pseudojets()
const {
1649 vector<PseudoJet> unclustered;
1650 for (
unsigned i = 0; i < _history.size() ; i++) {
1651 if ((_history[i].child == Invalid) && (_history[i].parent2 != BeamJet))
1652 unclustered.push_back(_jets[_history[i].jetp_index]);
1663bool ClusterSequence::contains(
const PseudoJet & jet)
const {
1674void ClusterSequence::_extract_tree_parents(
1676 valarray<bool> & extracted,
1677 const valarray<int> & lowest_constituent,
1678 vector<int> & unique_tree)
const {
1680 if (!extracted[position]) {
1681 int parent1 = _history[position].parent1;
1682 int parent2 = _history[position].parent2;
1685 if (parent1 >= 0 && parent2 >= 0) {
1686 if (lowest_constituent[parent1] > lowest_constituent[parent2])
1687 std::swap(parent1, parent2);
1690 if (parent1 >= 0 && !extracted[parent1])
1691 _extract_tree_parents(parent1,extracted,lowest_constituent,unique_tree);
1692 if (parent2 >= 0 && !extracted[parent2])
1693 _extract_tree_parents(parent2,extracted,lowest_constituent,unique_tree);
1696 unique_tree.push_back(position);
1697 extracted[position] =
true;
1706void ClusterSequence::_do_ij_recombination_step(
1707 const int jet_i,
const int jet_j,
1717 _jet_def.recombiner()->recombine(_jets[jet_i], _jets[jet_j], newjet);
1718 _jets.push_back(newjet);
1723 newjet_k = _jets.size()-1;
1726 int newstep_k = _history.size();
1728 _jets[newjet_k].set_cluster_hist_index(newstep_k);
1731 int hist_i = _jets[jet_i].cluster_hist_index();
1732 int hist_j = _jets[jet_j].cluster_hist_index();
1734 _add_step_to_history(min(hist_i, hist_j), max(hist_i,hist_j),
1747void ClusterSequence::_do_iB_recombination_step(
1748 const int jet_i,
const double diB) {
1750 _add_step_to_history(_jets[jet_i].cluster_hist_index(),BeamJet,
1767void ClusterSequence::_set_structure_shared_ptr(
PseudoJet & j) {
1771 _update_structure_use_count();
1776void ClusterSequence::_update_structure_use_count() {
1779 _structure_use_count_after_construction = _structure_shared_ptr.use_count();
1786void ClusterSequence::delete_self_when_unused() {
1794 int new_count = _structure_shared_ptr.use_count() - _structure_use_count_after_construction;
1795 if (new_count <= 0) {
1796 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");
1803 _structure_shared_ptr.set_count(new_count);
1805 _deletes_self_when_unused =
true;
1809FASTJET_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 recombined with another jet to form its child.
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.