32 #include "fastjet/internal/LazyTiling25.hh"
33 #include "fastjet/internal/TilingExtent.hh"
37 FASTJET_BEGIN_NAMESPACE
39 LazyTiling25::LazyTiling25(ClusterSequence & cs) :
40 _cs(cs), _jets(cs.jets())
47 _Rparam = cs.jet_def().R();
48 _R2 = _Rparam * _Rparam;
75 void LazyTiling25::_initialise_tiles() {
79 double default_size = max(0.1,_Rparam)/2;
80 _tile_size_eta = default_size;
84 _n_tiles_phi = max(5,
int(floor(twopi/default_size)));
85 _tile_size_phi = twopi / _n_tiles_phi;
87 #define _FASTJET_TILING25_USE_TILING_ANALYSIS_
88 #ifdef _FASTJET_TILING25_USE_TILING_ANALYSIS_
90 TilingExtent tiling_analysis(_cs);
91 _tiles_eta_min = tiling_analysis.minrap();
92 _tiles_eta_max = tiling_analysis.maxrap();
94 #else // not _FASTJET_TILING25_USE_TILING_ANALYSIS_
99 const double maxrap = 7.0;
102 for(
unsigned int i = 0; i < _jets.size(); i++) {
103 double eta = _jets[i].rap();
106 if (abs(eta) < maxrap) {
107 if (eta < _tiles_eta_min) {_tiles_eta_min = eta;}
108 if (eta > _tiles_eta_max) {_tiles_eta_max = eta;}
112 #endif // _FASTJET_TILING25_USE_TILING_ANALYSIS_
124 #define FASTJET_LAZY25_MIN3TILESY
125 #ifdef FASTJET_LAZY25_MIN3TILESY
126 if (_tiles_eta_max - _tiles_eta_min < 3*_tile_size_eta) {
131 _tile_size_eta = (_tiles_eta_max - _tiles_eta_min)/3;
136 _tiles_eta_max -= _tile_size_eta;
138 #endif //FASTJET_LAZY25_MIN3TILESY
139 _tiles_ieta_min = int(floor(_tiles_eta_min/_tile_size_eta));
140 _tiles_ieta_max = int(floor( _tiles_eta_max/_tile_size_eta));
141 _tiles_eta_min = _tiles_ieta_min * _tile_size_eta;
142 _tiles_eta_max = _tiles_ieta_max * _tile_size_eta;
143 #ifdef FASTJET_LAZY25_MIN3TILESY
146 _tile_half_size_eta = _tile_size_eta * 0.5;
147 _tile_half_size_phi = _tile_size_phi * 0.5;
151 vector<bool> use_periodic_delta_phi(_n_tiles_phi,
false);
152 if (_n_tiles_phi <= 5) {
153 fill(use_periodic_delta_phi.begin(), use_periodic_delta_phi.end(),
true);
155 use_periodic_delta_phi[0] =
true;
156 use_periodic_delta_phi[1] =
true;
157 use_periodic_delta_phi[_n_tiles_phi-2] =
true;
158 use_periodic_delta_phi[_n_tiles_phi-1] =
true;
162 _tiles.resize((_tiles_ieta_max-_tiles_ieta_min+1)*_n_tiles_phi);
165 for (
int ieta = _tiles_ieta_min; ieta <= _tiles_ieta_max; ieta++) {
166 for (
int iphi = 0; iphi < _n_tiles_phi; iphi++) {
167 Tile25 * tile = & _tiles[_tile_index(ieta,iphi)];
170 tile->begin_tiles[0] = tile;
171 Tile25 ** pptile = & (tile->begin_tiles[0]);
175 tile->surrounding_tiles = pptile;
176 if (ieta > _tiles_ieta_min) {
180 for (
int idphi = -2; idphi <=+2; idphi++) {
181 *pptile = & _tiles[_tile_index(ieta-1,iphi+idphi)];
186 if (ieta > _tiles_ieta_min + 1) {
190 for (
int idphi = -2; idphi <= +2; idphi++) {
191 *pptile = & _tiles[_tile_index(ieta-2,iphi+idphi)];
196 *pptile = & _tiles[_tile_index(ieta,iphi-1)];
198 *pptile = & _tiles[_tile_index(ieta,iphi-2)];
201 tile->RH_tiles = pptile;
202 *pptile = & _tiles[_tile_index(ieta,iphi+1)];
204 *pptile = & _tiles[_tile_index(ieta,iphi+2)];
207 if (ieta < _tiles_ieta_max) {
208 for (
int idphi = -2; idphi <= +2; idphi++) {
209 *pptile = & _tiles[_tile_index(ieta+1,iphi+idphi)];
214 if (ieta < _tiles_ieta_max - 1) {
215 for (
int idphi = -2; idphi <= +2; idphi++) {
216 *pptile = & _tiles[_tile_index(ieta+2,iphi+idphi)];
221 tile->end_tiles = pptile;
223 tile->tagged =
false;
225 tile->use_periodic_delta_phi = use_periodic_delta_phi[iphi];
227 tile->max_NN_dist = 0;
229 tile->eta_centre = (ieta-_tiles_ieta_min+0.5)*_tile_size_eta + _tiles_eta_min;
230 tile->phi_centre = (iphi+0.5)*_tile_size_phi;
238 int LazyTiling25::_tile_index(
const double eta,
const double phi)
const {
240 if (eta <= _tiles_eta_min) {ieta = 0;}
241 else if (eta >= _tiles_eta_max) {ieta = _tiles_ieta_max-_tiles_ieta_min;}
244 ieta = int(((eta - _tiles_eta_min) / _tile_size_eta));
246 if (ieta > _tiles_ieta_max-_tiles_ieta_min) {
247 ieta = _tiles_ieta_max-_tiles_ieta_min;}
253 iphi = int((phi+twopi)/_tile_size_phi) % _n_tiles_phi;
254 return (iphi + ieta * _n_tiles_phi);
260 inline void LazyTiling25::_tj_set_jetinfo( TiledJet *
const jet,
261 const int _jets_index) {
263 _bj_set_jetinfo<>(jet, _jets_index);
268 jet->tile_index = _tile_index(jet->eta, jet->phi);
271 Tile25 * tile = &_tiles[jet->tile_index];
272 jet->previous = NULL;
273 jet->next = tile->head;
274 if (jet->next != NULL) {jet->next->previous = jet;}
280 void LazyTiling25::_bj_remove_from_tiles(TiledJet *
const jet) {
281 Tile25 * tile = & _tiles[jet->tile_index];
283 if (jet->previous == NULL) {
286 tile->head = jet->next;
289 jet->previous->next = jet->next;
291 if (jet->next != NULL) {
293 jet->next->previous = jet->previous;
300 void LazyTiling25::_print_tiles(TiledJet * briefjets )
const {
301 for (vector<Tile25>::const_iterator tile = _tiles.begin();
302 tile < _tiles.end(); tile++) {
303 cout <<
"Tile " << tile - _tiles.begin()
304 <<
" at " << setw(10) << tile->eta_centre <<
"," << setw(10) << tile->phi_centre
307 for (TiledJet * jetI = tile->head; jetI != NULL; jetI = jetI->next) {
308 list.push_back(jetI-briefjets);
311 sort(list.begin(),list.end());
312 for (
unsigned int i = 0; i < list.size(); i++) {cout <<
" "<<list[i];}
325 void LazyTiling25::_add_neighbours_to_tile_union(
const int tile_index,
326 vector<int> & tile_union,
int & n_near_tiles)
const {
327 for (Tile25 *
const * near_tile = _tiles[tile_index].begin_tiles;
328 near_tile != _tiles[tile_index].end_tiles; near_tile++){
330 tile_union[n_near_tiles] = *near_tile - & _tiles[0];
340 inline void LazyTiling25::_add_untagged_neighbours_to_tile_union(
341 const int tile_index,
342 vector<int> & tile_union,
int & n_near_tiles) {
343 for (Tile25 ** near_tile = _tiles[tile_index].begin_tiles;
344 near_tile != _tiles[tile_index].end_tiles; near_tile++){
345 if (! (*near_tile)->tagged) {
346 (*near_tile)->tagged =
true;
348 tile_union[n_near_tiles] = *near_tile - & _tiles[0];
360 inline void LazyTiling25::_add_untagged_neighbours_to_tile_union_using_max_info(
361 const TiledJet * jet,
362 vector<int> & tile_union,
int & n_near_tiles) {
363 Tile25 & tile = _tiles[jet->tile_index];
365 for (Tile25 ** near_tile = tile.begin_tiles; near_tile != tile.end_tiles; near_tile++){
366 if ((*near_tile)->tagged)
continue;
367 double dist = _distance_to_tile(jet, *near_tile);
371 if (dist > (*near_tile)->max_NN_dist)
continue;
374 (*near_tile)->tagged =
true;
376 tile_union[n_near_tiles] = *near_tile - & _tiles[0];
394 inline double LazyTiling25::_distance_to_tile(
const TiledJet * bj,
const Tile25 * tile)
400 #endif // INSTRUMENT2
407 if (_tiles[bj->tile_index].eta_centre == tile->eta_centre) deta = 0;
409 else deta = std::abs(bj->eta - tile->eta_centre) - _tile_half_size_eta;
418 double dphi = std::abs(bj->phi - tile->phi_centre);
419 if (dphi > pi) dphi = twopi-dphi;
420 dphi -= _tile_half_size_phi;
422 if (dphi < 0) dphi = 0;
424 return dphi*dphi + deta*deta;
438 inline void LazyTiling25::_update_jetX_jetI_NN(TiledJet * jetX, TiledJet * jetI, vector<TiledJet *> & jets_for_minheap) {
439 double dist = _bj_dist(jetI,jetX);
440 if (dist < jetI->NN_dist) {
442 jetI->NN_dist = dist;
445 if (!jetI->minheap_update_needed()) {
446 jetI->label_minheap_update_needed();
447 jets_for_minheap.push_back(jetI);
451 if (dist < jetX->NN_dist) {
453 jetX->NN_dist = dist;
459 inline void LazyTiling25::_set_NN(TiledJet * jetI,
460 vector<TiledJet *> & jets_for_minheap) {
464 if (!jetI->minheap_update_needed()) {
465 jetI->label_minheap_update_needed();
466 jets_for_minheap.push_back(jetI);}
468 Tile25 * tile_ptr = &_tiles[jetI->tile_index];
470 for (Tile25 ** near_tile = tile_ptr->begin_tiles;
471 near_tile != tile_ptr->end_tiles; near_tile++) {
474 if (jetI->NN_dist < _distance_to_tile(jetI, *near_tile))
continue;
476 for (TiledJet * jetJ = (*near_tile)->head;
477 jetJ != NULL; jetJ = jetJ->next) {
478 double dist = _bj_dist(jetI,jetJ);
479 if (dist < jetI->NN_dist && jetJ != jetI) {
480 jetI->NN_dist = dist; jetI->NN = jetJ;
505 void LazyTiling25::run() {
509 int n = _jets.size();
512 TiledJet * briefjets =
new TiledJet[n];
513 TiledJet * jetA = briefjets, * jetB;
516 TiledJet oldB = briefjets[0];
520 vector<int> tile_union(3*25);
523 for (
int i = 0; i< n; i++) {
524 _tj_set_jetinfo(jetA, i);
528 TiledJet * head = briefjets;
531 vector<Tile25>::iterator tile;
532 for (tile = _tiles.begin(); tile != _tiles.end(); tile++) {
534 for (jetA = tile->head; jetA != NULL; jetA = jetA->next) {
535 for (jetB = tile->head; jetB != jetA; jetB = jetB->next) {
536 double dist = _bj_dist_not_periodic(jetA,jetB);
537 if (dist < jetA->NN_dist) {jetA->NN_dist = dist; jetA->NN = jetB;}
538 if (dist < jetB->NN_dist) {jetB->NN_dist = dist; jetB->NN = jetA;}
541 for (jetA = tile->head; jetA != NULL; jetA = jetA->next) {
542 if (jetA->NN_dist > tile->max_NN_dist) tile->max_NN_dist = jetA->NN_dist;
545 for (tile = _tiles.begin(); tile != _tiles.end(); tile++) {
546 if (tile->use_periodic_delta_phi) {
548 for (Tile25 ** RTile = tile->RH_tiles; RTile != tile->end_tiles; RTile++) {
549 for (jetA = tile->head; jetA != NULL; jetA = jetA->next) {
550 double dist_to_tile = _distance_to_tile(jetA, *RTile);
562 bool relevant_for_jetA = dist_to_tile <= jetA->NN_dist;
563 bool relevant_for_RTile = dist_to_tile <= (*RTile)->max_NN_dist;
564 if (relevant_for_jetA || relevant_for_RTile) {
565 for (jetB = (*RTile)->head; jetB != NULL; jetB = jetB->next) {
566 double dist = _bj_dist(jetA,jetB);
567 if (dist < jetA->NN_dist) {jetA->NN_dist = dist; jetA->NN = jetB;}
568 if (dist < jetB->NN_dist) {jetB->NN_dist = dist; jetB->NN = jetA;}
577 for (Tile25 ** RTile = tile->RH_tiles; RTile != tile->end_tiles; RTile++) {
578 for (jetA = tile->head; jetA != NULL; jetA = jetA->next) {
579 double dist_to_tile = _distance_to_tile(jetA, *RTile);
580 bool relevant_for_jetA = dist_to_tile <= jetA->NN_dist;
581 bool relevant_for_RTile = dist_to_tile <= (*RTile)->max_NN_dist;
582 if (relevant_for_jetA || relevant_for_RTile) {
583 for (jetB = (*RTile)->head; jetB != NULL; jetB = jetB->next) {
584 double dist = _bj_dist_not_periodic(jetA,jetB);
585 if (dist < jetA->NN_dist) {jetA->NN_dist = dist; jetA->NN = jetB;}
586 if (dist < jetB->NN_dist) {jetB->NN_dist = dist; jetB->NN = jetA;}
598 for (tile = _tiles.begin(); tile != _tiles.end(); tile++) {
599 tile->max_NN_dist = 0;
600 for (jetA = tile->head; jetA != NULL; jetA = jetA->next) {
601 if (jetA->NN_dist > tile->max_NN_dist) tile->max_NN_dist = jetA->NN_dist;
606 cout <<
"intermediate ncall, dtt = " << _ncall <<
" " << _ncall_dtt << endl;
607 #endif // INSTRUMENT2
615 vector<double> diJs(n);
616 for (
int i = 0; i < n; i++) {
617 diJs[i] = _bj_diJ(&briefjets[i]);
618 briefjets[i].label_minheap_update_done();
620 MinHeap minheap(diJs);
622 vector<TiledJet *> jets_for_minheap;
623 jets_for_minheap.reserve(n);
626 int history_location = n-1;
629 double diJ_min = minheap.minval() *_invR2;
630 jetA = head + minheap.minloc();
642 if (jetA < jetB) {std::swap(jetA,jetB);}
645 _cs.plugin_record_ij_recombination(jetA->_jets_index, jetB->_jets_index, diJ_min, nn);
648 _bj_remove_from_tiles(jetA);
650 _bj_remove_from_tiles(jetB);
651 _tj_set_jetinfo(jetB, nn);
656 _cs.plugin_record_iB_recombination(jetA->_jets_index, diJ_min);
657 _bj_remove_from_tiles(jetA);
661 minheap.remove(jetA-head);
663 int n_near_tiles = 0;
669 Tile25 & jetB_tile = _tiles[jetB->tile_index];
670 for (Tile25 ** near_tile = jetB_tile.begin_tiles;
671 near_tile != jetB_tile.end_tiles; near_tile++) {
673 double dist_to_tile = _distance_to_tile(jetB, *near_tile);
676 bool relevant_for_jetB = dist_to_tile <= jetB->NN_dist;
677 bool relevant_for_near_tile = dist_to_tile <= (*near_tile)->max_NN_dist;
678 bool relevant = relevant_for_jetB || relevant_for_near_tile;
679 if (! relevant)
continue;
682 tile_union[n_near_tiles] = *near_tile - & _tiles[0];
683 (*near_tile)->tagged =
true;
690 for (TiledJet * jetI = (*near_tile)->head; jetI != NULL; jetI = jetI->next) {
691 if (jetI->NN == jetA || jetI->NN == jetB) _set_NN(jetI, jets_for_minheap);
692 _update_jetX_jetI_NN(jetB, jetI, jets_for_minheap);
701 int n_done_tiles = n_near_tiles;
702 _add_untagged_neighbours_to_tile_union_using_max_info(jetA,
703 tile_union, n_near_tiles);
705 _add_untagged_neighbours_to_tile_union_using_max_info(&oldB,
706 tile_union,n_near_tiles);
707 jetB->label_minheap_update_needed();
708 jets_for_minheap.push_back(jetB);
713 for (
int itile = 0; itile < n_done_tiles; itile++) {
714 _tiles[tile_union[itile]].tagged =
false;
718 for (
int itile = n_done_tiles; itile < n_near_tiles; itile++) {
719 Tile25 * tile_ptr = &_tiles[tile_union[itile]];
720 tile_ptr->tagged =
false;
722 for (TiledJet * jetI = tile_ptr->head; jetI != NULL; jetI = jetI->next) {
724 if (jetI->NN == jetA || (jetI->NN == jetB && jetB != NULL)) {
725 _set_NN(jetI, jets_for_minheap);
732 while (jets_for_minheap.size() > 0) {
733 TiledJet * jetI = jets_for_minheap.back();
734 jets_for_minheap.pop_back();
735 minheap.update(jetI-head, _bj_diJ(jetI));
736 jetI->label_minheap_update_done();
739 Tile25 & tile_I = _tiles[jetI->tile_index];
740 if (tile_I.max_NN_dist < jetI->NN_dist) tile_I.max_NN_dist = jetI->NN_dist;
748 cout <<
"ncall, dtt = " << _ncall <<
" " << _ncall_dtt << endl;
749 #endif // INSTRUMENT2
753 FASTJET_END_NAMESPACE