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"
51 FASTJET_BEGIN_NAMESPACE
147 std::ostream * ClusterSequence::_fastjet_banner_ostr = &cout;
151 ClusterSequence::~ClusterSequence () {
154 if (_structure_shared_ptr()){
155 ClusterSequenceStructure* csi =
dynamic_cast<ClusterSequenceStructure*
>(_structure_shared_ptr());
161 csi->set_associated_cs(NULL);
169 if (_deletes_self_when_unused) {
170 _structure_shared_ptr.set_count(_structure_shared_ptr.use_count()
171 + _structure_use_count_after_construction);
177 void ClusterSequence::signal_imminent_self_deletion()
const {
191 assert(_deletes_self_when_unused);
192 _deletes_self_when_unused =
false;
207 void ClusterSequence::_initialise_and_run (
209 const bool & writeout_combinations) {
212 _decant_options(jet_def_in, writeout_combinations);
215 _initialise_and_run_no_decant();
219 void ClusterSequence::_initialise_and_run_no_decant () {
223 _fill_initial_history();
226 if (n_particles() == 0)
return;
229 if (_jet_algorithm == plugin_algorithm) {
231 _plugin_activated =
true;
233 _jet_def.plugin()->run_clustering( (*
this) );
234 _plugin_activated =
false;
235 _update_structure_use_count();
237 }
else if (_jet_algorithm == ee_kt_algorithm ||
238 _jet_algorithm == ee_genkt_algorithm) {
241 if (_jet_algorithm == ee_kt_algorithm) {
245 assert(_Rparam > 2.0);
259 _R2 = 2 * ( 3.0 + cos(_Rparam) );
261 _R2 = 2 * ( 1.0 - cos(_Rparam) );
265 _simple_N2_cluster_EEBriefJet();
267 }
else if (_jet_algorithm == undefined_jet_algorithm) {
268 throw Error(
"A ClusterSequence cannot be created with an uninitialised JetDefinition");
282 if (_strategy == Best) {
283 _strategy = _best_strategy();
288 }
else if (_strategy == BestFJ30) {
289 int N = _jets.size();
295 if (min(1.0,max(0.1,_Rparam)*3.3)*N <= 30) {
300 }
else if ((N > 16000/pow(_Rparam,1.15) && _jet_def.jet_algorithm() !=
antikt_algorithm)
301 || N > 35000/pow(_Rparam,1.15)) {
304 }
else if (N <= 450) {
315 if (_Rparam >= twopi) {
316 if ( _strategy == NlnN
317 || _strategy == NlnN3pi
318 || _strategy == NlnNCam
319 || _strategy == NlnNCam2pi2R
320 || _strategy == NlnNCam4pi) {
327 if (_jet_def.strategy() !=
Best && _strategy != _jet_def.strategy()) {
329 oss <<
"Cluster strategy " << strategy_string(_jet_def.strategy())
330 <<
" automatically changed to " << strategy_string()
331 <<
" because the former is not supported for R = " << _Rparam
333 _changed_strategy_warning.warn(oss.str());
343 if (_strategy == N2Plain) {
345 this->_simple_N2_cluster_BriefJet();
346 }
else if (_strategy == N2Tiled) {
347 this->_faster_tiled_N2_cluster();
348 }
else if (_strategy == N2MinHeapTiled) {
349 this->_minheap_faster_tiled_N2_cluster();
350 }
else if (_strategy == N2MHTLazy9Alt) {
353 _plugin_activated =
true;
354 LazyTiling9Alt tiling(*
this);
356 _plugin_activated =
false;
358 }
else if (_strategy == N2MHTLazy25) {
361 _plugin_activated =
true;
362 LazyTiling25 tiling(*
this);
364 _plugin_activated =
false;
366 }
else if (_strategy == N2MHTLazy9) {
369 _plugin_activated =
true;
370 LazyTiling9 tiling(*
this);
372 _plugin_activated =
false;
374 }
else if (_strategy == N2MHTLazy9AntiKtSeparateGhosts) {
378 _plugin_activated =
true;
379 LazyTiling9SeparateGhosts tiling(*
this);
381 _plugin_activated =
false;
383 throw Error(
"N2MHTLazy9AntiKtSeparateGhosts strategy not supported with FJCORE");
386 }
else if (_strategy == NlnN) {
387 this->_delaunay_cluster();
388 }
else if (_strategy == NlnNCam) {
389 this->_CP2DChan_cluster_2piMultD();
390 }
else if (_strategy == NlnN3pi || _strategy == NlnN4pi ) {
391 this->_delaunay_cluster();
392 }
else if (_strategy == N3Dumb ) {
393 this->_really_dumb_cluster();
394 }
else if (_strategy == N2PoorTiled) {
395 this->_tiled_N2_cluster();
396 }
else if (_strategy == NlnNCam4pi) {
397 this->_CP2DChan_cluster();
398 }
else if (_strategy == NlnNCam2pi2R) {
399 this->_CP2DChan_cluster_2pi2R();
402 err <<
"Unrecognised value for strategy: "<<_strategy;
403 throw Error(err.str());
410 bool ClusterSequence::_first_time =
true;
411 LimitedWarning ClusterSequence::_exclusive_warnings;
417 return "FastJet version "+string(fastjet_version);
423 void ClusterSequence::print_banner() {
425 if (!_first_time) {
return;}
429 ostream * ostr = _fastjet_banner_ostr;
432 (*ostr) <<
"#--------------------------------------------------------------------------\n";
433 (*ostr) <<
"# FastJet release " << fastjet_version << endl;
434 (*ostr) <<
"# M. Cacciari, G.P. Salam and G. Soyez \n";
435 (*ostr) <<
"# A software package for jet finding and analysis at colliders \n";
436 (*ostr) <<
"# http://fastjet.fr \n";
438 (*ostr) <<
"# Please cite EPJC72(2012)1896 [arXiv:1111.6097] if you use this package\n";
439 (*ostr) <<
"# for scientific work and optionally PLB641(2006)57 [hep-ph/0512210]. \n";
441 (*ostr) <<
"# FastJet is provided without warranty under the terms of the GNU GPLv2.\n";
442 (*ostr) <<
"# It uses T. Chan's closest pair algorithm, S. Fortune's Voronoi code";
444 (*ostr) <<
",\n# CGAL ";
448 (*ostr) <<
"and 3rd party plugin jet algorithms. See COPYING file for details.\n";
449 (*ostr) <<
"#--------------------------------------------------------------------------\n";
457 const bool & writeout_combinations) {
459 _jet_def = jet_def_in;
460 _writeout_combinations = writeout_combinations;
464 _decant_options_partial();
469 void ClusterSequence::_decant_options_partial() {
473 _jet_algorithm = _jet_def.jet_algorithm();
474 _Rparam = _jet_def.R(); _R2 = _Rparam*_Rparam; _invR2 = 1.0/_R2;
475 _strategy = _jet_def.strategy();
478 _plugin_activated =
false;
482 _update_structure_use_count();
488 void ClusterSequence::_fill_initial_history () {
493 _jets.reserve(_jets.size()*2);
494 _history.reserve(_jets.size()*2);
498 for (
int i = 0; i < static_cast<int>(_jets.size()) ; i++) {
500 element.parent1 = InexistentParent;
501 element.
parent2 = InexistentParent;
502 element.
child = Invalid;
507 _history.push_back(element);
510 _jet_def.recombiner()->preprocess(_jets[i]);
513 _jets[i].set_cluster_hist_index(i);
514 _set_structure_shared_ptr(_jets[i]);
517 _Qtot += _jets[i].E();
519 _initial_n = _jets.size();
520 _deletes_self_when_unused =
false;
525 string ClusterSequence::strategy_string (
Strategy strategy_in)
const {
527 switch(strategy_in) {
529 strategy =
"NlnN";
break;
531 strategy =
"NlnN3pi";
break;
533 strategy =
"NlnN4pi";
break;
535 strategy =
"N2Plain";
break;
537 strategy =
"N2Tiled";
break;
539 strategy =
"N2MinHeapTiled";
break;
541 strategy =
"N2PoorTiled";
break;
543 strategy =
"N2MHTLazy9";
break;
545 strategy =
"N2MHTLazy9Alt";
break;
547 strategy =
"N2MHTLazy25";
break;
549 strategy =
"N2MHTLazy9AntiKtSeparateGhosts";
break;
551 strategy =
"N3Dumb";
break;
553 strategy =
"NlnNCam4pi";
break;
555 strategy =
"NlnNCam2pi2R";
break;
557 strategy =
"NlnNCam";
break;
559 strategy =
"plugin strategy";
break;
561 strategy =
"Unrecognized";
567 double ClusterSequence::jet_scale_for_algorithm(
572 double kt2=jet.
kt2();
573 return kt2 > 1e-300 ? 1.0/kt2 : 1e300;
575 double kt2 = jet.
kt2();
576 double p = jet_def().extra_param();
577 if (p <= 0 && kt2 < 1e-300) kt2 = 1e-300;
580 double kt2 = jet.
kt2();
581 double lim = _jet_def.extra_param();
582 if (kt2 < lim*lim && kt2 != 0.0) {
585 }
else {
throw Error(
"Unrecognised jet algorithm");}
604 int N = _jets.size();
607 double bounded_R = max(_Rparam, 0.1);
611 if (N <= 30 || N <= 39.0/(bounded_R + 0.6)) {
622 const static _Parabola N_Tiled_to_MHT_lowR (-45.4947,54.3528,44.6283);
623 const static _Parabola L_MHT_to_MHTLazy9_lowR (0.677807,-1.05006,10.6994);
624 const static _Parabola L_MHTLazy9_to_MHTLazy25_akt_lowR(0.169967,-0.512589,12.1572);
625 const static _Parabola L_MHTLazy9_to_MHTLazy25_kt_lowR (0.16237,-0.484612,12.3373);
626 const static _Parabola L_MHTLazy9_to_MHTLazy25_cam_lowR = L_MHTLazy9_to_MHTLazy25_kt_lowR;
627 const static _Parabola L_MHTLazy25_to_NlnN_akt_lowR (0.0472051,-0.22043,15.9196);
628 const static _Parabola L_MHTLazy25_to_NlnN_kt_lowR (0.118609,-0.326811,14.8287);
629 const static _Parabola L_MHTLazy25_to_NlnN_cam_lowR (0.10119,-0.295748,14.3924);
631 const static _Line L_Tiled_to_MHTLazy9_medR (-1.31304,7.29621);
632 const static _Parabola L_MHTLazy9_to_MHTLazy25_akt_medR = L_MHTLazy9_to_MHTLazy25_akt_lowR;
633 const static _Parabola L_MHTLazy9_to_MHTLazy25_kt_medR = L_MHTLazy9_to_MHTLazy25_kt_lowR;
634 const static _Parabola L_MHTLazy9_to_MHTLazy25_cam_medR = L_MHTLazy9_to_MHTLazy25_cam_lowR;
635 const static _Parabola L_MHTLazy25_to_NlnN_akt_medR = L_MHTLazy25_to_NlnN_akt_lowR;
636 const static _Parabola L_MHTLazy25_to_NlnN_kt_medR = L_MHTLazy25_to_NlnN_kt_lowR;
637 const static _Parabola L_MHTLazy25_to_NlnN_cam_medR = L_MHTLazy25_to_NlnN_cam_lowR;
639 const static double N_Plain_to_MHTLazy9_largeR = 75;
640 const static double N_MHTLazy9_to_MHTLazy25_akt_largeR = 700;
641 const static double N_MHTLazy9_to_MHTLazy25_kt_largeR = 1000;
642 const static double N_MHTLazy9_to_MHTLazy25_cam_largeR = 1000;
643 const static double N_MHTLazy25_to_NlnN_akt_largeR = 100000;
644 const static double N_MHTLazy25_to_NlnN_kt_largeR = 40000;
645 const static double N_MHTLazy25_to_NlnN_cam_largeR = 15000;
656 double p = jet_def().extra_param();
664 jet_algorithm = _jet_algorithm;
667 if (bounded_R < 0.65) {
669 if (N < N_Tiled_to_MHT_lowR(bounded_R))
return N2Tiled;
670 double logN = log(
double(N));
671 if (logN < L_MHT_to_MHTLazy9_lowR(bounded_R))
return N2MinHeapTiled;
674 if (logN < L_MHTLazy9_to_MHTLazy25_akt_lowR(bounded_R))
return N2MHTLazy9;
675 else if (logN < L_MHTLazy25_to_NlnN_akt_lowR(bounded_R))
return N2MHTLazy25;
678 if (logN < L_MHTLazy9_to_MHTLazy25_kt_lowR(bounded_R))
return N2MHTLazy9;
679 else if (logN < L_MHTLazy25_to_NlnN_kt_lowR(bounded_R))
return N2MHTLazy25;
682 if (logN < L_MHTLazy9_to_MHTLazy25_cam_lowR(bounded_R))
return N2MHTLazy9;
683 else if (logN < L_MHTLazy25_to_NlnN_cam_lowR(bounded_R))
return N2MHTLazy25;
687 }
else if (bounded_R < 0.5*pi) {
689 double logN = log(
double(N));
690 if (logN < L_Tiled_to_MHTLazy9_medR(bounded_R))
return N2Tiled;
693 if (logN < L_MHTLazy9_to_MHTLazy25_akt_medR(bounded_R))
return N2MHTLazy9;
694 else if (logN < L_MHTLazy25_to_NlnN_akt_medR(bounded_R))
return N2MHTLazy25;
697 if (logN < L_MHTLazy9_to_MHTLazy25_kt_medR(bounded_R))
return N2MHTLazy9;
698 else if (logN < L_MHTLazy25_to_NlnN_kt_medR(bounded_R))
return N2MHTLazy25;
701 if (logN < L_MHTLazy9_to_MHTLazy25_cam_medR(bounded_R))
return N2MHTLazy9;
702 else if (logN < L_MHTLazy25_to_NlnN_cam_medR(bounded_R))
return N2MHTLazy25;
708 if (N < N_Plain_to_MHTLazy9_largeR)
return N2Plain;
711 if (N < N_MHTLazy9_to_MHTLazy25_akt_largeR)
return N2MHTLazy9;
712 else if (N < N_MHTLazy25_to_NlnN_akt_largeR)
return N2MHTLazy25;
715 if (N < N_MHTLazy9_to_MHTLazy25_kt_largeR)
return N2MHTLazy9;
716 else if (N < N_MHTLazy25_to_NlnN_kt_largeR)
return N2MHTLazy25;
719 if (N < N_MHTLazy9_to_MHTLazy25_cam_largeR)
return N2MHTLazy9;
720 else if (N < N_MHTLazy25_to_NlnN_cam_largeR)
return N2MHTLazy25;
726 bool code_should_never_reach_here =
false;
727 assert(code_should_never_reach_here);
790 if (will_delete_self_when_unused())
791 throw(
Error(
"cannot use CS::transfer_from_sequence after a call to delete_self_when_unused()"));
794 _jet_def = from_seq._jet_def ;
795 _writeout_combinations = from_seq._writeout_combinations ;
796 _initial_n = from_seq._initial_n ;
797 _Rparam = from_seq._Rparam ;
799 _invR2 = from_seq._invR2 ;
800 _strategy = from_seq._strategy ;
801 _jet_algorithm = from_seq._jet_algorithm ;
802 _plugin_activated = from_seq._plugin_activated ;
808 _jets = (*action_on_jets)(from_seq.
_jets);
810 _jets = from_seq.
_jets;
814 _extras = from_seq._extras;
817 if (_structure_shared_ptr()) {
821 if (_deletes_self_when_unused)
throw Error(
"transfer_from_sequence cannot be used for a cluster sequence that deletes self when unused");
832 _update_structure_use_count();
834 for (
unsigned int i=0; i<_jets.size(); i++){
837 _jets[i].set_cluster_hist_index(from_seq.
_jets[i].cluster_hist_index());
840 _set_structure_shared_ptr(_jets[i]);
848 void ClusterSequence::plugin_record_ij_recombination(
849 int jet_i,
int jet_j,
double dij,
850 const PseudoJet & newjet,
int & newjet_k) {
852 plugin_record_ij_recombination(jet_i, jet_j, dij, newjet_k);
855 int tmp_index = _jets[newjet_k].cluster_hist_index();
856 _jets[newjet_k] = newjet;
858 _set_structure_shared_ptr(_jets[newjet_k]);
864 vector<PseudoJet> ClusterSequence::inclusive_jets (
const double ptmin)
const{
865 double dcut = ptmin*ptmin;
866 int i = _history.size() - 1;
867 vector<PseudoJet> jets_local;
873 if (_history[i].max_dij_so_far < dcut) {
break;}
874 if (_history[i].parent2 == BeamJet && _history[i].dij >= dcut) {
876 int parent1 = _history[i].parent1;
877 jets_local.push_back(_jets[_history[parent1].jetp_index]);}
885 if (_history[i].parent2 != BeamJet) {
break;}
886 int parent1 = _history[i].parent1;
887 const PseudoJet & jet = _jets[_history[parent1].jetp_index];
888 if (jet.
perp2() >= dcut) {jets_local.push_back(jet);}
901 if (_history[i].parent2 == BeamJet) {
902 int parent1 = _history[i].parent1;
903 const PseudoJet & jet = _jets[_history[parent1].jetp_index];
904 if (jet.
perp2() >= dcut) {jets_local.push_back(jet);}
908 }
else {
throw Error(
"cs::inclusive_jets(...): Unrecognized jet algorithm");}
916 int ClusterSequence::n_exclusive_jets (
const double dcut)
const {
920 int i = _history.size() - 1;
922 if (_history[i].max_dij_so_far <= dcut) {
break;}
925 int stop_point = i + 1;
928 int njets = 2*_initial_n - stop_point;
935 vector<PseudoJet> ClusterSequence::exclusive_jets (
const double dcut)
const {
936 int njets = n_exclusive_jets(dcut);
937 return exclusive_jets(njets);
944 vector<PseudoJet> ClusterSequence::exclusive_jets (
const int njets)
const {
948 if (njets > _initial_n) {
950 err <<
"Requested " << njets <<
" exclusive jets, but there were only "
951 << _initial_n <<
" particles in the event";
952 throw Error(err.str());
955 return exclusive_jets_up_to(njets);
961 vector<PseudoJet> ClusterSequence::exclusive_jets_up_to (
const int njets)
const {
973 (_jet_def.extra_param() <0)) &&
975 (!_jet_def.plugin()->exclusive_sequence_meaningful()))) {
976 _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.");
983 int stop_point = 2*_initial_n - njets;
985 if (stop_point < _initial_n) stop_point = _initial_n;
989 if (2*_initial_n != static_cast<int>(_history.size())) {
991 err <<
"2*_initial_n != _history.size() -- this endangers internal assumptions!\n";
992 throw Error(err.str());
1001 vector<PseudoJet> jets_local;
1002 for (
unsigned int i = stop_point; i < _history.size(); i++) {
1003 int parent1 = _history[i].parent1;
1004 if (parent1 < stop_point) {
1005 jets_local.push_back(_jets[_history[parent1].jetp_index]);
1007 int parent2 = _history[i].parent2;
1008 if (parent2 < stop_point && parent2 > 0) {
1009 jets_local.push_back(_jets[_history[parent2].jetp_index]);
1015 if (
int(jets_local.size()) != min(_initial_n, njets)) {
1017 err <<
"ClusterSequence::exclusive_jets: size of returned vector ("
1018 <<jets_local.size()<<
") does not coincide with requested number of jets ("
1020 throw Error(err.str());
1029 double ClusterSequence::exclusive_dmerge (
const int njets)
const {
1031 if (njets >= _initial_n) {
return 0.0;}
1032 return _history[2*_initial_n-njets-1].dij;
1041 double ClusterSequence::exclusive_dmerge_max (
const int njets)
const {
1043 if (njets >= _initial_n) {
return 0.0;}
1044 return _history[2*_initial_n-njets-1].max_dij_so_far;
1052 std::vector<PseudoJet> ClusterSequence::exclusive_subjets
1055 set<const history_element*> subhist;
1059 get_subhist_set(subhist, jet, dcut, 0);
1062 vector<PseudoJet> subjets;
1063 subjets.reserve(subhist.size());
1064 for (set<const history_element*>::iterator elem = subhist.begin();
1065 elem != subhist.end(); elem++) {
1066 subjets.push_back(_jets[(*elem)->jetp_index]);
1075 int ClusterSequence::n_exclusive_subjets(
const PseudoJet & jet,
1076 const double dcut)
const {
1077 set<const history_element*> subhist;
1080 get_subhist_set(subhist, jet, dcut, 0);
1081 return subhist.size();
1088 std::vector<PseudoJet> ClusterSequence::exclusive_subjets
1090 vector<PseudoJet> subjets = exclusive_subjets_up_to(jet, nsub);
1091 if (
int(subjets.size()) < nsub) {
1093 err <<
"Requested " << nsub <<
" exclusive subjets, but there were only "
1094 << subjets.size() <<
" particles in the jet";
1095 throw Error(err.str());
1105 std::vector<PseudoJet> ClusterSequence::exclusive_subjets_up_to
1108 set<const history_element*> subhist;
1111 vector<PseudoJet> subjets;
1112 if (nsub < 0)
throw Error(
"Requested a negative number of subjets. This is nonsensical.");
1113 if (nsub == 0)
return subjets;
1117 get_subhist_set(subhist, jet, -1.0, nsub);
1120 subjets.reserve(subhist.size());
1121 for (set<const history_element*>::iterator elem = subhist.begin();
1122 elem != subhist.end(); elem++) {
1123 subjets.push_back(_jets[(*elem)->jetp_index]);
1134 double ClusterSequence::exclusive_subdmerge(
const PseudoJet & jet,
int nsub)
const {
1135 set<const history_element*> subhist;
1139 get_subhist_set(subhist, jet, -1.0, nsub);
1141 set<const history_element*>::iterator highest = subhist.end();
1145 return (*highest)->dij;
1155 double ClusterSequence::exclusive_subdmerge_max(
const PseudoJet & jet,
int nsub)
const {
1157 set<const history_element*> subhist;
1161 get_subhist_set(subhist, jet, -1.0, nsub);
1163 set<const history_element*>::iterator highest = subhist.end();
1167 return (*highest)->max_dij_so_far;
1179 void ClusterSequence::get_subhist_set(set<const history_element*> & subhist,
1181 double dcut,
int maxjet)
const {
1182 assert(contains(jet));
1192 set<const history_element*>::iterator highest = subhist.end();
1193 assert (highest != subhist.begin());
1197 if (njet == maxjet)
break;
1199 if (elem->parent1 < 0)
break;
1204 subhist.erase(highest);
1205 subhist.insert(&(_history[elem->parent1]));
1206 subhist.insert(&(_history[elem->
parent2]));
1213 bool ClusterSequence::object_in_jet(
const PseudoJet &
object,
1218 assert(contains(
object) && contains(jet));
1220 const PseudoJet * this_object = &object;
1225 }
else if (has_child(*this_object, childp)) {
1226 this_object = childp;
1246 assert ((hist.parent1 >= 0 && hist.
parent2 >= 0) ||
1247 (hist.parent1 < 0 && hist.
parent2 < 0));
1249 if (hist.parent1 < 0) {
1254 parent1 = _jets[_history[hist.parent1].jetp_index];
1255 parent2 = _jets[_history[hist.
parent2].jetp_index];
1257 if (parent1.
perp2() < parent2.
perp2()) std::swap(parent1,parent2);
1277 bool res = has_child(jet, childp);
1294 if (hist.
child >= 0 && _history[hist.
child].jetp_index >= 0) {
1295 childp = &(_jets[_history[hist.
child].jetp_index]);
1315 if (hist.
child >= 0 && _history[hist.
child].parent2 >= 0) {
1320 partner = _jets[_history[child_hist.
parent2].jetp_index];
1323 partner = _jets[_history[child_hist.parent1].jetp_index];
1335 vector<PseudoJet> ClusterSequence::constituents (
const PseudoJet & jet)
const {
1336 vector<PseudoJet> subjets;
1337 add_constituents(jet, subjets);
1350 void ClusterSequence::print_jets_for_root(
const std::vector<PseudoJet> & jets_in,
1351 ostream & ostr)
const {
1352 for (
unsigned i = 0; i < jets_in.size(); i++) {
1354 << jets_in[i].px() <<
" "
1355 << jets_in[i].py() <<
" "
1356 << jets_in[i].pz() <<
" "
1357 << jets_in[i].E() << endl;
1358 vector<PseudoJet> cst = constituents(jets_in[i]);
1359 for (
unsigned j = 0; j < cst.size() ; j++) {
1360 ostr <<
" " << j <<
" "
1361 << cst[j].rap() <<
" "
1362 << cst[j].phi() <<
" "
1363 << cst[j].perp() << endl;
1365 ostr <<
"#END" << endl;
1369 void ClusterSequence::print_jets_for_root(
const std::vector<PseudoJet> & jets_in,
1370 const std::string & filename,
1371 const std::string & comment )
const {
1372 std::ofstream ostr(filename.c_str());
1373 if (comment !=
"") ostr <<
"# " << comment << endl;
1374 print_jets_for_root(jets_in, ostr);
1395 vector<int> ClusterSequence::particle_jet_indices(
1396 const vector<PseudoJet> & jets_in)
const {
1398 vector<int> indices(n_particles());
1401 for (
unsigned ipart = 0; ipart < n_particles(); ipart++)
1402 indices[ipart] = -1;
1406 for (
unsigned ijet = 0; ijet < jets_in.size(); ijet++) {
1408 vector<PseudoJet> jet_constituents(constituents(jets_in[ijet]));
1410 for (
unsigned ip = 0; ip < jet_constituents.size(); ip++) {
1414 unsigned iclust = jet_constituents[ip].cluster_hist_index();
1415 unsigned ipart = history()[iclust].jetp_index;
1416 indices[ipart] = ijet;
1426 void ClusterSequence::add_constituents (
1427 const PseudoJet & jet, vector<PseudoJet> & subjet_vector)
const {
1430 int parent1 = _history[i].parent1;
1431 int parent2 = _history[i].parent2;
1433 if (parent1 == InexistentParent) {
1439 subjet_vector.push_back(_jets[i]);
1444 add_constituents(_jets[_history[parent1].jetp_index], subjet_vector);
1447 if (parent2 != BeamJet) {
1448 add_constituents(_jets[_history[parent2].jetp_index], subjet_vector);
1456 void ClusterSequence::_add_step_to_history (
1457 const int step_number,
const int parent1,
1458 const int parent2,
const int jetp_index,
1461 history_element element;
1462 element.parent1 = parent1;
1463 element.parent2 = parent2;
1464 element.jetp_index = jetp_index;
1465 element.child = Invalid;
1467 element.max_dij_so_far = max(dij,_history[_history.size()-1].max_dij_so_far);
1468 _history.push_back(element);
1470 int local_step = _history.size()-1;
1471 assert(local_step == step_number);
1480 assert(parent1 >= 0);
1481 if (_history[parent1].child != Invalid){
1482 throw InternalError(
"trying to recomine an object that has previsously been recombined");
1484 _history[parent1].child = local_step;
1486 if (_history[parent2].child != Invalid){
1487 throw InternalError(
"trying to recomine an object that has previsously been recombined");
1489 _history[parent2].child = local_step;
1493 if (jetp_index != Invalid) {
1494 assert(jetp_index >= 0);
1496 _jets[jetp_index].set_cluster_hist_index(local_step);
1497 _set_structure_shared_ptr(_jets[jetp_index]);
1500 if (_writeout_combinations) {
1501 cout << local_step <<
": "
1502 << parent1 <<
" with " << parent2
1503 <<
"; y = "<< dij<<endl;
1516 vector<int> ClusterSequence::unique_history_order()
const {
1522 valarray<int> lowest_constituent(_history.size());
1523 int hist_n = _history.size();
1524 lowest_constituent = hist_n;
1525 for (
int i = 0; i < hist_n; i++) {
1527 lowest_constituent[i] = min(lowest_constituent[i],i);
1529 if (_history[i].child > 0) lowest_constituent[_history[i].child]
1530 = min(lowest_constituent[_history[i].child],lowest_constituent[i]);
1534 valarray<bool> extracted(_history.size()); extracted =
false;
1535 vector<int> unique_tree;
1536 unique_tree.reserve(_history.size());
1539 for (
unsigned i = 0; i < n_particles(); i++) {
1540 if (!extracted[i]) {
1541 unique_tree.push_back(i);
1542 extracted[i] =
true;
1543 _extract_tree_children(i, extracted, lowest_constituent, unique_tree);
1552 void ClusterSequence::_extract_tree_children(
1554 valarray<bool> & extracted,
1555 const valarray<int> & lowest_constituent,
1556 vector<int> & unique_tree)
const {
1557 if (!extracted[position]) {
1560 _extract_tree_parents(position,extracted,lowest_constituent,unique_tree);
1564 int child = _history[position].child;
1565 if (child >= 0) _extract_tree_children(child,extracted,lowest_constituent,unique_tree);
1571 vector<PseudoJet> ClusterSequence::unclustered_particles()
const {
1572 vector<PseudoJet> unclustered;
1573 for (
unsigned i = 0; i < n_particles() ; i++) {
1574 if (_history[i].child == Invalid)
1575 unclustered.push_back(_jets[_history[i].jetp_index]);
1585 vector<PseudoJet> ClusterSequence::childless_pseudojets()
const {
1586 vector<PseudoJet> unclustered;
1587 for (
unsigned i = 0; i < _history.size() ; i++) {
1588 if ((_history[i].child == Invalid) && (_history[i].parent2 != BeamJet))
1589 unclustered.push_back(_jets[_history[i].jetp_index]);
1600 bool ClusterSequence::contains(
const PseudoJet & jet)
const {
1611 void ClusterSequence::_extract_tree_parents(
1613 valarray<bool> & extracted,
1614 const valarray<int> & lowest_constituent,
1615 vector<int> & unique_tree)
const {
1617 if (!extracted[position]) {
1618 int parent1 = _history[position].parent1;
1619 int parent2 = _history[position].parent2;
1622 if (parent1 >= 0 && parent2 >= 0) {
1623 if (lowest_constituent[parent1] > lowest_constituent[parent2])
1624 std::swap(parent1, parent2);
1627 if (parent1 >= 0 && !extracted[parent1])
1628 _extract_tree_parents(parent1,extracted,lowest_constituent,unique_tree);
1629 if (parent2 >= 0 && !extracted[parent2])
1630 _extract_tree_parents(parent2,extracted,lowest_constituent,unique_tree);
1633 unique_tree.push_back(position);
1634 extracted[position] =
true;
1643 void ClusterSequence::_do_ij_recombination_step(
1644 const int jet_i,
const int jet_j,
1654 _jet_def.recombiner()->recombine(_jets[jet_i], _jets[jet_j], newjet);
1655 _jets.push_back(newjet);
1660 newjet_k = _jets.size()-1;
1663 int newstep_k = _history.size();
1665 _jets[newjet_k].set_cluster_hist_index(newstep_k);
1668 int hist_i = _jets[jet_i].cluster_hist_index();
1669 int hist_j = _jets[jet_j].cluster_hist_index();
1671 _add_step_to_history(newstep_k, min(hist_i, hist_j), max(hist_i,hist_j),
1680 void ClusterSequence::_do_iB_recombination_step(
1681 const int jet_i,
const double diB) {
1683 int newstep_k = _history.size();
1686 _add_step_to_history(newstep_k,_jets[jet_i].cluster_hist_index(),BeamJet,
1697 void ClusterSequence::_set_structure_shared_ptr(
PseudoJet & j) {
1701 _update_structure_use_count();
1706 void ClusterSequence::_update_structure_use_count() {
1709 _structure_use_count_after_construction = _structure_shared_ptr.use_count();
1716 void ClusterSequence::delete_self_when_unused() {
1724 int new_count = _structure_shared_ptr.use_count() - _structure_use_count_after_construction;
1725 if (new_count <= 0) {
1726 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");
1729 _structure_shared_ptr.set_count(new_count);
1730 _deletes_self_when_unused =
true;
1734 FASTJET_END_NAMESPACE
a version of cambridge with a special distance measure for particles whose pt is < extra_param(); thi...
like the k_t but with distance measures dij = min(kti^{2p},ktj^{2p}) Delta R_{ij}^2 / R^2 diB = 1/kti...
double kt2() const
returns the squared transverse momentum
Chan's closest pair method (in a variant with 4pi coverage), for use exclusively with the Cambridge a...
best of the NlnN variants – best overall for N>10^4.
Similar to N2MHTLazy9, but uses tiles of size R/2 and a 5x5 tile grid around the particle.
legacy N ln N using 4pi coverage of cylinder
int cluster_hist_index() const
return the cluster_hist_index, intended to be used by clustering routines.
the e+e- genkt algorithm (R > 2 and p=1 gives ee_kt)
fastest from about 50..500
the longitudinally invariant kt algorithm
int jetp_index
index in _history where the current jet is recombined with another jet to form its child...
Contains any information related to the clustering that should be directly accessible to PseudoJet...
std::vector< PseudoJet > _jets
This contains the physical PseudoJets; for each PseudoJet one can find the corresponding position in ...
Chan's closest pair method (in a variant with 2pi+minimal extra variant), for use exclusively with th...
std::vector< history_element > _history
this vector will contain the branching history; for each stage, _history[i].jetp_index indicates wher...
class to provide facilities for giving warnings up to some maximum number of times and to provide glo...
any plugin algorithm supplied by the user
int parent2
index in _history where first parent of this jet was created (InexistentParent if this jet is an orig...
Like N2MHTLazy9 in a number of respects, but does not calculate ghost-ghost distances and so does not...
the plugin has been used...
virtual void set_associated_cs(const ClusterSequence *new_cs)
set the associated csw
double dij
index in the _jets vector where we will find the
void set_cluster_hist_index(const int index)
set the cluster_hist_index, intended to be used by clustering routines.
bool has_valid_cluster_sequence() const
returns true if this PseudoJet has an associated and still valid(ated) ClusterSequence.
worse even than the usual N^3 algorithms
Like to N2MHTLazy9 but uses slightly different optimizations, e.g.
base class corresponding to errors that can be thrown by FastJet
string fastjet_version_string()
return a string containing information about the release
int child
index in _history where second parent of this jet was created (InexistentParent if this jet is an ori...
automatic selection of the best (based on N), including the LazyTiled strategies that are new to FJ3...
class corresponding to critical internal errors
a single element in the clustering history
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...
const ClusterSequence * associated_cluster_sequence() const
get a (const) pointer to the parent ClusterSequence (NULL if inexistent)
void set_structure_shared_ptr(const SharedPtr< PseudoJetStructureBase > &structure)
set the associated structure
the longitudinally invariant variant of the cambridge algorithm (aka Aachen algoithm).
JetAlgorithm
the various families of jet-clustering algorithm
only looks into a neighbouring tile for a particle's nearest neighbour (NN) if that particle's in-til...
Strategy
the various options for the algorithmic strategy to adopt in clustering events with kt and cambridge ...
Class to contain pseudojets, including minimal information of use to jet-clustering routines...
double max_dij_so_far
the distance corresponding to the recombination at this stage of the clustering.
legacy N ln N using 3pi coverage of cylinder.
faster that N2Tiled above about 500 particles; differs from it by retainig the di(closest j) distance...
Chan's closest pair method (in a variant with 2pi+2R coverage), for use exclusively with the Cambridg...
double perp2() const
returns the squared transverse momentum
class that is intended to hold a full definition of the jet clusterer