32 #include "fastjet/internal/LazyTiling9Alt.hh"
33 #include "fastjet/internal/TilingExtent.hh"
36 FASTJET_BEGIN_NAMESPACE
38 LazyTiling9Alt::LazyTiling9Alt(ClusterSequence & cs) :
39 _cs(cs), _jets(cs.jets())
42 _Rparam = cs.jet_def().R();
43 _R2 = _Rparam * _Rparam;
68 void LazyTiling9Alt::_initialise_tiles() {
72 double default_size = max(0.1,_Rparam);
73 _tile_size_eta = default_size;
77 _n_tiles_phi = max(3,
int(floor(twopi/default_size)));
78 _tile_size_phi = twopi / _n_tiles_phi;
84 const double maxrap = 7.0;
87 for(
unsigned int i = 0; i < _jets.size(); i++) {
88 double eta = _jets[i].rap();
91 if (abs(eta) < maxrap) {
92 if (eta < _tiles_eta_min) {_tiles_eta_min = eta;}
93 if (eta > _tiles_eta_max) {_tiles_eta_max = eta;}
98 _tiles_ieta_min = int(floor(_tiles_eta_min/_tile_size_eta));
99 _tiles_ieta_max = int(floor( _tiles_eta_max/_tile_size_eta));
100 _tiles_eta_min = _tiles_ieta_min * _tile_size_eta;
101 _tiles_eta_max = _tiles_ieta_max * _tile_size_eta;
103 _tile_half_size_eta = _tile_size_eta * 0.5;
104 _tile_half_size_phi = _tile_size_phi * 0.5;
108 vector<bool> use_periodic_delta_phi(_n_tiles_phi,
false);
109 if (_n_tiles_phi <= 3) {
110 fill(use_periodic_delta_phi.begin(), use_periodic_delta_phi.end(),
true);
112 use_periodic_delta_phi[0] =
true;
113 use_periodic_delta_phi[_n_tiles_phi-1] =
true;
117 _tiles.resize((_tiles_ieta_max-_tiles_ieta_min+1)*_n_tiles_phi);
120 for (
int ieta = _tiles_ieta_min; ieta <= _tiles_ieta_max; ieta++) {
121 for (
int iphi = 0; iphi < _n_tiles_phi; iphi++) {
122 Tile * tile = & _tiles[_tile_index(ieta,iphi)];
125 tile->begin_tiles[0] = Tile::TileFnPair(tile,&Tile::distance_to_centre);
126 Tile::TileFnPair * pptile = & (tile->begin_tiles[0]);
130 tile->surrounding_tiles = pptile;
131 if (ieta > _tiles_ieta_min) {
139 *pptile = Tile::TileFnPair(& _tiles[_tile_index(ieta-1,iphi-1)],
140 &Tile::distance_to_left_bottom);
142 *pptile = Tile::TileFnPair(& _tiles[_tile_index(ieta-1,iphi)],
143 &Tile::distance_to_left);
145 *pptile = Tile::TileFnPair(& _tiles[_tile_index(ieta-1,iphi+1)],
146 &Tile::distance_to_left_top);
150 *pptile = Tile::TileFnPair(& _tiles[_tile_index(ieta,iphi-1)],
151 &Tile::distance_to_bottom);
154 tile->RH_tiles = pptile;
155 *pptile = Tile::TileFnPair(& _tiles[_tile_index(ieta,iphi+1)],
156 &Tile::distance_to_top);
159 if (ieta < _tiles_ieta_max) {
164 *pptile = Tile::TileFnPair(& _tiles[_tile_index(ieta+1,iphi-1)],
165 &Tile::distance_to_right_bottom);
167 *pptile = Tile::TileFnPair(& _tiles[_tile_index(ieta+1,iphi)],
168 &Tile::distance_to_right);
170 *pptile = Tile::TileFnPair(& _tiles[_tile_index(ieta+1,iphi+1)],
171 &Tile::distance_to_right_top);
175 tile->end_tiles = pptile;
177 tile->tagged =
false;
179 tile->use_periodic_delta_phi = use_periodic_delta_phi[iphi];
181 tile->max_NN_dist = 0;
183 tile->eta_min = ieta*_tile_size_eta;
184 tile->eta_max = (ieta+1)*_tile_size_eta;
185 tile->phi_min = iphi*_tile_size_phi;
186 tile->phi_max = (iphi+1)*_tile_size_phi;
194 int LazyTiling9Alt::_tile_index(
const double eta,
const double phi)
const {
196 if (eta <= _tiles_eta_min) {ieta = 0;}
197 else if (eta >= _tiles_eta_max) {ieta = _tiles_ieta_max-_tiles_ieta_min;}
200 ieta = int(((eta - _tiles_eta_min) / _tile_size_eta));
202 if (ieta > _tiles_ieta_max-_tiles_ieta_min) {
203 ieta = _tiles_ieta_max-_tiles_ieta_min;}
209 iphi = int((phi+twopi)/_tile_size_phi) % _n_tiles_phi;
210 return (iphi + ieta * _n_tiles_phi);
216 inline void LazyTiling9Alt::_tj_set_jetinfo( TiledJet *
const jet,
217 const int _jets_index) {
219 _bj_set_jetinfo<>(jet, _jets_index);
224 jet->tile_index = _tile_index(jet->eta, jet->phi);
227 Tile * tile = &_tiles[jet->tile_index];
228 jet->previous = NULL;
229 jet->next = tile->head;
230 if (jet->next != NULL) {jet->next->previous = jet;}
236 void LazyTiling9Alt::_bj_remove_from_tiles(TiledJet *
const jet) {
237 Tile * tile = & _tiles[jet->tile_index];
239 if (jet->previous == NULL) {
242 tile->head = jet->next;
245 jet->previous->next = jet->next;
247 if (jet->next != NULL) {
249 jet->next->previous = jet->previous;
256 void LazyTiling9Alt::_print_tiles(TiledJet * briefjets )
const {
257 for (vector<Tile>::const_iterator tile = _tiles.begin();
258 tile < _tiles.end(); tile++) {
259 cout <<
"Tile " << tile - _tiles.begin()<<
" = ";
261 for (TiledJet * jetI = tile->head; jetI != NULL; jetI = jetI->next) {
262 list.push_back(jetI-briefjets);
265 sort(list.begin(),list.end());
266 for (
unsigned int i = 0; i < list.size(); i++) {cout <<
" "<<list[i];}
279 void LazyTiling9Alt::_add_neighbours_to_tile_union(
const int tile_index,
280 vector<int> & tile_union,
int & n_near_tiles)
const {
281 for (Tile::TileFnPair
const * near_tile = _tiles[tile_index].begin_tiles;
282 near_tile != _tiles[tile_index].end_tiles; near_tile++){
284 tile_union[n_near_tiles] = near_tile->first - & _tiles[0];
294 inline void LazyTiling9Alt::_add_untagged_neighbours_to_tile_union(
295 const int tile_index,
296 vector<int> & tile_union,
int & n_near_tiles) {
297 for (Tile::TileFnPair * near_tile = _tiles[tile_index].begin_tiles;
298 near_tile != _tiles[tile_index].end_tiles; near_tile++){
299 if (! (near_tile->first)->tagged) {
300 (near_tile->first)->tagged =
true;
302 tile_union[n_near_tiles] = near_tile->first - & _tiles[0];
314 inline void LazyTiling9Alt::_add_untagged_neighbours_to_tile_union_using_max_info(
315 const TiledJet * jet,
316 vector<int> & tile_union,
int & n_near_tiles) {
317 Tile & tile = _tiles[jet->tile_index];
319 for (Tile::TileFnPair * near_tile = tile.begin_tiles; near_tile != tile.end_tiles; near_tile++){
320 if ((near_tile->first)->tagged)
continue;
321 double dist = (tile.*(near_tile->second))(jet);
325 if (dist > (near_tile->first)->max_NN_dist)
continue;
328 (near_tile->first)->tagged =
true;
330 tile_union[n_near_tiles] = near_tile->first - & _tiles[0];
336 ostream &
operator<<(ostream & ostr,
const TiledJet & jet) {
337 ostr <<
"j" << setw(3) << jet._jets_index <<
":pt2,rap,phi=" ; ostr.flush();
338 ostr << jet.kt2 <<
","; ostr.flush();
339 ostr << jet.eta <<
","; ostr.flush();
340 ostr << jet.phi; ostr.flush();
341 ostr <<
", tile=" << jet.tile_index; ostr.flush();
389 inline void LazyTiling9Alt::_update_jetX_jetI_NN(TiledJet * jetX, TiledJet * jetI, vector<TiledJet *> & jets_for_minheap) {
390 double dist = _bj_dist(jetI,jetX);
391 if (dist < jetI->NN_dist) {
393 jetI->NN_dist = dist;
396 if (!jetI->minheap_update_needed()) {
397 jetI->label_minheap_update_needed();
398 jets_for_minheap.push_back(jetI);
402 if (dist < jetX->NN_dist) {
404 jetX->NN_dist = dist;
410 inline void LazyTiling9Alt::_set_NN(TiledJet * jetI,
411 vector<TiledJet *> & jets_for_minheap) {
415 if (!jetI->minheap_update_needed()) {
416 jetI->label_minheap_update_needed();
417 jets_for_minheap.push_back(jetI);}
419 Tile * tile_ptr = &_tiles[jetI->tile_index];
421 for (Tile::TileFnPair * near_tile = tile_ptr->begin_tiles;
422 near_tile != tile_ptr->end_tiles; near_tile++) {
425 if (jetI->NN_dist < (tile_ptr->*(near_tile->second))(jetI))
continue;
427 for (TiledJet * jetJ = (near_tile->first)->head;
428 jetJ != NULL; jetJ = jetJ->next) {
429 double dist = _bj_dist(jetI,jetJ);
430 if (dist < jetI->NN_dist && jetJ != jetI) {
431 jetI->NN_dist = dist; jetI->NN = jetJ;
456 void LazyTiling9Alt::run() {
460 int n = _jets.size();
461 TiledJet * briefjets =
new TiledJet[n];
462 TiledJet * jetA = briefjets, * jetB;
468 vector<int> tile_union(3*n_tile_neighbours);
471 for (
int i = 0; i< n; i++) {
472 _tj_set_jetinfo(jetA, i);
476 TiledJet * head = briefjets;
479 vector<Tile>::iterator tile;
480 for (tile = _tiles.begin(); tile != _tiles.end(); tile++) {
482 for (jetA = tile->head; jetA != NULL; jetA = jetA->next) {
483 for (jetB = tile->head; jetB != jetA; jetB = jetB->next) {
484 double dist = _bj_dist_not_periodic(jetA,jetB);
485 if (dist < jetA->NN_dist) {jetA->NN_dist = dist; jetA->NN = jetB;}
486 if (dist < jetB->NN_dist) {jetB->NN_dist = dist; jetB->NN = jetA;}
489 for (jetA = tile->head; jetA != NULL; jetA = jetA->next) {
490 if (jetA->NN_dist > tile->max_NN_dist) tile->max_NN_dist = jetA->NN_dist;
493 for (tile = _tiles.begin(); tile != _tiles.end(); tile++) {
494 if (tile->use_periodic_delta_phi) {
496 for (Tile::TileFnPair * RTileFnPair = tile->RH_tiles;
497 RTileFnPair != tile->end_tiles; RTileFnPair++) {
498 Tile *RTile = RTileFnPair->first;
499 for (jetA = tile->head; jetA != NULL; jetA = jetA->next) {
500 double dist_to_tile = ((*tile).*(RTileFnPair->second))(jetA);
512 bool relevant_for_jetA = dist_to_tile <= jetA->NN_dist;
513 bool relevant_for_RTile = dist_to_tile <= RTile->max_NN_dist;
514 if (relevant_for_jetA || relevant_for_RTile) {
515 for (jetB = RTile->head; jetB != NULL; jetB = jetB->next) {
516 double dist = _bj_dist(jetA,jetB);
517 if (dist < jetA->NN_dist) {jetA->NN_dist = dist; jetA->NN = jetB;}
518 if (dist < jetB->NN_dist) {jetB->NN_dist = dist; jetB->NN = jetA;}
527 for (Tile::TileFnPair* RTileFnPair = tile->RH_tiles;
528 RTileFnPair != tile->end_tiles; RTileFnPair++) {
529 Tile *RTile = RTileFnPair->first;
530 for (jetA = tile->head; jetA != NULL; jetA = jetA->next) {
531 double dist_to_tile = ((*tile).*(RTileFnPair->second))(jetA);
532 bool relevant_for_jetA = dist_to_tile <= jetA->NN_dist;
533 bool relevant_for_RTile = dist_to_tile <= RTile->max_NN_dist;
534 if (relevant_for_jetA || relevant_for_RTile) {
535 for (jetB = RTile->head; jetB != NULL; 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;}
550 for (tile = _tiles.begin(); tile != _tiles.end(); tile++) {
551 tile->max_NN_dist = 0;
552 for (jetA = tile->head; jetA != NULL; jetA = jetA->next) {
553 if (jetA->NN_dist > tile->max_NN_dist) tile->max_NN_dist = jetA->NN_dist;
558 vector<double> diJs(n);
559 for (
int i = 0; i < n; i++) {
560 diJs[i] = _bj_diJ(&briefjets[i]);
561 briefjets[i].label_minheap_update_done();
563 MinHeap minheap(diJs);
565 vector<TiledJet *> jets_for_minheap;
566 jets_for_minheap.reserve(n);
569 int history_location = n-1;
572 double diJ_min = minheap.minval() *_invR2;
573 jetA = head + minheap.minloc();
585 if (jetA < jetB) {std::swap(jetA,jetB);}
588 _cs.plugin_record_ij_recombination(jetA->_jets_index, jetB->_jets_index, diJ_min, nn);
591 _bj_remove_from_tiles(jetA);
593 _bj_remove_from_tiles(jetB);
594 _tj_set_jetinfo(jetB, nn);
599 _cs.plugin_record_iB_recombination(jetA->_jets_index, diJ_min);
600 _bj_remove_from_tiles(jetA);
604 minheap.remove(jetA-head);
610 int n_near_tiles = 0;
611 _add_untagged_neighbours_to_tile_union_using_max_info(jetA,
612 tile_union, n_near_tiles);
614 _add_untagged_neighbours_to_tile_union_using_max_info(&oldB,
615 tile_union,n_near_tiles);
616 jetB->label_minheap_update_needed();
617 jets_for_minheap.push_back(jetB);
626 Tile & jetB_tile = _tiles[jetB->tile_index];
627 for (Tile::TileFnPair * near_tile_fn_pair = jetB_tile.begin_tiles;
628 near_tile_fn_pair != jetB_tile.end_tiles; near_tile_fn_pair++) {
629 Tile * near_tile = near_tile_fn_pair->first;
631 double dist_to_tile = (jetB_tile.*(near_tile_fn_pair->second))(jetB);
634 bool relevant_for_jetB = dist_to_tile <= jetB->NN_dist;
635 bool relevant_for_near_tile = dist_to_tile <= near_tile->max_NN_dist;
636 bool relevant = relevant_for_jetB || relevant_for_near_tile;
641 if (near_tile->tagged) {
642 for (TiledJet * jetI = near_tile->head; jetI != NULL; jetI = jetI->next) {
643 if (jetI->NN == jetA || jetI->NN == jetB) _set_NN(jetI, jets_for_minheap);
644 _update_jetX_jetI_NN(jetB, jetI, jets_for_minheap);
646 near_tile->tagged =
false;
648 for (TiledJet * jetI = near_tile->head; jetI != NULL; jetI = jetI->next) {
649 _update_jetX_jetI_NN(jetB, jetI, jets_for_minheap);
691 for (
int itile = 0; itile < n_near_tiles; itile++) {
692 Tile * tile_ptr = &_tiles[tile_union[itile]];
693 if (!tile_ptr->tagged)
continue;
694 tile_ptr->tagged =
false;
696 for (TiledJet * jetI = tile_ptr->head; jetI != NULL; jetI = jetI->next) {
698 if (jetI->NN == jetA || (jetI->NN == jetB && jetB != NULL)) {
699 _set_NN(jetI, jets_for_minheap);
706 while (jets_for_minheap.size() > 0) {
707 TiledJet * jetI = jets_for_minheap.back();
708 jets_for_minheap.pop_back();
709 minheap.update(jetI-head, _bj_diJ(jetI));
710 jetI->label_minheap_update_done();
713 Tile & tile_I = _tiles[jetI->tile_index];
714 if (tile_I.max_NN_dist < jetI->NN_dist) tile_I.max_NN_dist = jetI->NN_dist;
724 FASTJET_END_NAMESPACE
ostream & operator<<(ostream &, PseudoJet &)
does the actual work for printing out a jet