openCARP
Doxygen code documentation for the open cardiac electrophysiology simulator openCARP
kdpart.hpp
Go to the documentation of this file.
1 // ----------------------------------------------------------------------------
2 // openCARP is an open cardiac electrophysiology simulator.
3 //
4 // Copyright (C) 2020 openCARP project
5 //
6 // This program is licensed under the openCARP Academic Public License (APL)
7 // v1.0: You can use and redistribute it and/or modify it in non-commercial
8 // academic environments under the terms of APL as published by the openCARP
9 // project v1.0, or (at your option) any later version. Commercial use requires
10 // a commercial license (info@opencarp.org).
11 //
12 // This program is distributed without any warranty; see the openCARP APL for
13 // more details.
14 //
15 // You should have received a copy of the openCARP APL along with this program
16 // and can find it online: http://www.opencarp.org/license
17 // ----------------------------------------------------------------------------
18 
27 #ifndef _KDPART_HPP
28 #define _KDPART_HPP
29 
30 #include <algorithm>
31 #include <iostream>
32 #include <math.h>
33 #include <numeric>
34 #include <vector>
35 
36 #ifdef KDPART_MPI
37 #include <mpi.h>
38 #endif
39 
40 namespace kdpart {
41 
43 template<class S>
44 struct vec3 {
45  S x = S(), y = S(), z = S();
46 
47  void get(const S* p) {
48  x = p[0], y = p[1], z = p[2];
49  }
50  void set(S* p) {
51  p[0] = x, p[1] = y, p[2] = z;
52  }
53 };
54 
56 template<class S>
57 struct bbox {
59  kdpart::vec3<S> bounds[2];
60 };
61 
63 enum axis {
64  X = 0, Y, Z, UNSET
65 };
66 
68 template<class T, class S>
69 struct elem {
71  T eidx = -1;
72 };
73 
75 template<class T, class S>
76 struct partition {
77  T pidx = T(-1);
78  T cnt = T(-1);
80  std::vector<elem<T,S> > * elems = NULL;
81 
83  if(elems) delete elems;
84  }
85 };
86 
91 inline bool is_power_of_two(double val)
92 {
93  double log_two = log(val) / log(2.0);
94  int log_two_int = log_two;
95 
96  return (fabs(double(log_two_int) - log_two) < 1e-6);
97 }
98 
99 
100 
107 template<class T, class S>
108 struct mixed_pair {
109  T v1;
110  S v2;
111 };
112 
114 template<class T, class S>
115 bool operator< (const mixed_pair<T,S> & lhs, const mixed_pair<T,S> & rhs)
116 {
117  return lhs.v1 < rhs.v1;
118 }
119 
131 template<typename V, typename W>
132 V clamp(const V val, const W start, const W end) {
133  if(val < start) return start;
134  if(val > end) return end;
135  return val;
136 }
137 
139 template<class T>
140 inline void dsp_from_cnt(const std::vector<T> & cnt, std::vector<T> & dsp)
141 {
142  dsp.resize(cnt.size()+1);
143  dsp[0] = 0;
144  for(size_t i=0; i<cnt.size(); i++) dsp[i+1] = dsp[i] + cnt[i];
145 }
146 
148 template<class T>
149 inline void cnt_from_dsp(const std::vector<T> & dsp, std::vector<T> & cnt)
150 {
151  cnt.resize(dsp.size() - 1);
152  for(size_t i=0; i<dsp.size()-1; i++) cnt[i] = dsp[i+1] - dsp[i];
153 }
154 
155 template<class V, class W>
156 void sort_copy(std::vector<V> & v1, std::vector<W> & v2)
157 {
158  assert(v1.size() == v2.size());
159 
160  std::vector<mixed_pair<V,W> > pair_array(v1.size());
161 
162  for(size_t i=0; i<v1.size(); i++) {
163  pair_array[i].v1 = v1[i];
164  pair_array[i].v2 = v2[i];
165  }
166 
167  std::sort(pair_array.begin(), pair_array.end());
168 
169  for(size_t i=0; i<v1.size(); i++) {
170  v1[i] = pair_array[i].v1;
171  v2[i] = pair_array[i].v2;
172  }
173 }
174 
175 
185 #define KD_ORDER_INC 5.0
186 
187 #define KD_MIN_SIZE 64
188 
189 #ifdef KDPART_MPI
190 template<class T, class S>
192 {
193  private:
194  // the kdtree_partitioner members
195  std::vector<kdpart::partition<T,S> > _layout;
196  T _cur_part_number;
197  MPI_Comm _comm;
198 
200  inline kdpart::bbox<S> get_bbox(const std::vector<elem<T,S> > & elems)
201  {
202  double minmax[6], minmax_red[6];
203 
204  minmax[0] = 1e100;
205  minmax[1] = 1e100;
206  minmax[2] = 1e100;
207  minmax[3] = -1e100;
208  minmax[4] = -1e100;
209  minmax[5] = -1e100;
210 
211  for(size_t i=0; i<elems.size(); i++) {
212  kdpart::vec3<S> p = elems[i].ctr;
213 
214  if(minmax[0] > p.x) minmax[0] = p.x;
215  if(minmax[1] > p.y) minmax[1] = p.y;
216  if(minmax[2] > p.z) minmax[2] = p.z;
217  if(minmax[3] < p.x) minmax[3] = p.x;
218  if(minmax[4] < p.y) minmax[4] = p.y;
219  if(minmax[5] < p.z) minmax[5] = p.z;
220  }
221 
222  MPI_Allreduce(minmax, minmax_red, 3, MPI_DOUBLE, MPI_MIN, _comm);
223  MPI_Allreduce(minmax+3, minmax_red+3, 3, MPI_DOUBLE, MPI_MAX, _comm);
224 
225  kdpart::bbox<S> box;
226  kdpart::vec3<S> & min = box.bounds[0];
227  kdpart::vec3<S> & max = box.bounds[1];
228 
229  min.x = minmax_red[0];
230  min.y = minmax_red[1];
231  min.z = minmax_red[2];
232  max.x = minmax_red[3];
233  max.y = minmax_red[4];
234  max.z = minmax_red[5];
235 
236  return box;
237  }
238 
240  inline kdpart::axis get_longest_axis(const kdpart::bbox<S> & box)
241  {
242  const kdpart::vec3<S> & min = box.bounds[0];
243  const kdpart::vec3<S> & max = box.bounds[1];
244 
245  S x = fabs(min.x - max.x);
246  S y = fabs(min.y - max.y);
247  S z = fabs(min.z - max.z);
248 
249  return x > y && x > z ? X : y > z ? Y : Z;
250  }
251 
259  inline void print_layout(const T split_pos = -1)
260  {
261  int rank; MPI_Comm_rank(_comm, &rank);
262 
263  if(rank != 0) return;
264 
265  if(split_pos > -1) {
266  T idx = 0;
267  while(idx < split_pos) {
268  printf("%d : %d \n", int(_layout[idx].pidx), int(_layout[idx].cnt));
269  idx++;
270  }
271 
272  printf("----\n");
273  printf("%d : %d \n", int(_layout[idx].pidx), int(_layout[idx].cnt));
274  printf("%d : %d \n", int(_layout[idx+1].pidx), int(_layout[idx+1].cnt));
275  printf("----\n");
276  idx += 2;
277 
278  while(size_t(idx) < _layout.size() && _layout[idx].pidx > -1) {
279  printf("%d : %d \n", int(_layout[idx].pidx), int(_layout[idx].cnt));
280  idx++;
281  }
282  printf("\n");
283  }
284  else {
285  T idx = 0;
286  while(size_t(idx) < _layout.size() && _layout[idx].pidx > -1) {
287  printf("%d : %d \n", int(_layout[idx].pidx), int(_layout[idx].cnt));
288  idx++;
289  }
290  printf("\n");
291  }
292  }
293 
294  inline void get_parallel_median_split(std::vector<mixed_pair<S,T>> & vals,
295  double min, double max,
296  std::vector<bool> & on_left_side)
297  {
298  int size, rank;
299  MPI_Comm_size(_comm, &size);
300  MPI_Comm_rank(_comm, &rank);
301 
302  on_left_side.resize(vals.size());
303  std::sort(vals.begin(), vals.end());
304 
305  size_t nelem = vals.size();
306  double bucket_size = (max*1.05 - min) / double(size);
307  std::vector<int> buckets(size_t(size), int(0));
308  for(size_t i=0; i<nelem; i++) {
309  int idx = (vals[i].v1 - min) / bucket_size;
310  idx = kdpart::clamp(idx, 0, size-1);
311  buckets[idx]++;
312  }
313 
314  std::vector<int> glob_buckets(size_t(size), int(0));
315  MPI_Allreduce(buckets.data(), glob_buckets.data(), size, MPI_INT, MPI_SUM, _comm);
316 
317  // compute the global number of values and consequently the half number
318  int glob_sum = std::accumulate(glob_buckets.begin(), glob_buckets.end(), 0);
319  int lhalf = (glob_sum + 1) / 2;
320  // figure out the bucket index we have to go through to get the median value
321  int dsp = 0, bucket_idx = 0;
322  while(bucket_idx < size && (dsp + glob_buckets[bucket_idx]) <= lhalf) {
323  dsp += glob_buckets[bucket_idx];
324  bucket_idx++;
325  }
326 
327  // sanity check. we could clamp here, but its better to throw an error since
328  // an illegal index should not occur.
329  if(bucket_idx < 0 || bucket_idx >= size) {
330  fprintf(stderr, "Error: Illegal bucket index !!\n");
331  }
332 
333  // containers for the values we compute the median on
334  std::vector<double> loc_val_bucket, glob_val_bucket;
335  // containers for the split decision
336  std::vector<short> global_split_bucket, split_bucket(buckets[bucket_idx]);
337 
338  loc_val_bucket.resize(buckets[bucket_idx]);
339  for(size_t i=0, widx=0; i<nelem; i++) {
340  int idx = (vals[i].v1 - min) / bucket_size;
341  idx = kdpart::clamp(idx, 0, size-1);
342 
343  if(idx == bucket_idx)
344  loc_val_bucket[widx++] = vals[i].v1;
345  }
346 
347  std::vector<int> rcnt(size), rdsp(size);
348  MPI_Gather(&buckets[bucket_idx], 1, MPI_INT, rcnt.data(), 1, MPI_INT, bucket_idx, _comm);
349  kdpart::dsp_from_cnt(rcnt, rdsp);
350 
351  if(rank == bucket_idx) {
352  glob_val_bucket.resize(glob_buckets[bucket_idx]);
353  global_split_bucket.resize(glob_buckets[bucket_idx]);
354  }
355 
356  // the elements of the split bucket are communicated to the associated rank
357  MPI_Gatherv(loc_val_bucket.data(), buckets[bucket_idx], MPI_DOUBLE,
358  glob_val_bucket.data(), rcnt.data(), rdsp.data(), MPI_DOUBLE,
359  bucket_idx, _comm);
360 
361 
362  // the rank owning the bucket where the median split will take place is categorizing
363  // his elements into left and right
364  if(rank == bucket_idx) {
365  size_t bsize = glob_val_bucket.size();
366 
367  std::vector<int> bucket_perm(bsize);
368  for(size_t i=0; i<bsize; i++) bucket_perm[i] = i;
369 
370  sort_copy(glob_val_bucket, bucket_perm);
371  int loc_half_idx = lhalf - dsp;
372 
373  // sanity check
374  if(loc_half_idx < 0 || loc_half_idx >= int(glob_val_bucket.size())) {
375  fprintf(stderr, "Error: Illegal local val index!!\n");
376  }
377 
378  for(int i=0; i <= loc_half_idx; i++)
379  global_split_bucket[bucket_perm[i]] = 1;
380 
381  for(int i=loc_half_idx+1; i < int(bsize); i++)
382  global_split_bucket[bucket_perm[i]] = 0;
383  }
384 
385  MPI_Scatterv(global_split_bucket.data(), rcnt.data(), rdsp.data(), MPI_SHORT,
386  split_bucket.data(), buckets[bucket_idx], MPI_SHORT,
387  bucket_idx, _comm);
388 
389  for(size_t i=0, ridx=0; i<nelem; i++) {
390  int pidx = vals[i].v2;
391  int idx = (vals[i].v1 - min) / bucket_size;
392  idx = kdpart::clamp(idx, 0, size-1);
393 
394  if(idx < bucket_idx)
395  on_left_side[pidx] = true;
396  else if(idx == bucket_idx)
397  on_left_side[pidx] = split_bucket[ridx++] == 1;
398  else
399  on_left_side[pidx] = false;
400  }
401  }
402 
410  inline void median_split(kdpart::partition<T,S> parent,
411  kdpart::partition<T,S> & lchild,
412  kdpart::partition<T,S> & rchild)
413  {
414  kdpart::axis longest_axis = get_longest_axis(parent.box);
415  size_t nelem = parent.elems->size();
416 
417  std::vector<mixed_pair<S,T>> vals(nelem);
418 
419  double min = 0, max = 0;
420  switch(longest_axis) {
421  case X: {
422  min = parent.box.bounds[0].x, max = parent.box.bounds[1].x;
423  for(size_t i=0; i<nelem; i++) {
424  const kdpart::vec3<S> & p = (*parent.elems)[i].ctr;
425  vals[i].v1 = p.x;
426  vals[i].v2 = i;
427  }
428  break;
429  }
430  case Y: {
431  min = parent.box.bounds[0].y, max = parent.box.bounds[1].y;
432  for(size_t i=0; i<nelem; i++) {
433  const kdpart::vec3<S> & p = (*parent.elems)[i].ctr;
434  vals[i].v1 = p.y;
435  vals[i].v2 = i;
436  }
437  break;
438  }
439  case Z: {
440  min = parent.box.bounds[0].z, max = parent.box.bounds[1].z;
441  for(size_t i=0; i<nelem; i++) {
442  const kdpart::vec3<S> & p = (*parent.elems)[i].ctr;
443  vals[i].v1 = p.z;
444  vals[i].v2 = i;
445  }
446  break;
447  }
448  case UNSET: break;
449  }
450 
451  // we compute the median value in parallel
452  std::vector<bool> on_left;
453  get_parallel_median_split(vals, min, max, on_left);
454 
455  // then we split the local elements based on the median
456  size_t left_size = 0, right_size = 0;
457 
458  for(size_t i=0; i<nelem; i++) {
459  if(on_left[i]) left_size++;
460  else right_size++;
461  }
462 
463  lchild.elems = new std::vector<elem<T,S> >(left_size), rchild.elems = new std::vector<elem<T,S> >(right_size);
464 
465  left_size = 0, right_size = 0;
466  for(size_t i=0; i<nelem; i++) {
467  if(on_left[i]) (*lchild.elems)[left_size++ ] = (*parent.elems)[i];
468  else (*rchild.elems)[right_size++] = (*parent.elems)[i];
469  }
470 
471  // remove elements of parent partition
472  delete parent.elems;
473  parent.elems = NULL;
474 
475  int buff[2] = {int(lchild.elems->size()), int(rchild.elems->size())};
476  MPI_Allreduce(MPI_IN_PLACE, buff, 2, MPI_INT, MPI_SUM, _comm);
477  lchild.cnt = buff[0], rchild.cnt = buff[1];
478 
479  if(lchild.cnt == 0 || rchild.cnt == 0) {
480  fprintf(stderr, "Error: Empty partitioning!!\n");
481  }
482 
483  lchild.box = get_bbox(*lchild.elems);
484  rchild.box = get_bbox(*rchild.elems);
485  }
486 
492  inline void update_layout(const T split_pos) {
493  assert(size_t(split_pos) < _layout.size());
494 
495  T idx_at_split = _layout[split_pos].pidx;
496  T start = _cur_part_number - 1, stop = split_pos + 1;
497 
498  // increment partition index for partitions after the split
499  for(T i = start; i > stop; i--) {
500  _layout[i] = _layout[i-1]; // copy partition
501  _layout[i].pidx++; // increment partition index
502  }
503 
504  // set partition indices for the split
505  median_split(_layout[split_pos], _layout[split_pos], _layout[split_pos+1]);
506  _layout[split_pos].pidx = idx_at_split;
507  _layout[split_pos+1].pidx = idx_at_split+1;
508 
509  // print layout
510  // print_layout(split_pos);
511  }
512 
518  inline T get_split_pos()
519  {
520  T idx = 0, max = _layout[0].cnt, maxidx = 0;
521 
522  while(idx < _cur_part_number && _layout[idx].cnt > -1) {
523  if(max < _layout[idx].cnt) {
524  max = _layout[idx].cnt;
525  maxidx = idx;
526  }
527  idx++;
528  }
529 
530  return maxidx;
531  }
532 
533  public:
534  inline void operator() (const MPI_Comm comm, const std::vector<S> & ctr, const int req_part,
535  std::vector<T> & part_vec)
536  {
537  _comm = comm;
538 
539  assert(ctr.size() % 3 == 0); // coord components need to be a multiple of 3
540 
541  // get domain sizes
542  long int l_numelem = ctr.size() / 3, g_numelem;
543  MPI_Allreduce(&l_numelem, &g_numelem, 1, MPI_LONG, MPI_SUM, _comm);
544  assert(g_numelem > 0);
545 
546  int size, rank;
547  MPI_Comm_size(_comm, &size); MPI_Comm_rank(_comm, &rank);
548 
549  T npart = req_part;
550 
551  bool redistribute_remainder = false;
552  T npart_old = npart;
553 
561  if(!kdpart::is_power_of_two(npart)) {
562  npart = (log(double(npart)) / log(2.0)) + KD_ORDER_INC;
563  npart = pow(2, npart);
564 
565  // we can always afford to split at least into KD_MIN_SIZE parts. if the initial npart was very small,
566  // the computed new npart is still below KD_MIN_SIZE. so we set it explicitly
567  if(npart < KD_MIN_SIZE && g_numelem > KD_MIN_SIZE) npart = KD_MIN_SIZE;
568 
569  redistribute_remainder = true;
570  }
571 
572  assert(g_numelem > (long int)npart);
573 
574  // initialize the first partition
575  _cur_part_number = 1;
576  _layout.resize(size_t(npart));
577  _layout[0].pidx = 0;
578  _layout[0].elems = new std::vector<kdpart::elem<T,S> >(l_numelem);
579  _layout[0].cnt = g_numelem;
580 
581  std::vector<kdpart::elem<T,S> > & elems = *_layout[0].elems;
582 
583  // convert the elements into an array of center-points
584  for(long int i=0; i<l_numelem; i++) {
585  elems[i].ctr.get(ctr.data() + i*3); // get elem center coord
586  elems[i].eidx = i; // get elem index
587  }
588 
589  _layout[0].box = get_bbox(elems);
590 
591  /*
592  * Main splitting loop: We iterate until we have enough partitions. In each
593  * iteration, we split the largest partition into two. Thus we can get to any
594  * number of partitions.
595  *
596  */
597  while(_cur_part_number < npart) {
598  T split_pos = get_split_pos();
599  _cur_part_number++;
600 
601  update_layout(split_pos);
602  }
603 
604  // in case the requested number of partitions was not a power of two,
605  // we have computed more partitions than we have processes and we need to re-
606  // index the computed partitions into [0, size]
607  if(redistribute_remainder) {
608  // number of partitions we will at least assign to a process
609  T base_size = npart / npart_old;
610  // number of partitions we will assign to a subset of processes to reduce
611  // the remainder (npart / npart_old)
612  T extended_size = base_size + 1;
613  T remainder = npart % npart_old;
614 
615  // 'remainder' many processes get chunks of size 'extended_size'
616  for(T i=0; i<remainder; i++) {
617  for(T j=0; j<extended_size; j++)
618  _layout[i*extended_size+j].pidx = i;
619  }
620 
621  // '_layout.size() - remainder' many processes get chunks of size 'base_size'
622  for(T i=remainder*extended_size, pidx=remainder; i<T(_layout.size()); i+=base_size, pidx++) {
623  for(T j=0; j<base_size; j++)
624  _layout[i+j].pidx = pidx;
625  }
626 
627  // print_layout();
628  }
629 
630  // assign partition index to the individual elements
631  part_vec.assign(l_numelem, T(-1));
632  for(const kdpart::partition<T,S> & p : _layout) {
633  for(const kdpart::elem<T,S> & e : (*p.elems)) {
634  part_vec[e.eidx] = p.pidx;
635  }
636  }
637  }
638 };
639 #endif
640 
641 template<class T, class S>
643 
644  private:
645  // the kdpart members
646  std::vector<kdpart::partition<T,S> > _layout;
647  T _cur_part_number;
648 
650  inline kdpart::bbox<S> get_bbox(const std::vector<elem<T,S> > & elems)
651  {
652  kdpart::bbox<S> box;
653  kdpart::vec3<S> & min = box.bounds[0];
654  kdpart::vec3<S> & max = box.bounds[1];
655  kdpart::vec3<S> p = elems[0].ctr;
656 
657  min.x = p.x, min.y = p.y, min.z = p.z;
658  max.x = p.x, max.y = p.y, max.z = p.z;
659 
660  for(size_t i=1; i<elems.size(); i++) {
661  p = elems[i].ctr;
662 
663  if(min.x > p.x) min.x = p.x;
664  if(min.y > p.y) min.y = p.y;
665  if(min.z > p.z) min.z = p.z;
666  if(max.x < p.x) max.x = p.x;
667  if(max.y < p.y) max.y = p.y;
668  if(max.z < p.z) max.z = p.z;
669  }
670 
671  return box;
672  }
673 
675  inline kdpart::axis get_longest_axis(const kdpart::bbox<S> & box)
676  {
677  const kdpart::vec3<S> & min = box.bounds[0];
678  const kdpart::vec3<S> & max = box.bounds[1];
679 
680  S x = fabs(min.x - max.x);
681  S y = fabs(min.y - max.y);
682  S z = fabs(min.z - max.z);
683 
684  return x > y && x > z ? X : y > z ? Y : Z;
685  }
686 
694  inline void print_layout(const T split_pos = -1)
695  {
696  if(split_pos > -1) {
697  T idx = 0;
698  while(idx < split_pos) {
699  printf("%d : %d \n", int(_layout[idx].pidx), int(_layout[idx].cnt));
700  idx++;
701  }
702 
703  printf("----\n");
704  printf("%d : %d \n", int(_layout[idx].pidx), int(_layout[idx].cnt));
705  printf("%d : %d \n", int(_layout[idx+1].pidx), int(_layout[idx+1].cnt));
706  printf("----\n");
707  idx += 2;
708 
709  while(size_t(idx) < _layout.size() && _layout[idx].pidx > -1) {
710  printf("%d : %d \n", int(_layout[idx].pidx), int(_layout[idx].cnt));
711  idx++;
712  }
713  printf("\n");
714  }
715  else {
716  T idx = 0;
717  while(size_t(idx) < _layout.size() && _layout[idx].pidx > -1) {
718  printf("%d : %d \n", int(_layout[idx].pidx), int(_layout[idx].cnt));
719  idx++;
720  }
721  printf("\n");
722  }
723  }
724 
732  inline void median_split(kdpart::partition<T,S> parent,
734  {
735 
736  kdpart::axis longest_axis = get_longest_axis(parent.box);
737  size_t nelem = parent.elems->size();
738 
739  // put the (coord, idx) pairs into a vector and sort them
740  std::vector<mixed_pair<S,T> > pairs(nelem);
741  for(size_t i=0; i<nelem; i++) {
742  const kdpart::vec3<S> & p = (*parent.elems)[i].ctr;
743  S val = S();
744  switch(longest_axis) {
745  case X: val = p.x; break;
746  case Y: val = p.y; break;
747  case Z: val = p.z; break;
748  case UNSET: break;
749  }
750 
751  pairs[i] = {val, T(i)};
752  }
753  std::sort(pairs.begin(), pairs.end());
754 
755  // we then copy the first half into left children and the other half into
756  // right children
757  size_t lnum = (nelem + 1) / 2, rnum = nelem - lnum;
758  lchild.elems = new std::vector<elem<T,S> >(lnum), rchild.elems = new std::vector<elem<T,S> >(rnum);
759 
760  for(size_t i=0; i<lnum; i++) {
761  T lidx = pairs[i].v2;
762  (*lchild.elems)[i] = (*parent.elems)[lidx];
763  }
764 
765  for(size_t i=0; i<rnum; i++) {
766  T lidx = pairs[lnum + i].v2;
767  (*rchild.elems)[i] = (*parent.elems)[lidx];
768  }
769 
770  // remove elements of parent partition
771  delete parent.elems;
772  parent.elems = NULL;
773 
774  lchild.cnt = lchild.elems->size();
775  rchild.cnt = rchild.elems->size();
776  lchild.box = get_bbox(*lchild.elems);
777  rchild.box = get_bbox(*rchild.elems);
778  }
779 
785  inline void update_layout(const T split_pos) {
786  assert(size_t(split_pos) < _layout.size());
787 
788  T idx_at_split = _layout[split_pos].pidx;
789  T start = _cur_part_number - 1, stop = split_pos + 1;
790 
791  // increment partition index for partitions after the split
792  for(T i = start; i > stop; i--) {
793  _layout[i] = _layout[i-1]; // copy partition
794  _layout[i].pidx++; // increment partition index
795  }
796 
797  // set partition indices for the split
798  median_split(_layout[split_pos], _layout[split_pos], _layout[split_pos+1]);
799  _layout[split_pos].pidx = idx_at_split;
800  _layout[split_pos+1].pidx = idx_at_split+1;
801 
802  // print layout
803  // print_layout(split_pos);
804  }
805 
811  inline T get_split_pos()
812  {
813  T idx = 0, max = _layout[0].cnt, maxidx = 0;
814 
815  while(idx < _cur_part_number && _layout[idx].cnt > -1) {
816  if(max < _layout[idx].cnt) {
817  max = _layout[idx].cnt;
818  maxidx = idx;
819  }
820  idx++;
821  }
822 
823  return maxidx;
824  }
825 
826  public:
827  inline void operator() (const std::vector<S> & ctr, T npart, std::vector<T> & part_vec)
828  {
829  assert(ctr.size() % 3 == 0);
830 
831  size_t nelem = ctr.size() / 3;
832  assert(nelem > size_t(npart));
833 
834  bool redistribute_remainder = false;
835  T npart_old = npart;
836 
844  if(!is_power_of_two(npart)) {
845  npart = (log(double(npart)) / log(2.0)) + 5.0;
846  npart = pow(2, npart);
847 
848  // we can always afford to split at least into KD_MIN_SIZE parts. if the initial npart was very small,
849  // the computed new npart is still below KD_MIN_SIZE. so we set it explicitly
850  if(npart < KD_MIN_SIZE && nelem > KD_MIN_SIZE) npart = KD_MIN_SIZE;
851 
852  redistribute_remainder = true;
853  }
854 
855  // initialize _layout and _cnt
856  _layout.resize(size_t(npart));
857  part_vec.assign(nelem, T(-1));
858 
859  _cur_part_number = 1;
860  _layout[0].pidx = 0;
861  _layout[0].elems = new std::vector<kdpart::elem<T,S> >(nelem);
862  std::vector<kdpart::elem<T,S> > & elems = *_layout[0].elems;
863 
864  for(size_t i=0; i<nelem; i++) {
865  elems[i].ctr.get(ctr.data() + i*3); // get elem center coord
866  elems[i].eidx = i; // get elem index
867  }
868 
869  _layout[0].box = get_bbox(elems);
870 
871  // splitting loop
872  while(_cur_part_number < npart) {
873  T split_pos = get_split_pos();
874  _cur_part_number++;
875 
876  update_layout(split_pos);
877  }
878 
879  if(redistribute_remainder) {
880  T base_size = npart / npart_old, extended_size = base_size + 1;
881  T remainder = npart % npart_old;
882 
883  for(T i=0; i<remainder; i++) {
884  for(T j=0; j<extended_size; j++)
885  _layout[i*extended_size+j].pidx = i;
886  }
887 
888  // the rest gets reindexed ascendingly, starting with index "remainder"
889  for(T i=remainder*extended_size, pidx=remainder; i<T(_layout.size()); i+=base_size, pidx++) {
890  for(T j=0; j<base_size; j++)
891  _layout[i+j].pidx = pidx;
892  }
893 
894  // print_layout();
895  }
896 
897  // assign partition index to the individual elements
898  for(const kdpart::partition<T,S> & p : _layout) {
899  for(const kdpart::elem<T,S> & e : (*p.elems)) {
900  part_vec[e.eidx] = p.pidx;
901  }
902  }
903  }
904 };
905 
906 }
907 
908 #endif
constexpr T min(T a, T b)
Definition: ion_type.h:33
#define KD_ORDER_INC
The value we add to the power of two if we need higher partitioning resolution.
Definition: kdpart.hpp:185
constexpr T max(T a, T b)
Definition: ion_type.h:31
Combined floating point and integer pair.
Definition: kdpart.hpp:108
minimalistic internal point struct
Definition: kdpart.hpp:44
void cnt_from_dsp(const std::vector< T > &dsp, std::vector< T > &cnt)
Compute counts from displacements.
Definition: kdpart.hpp:149
element definition
Definition: kdpart.hpp:69
std::vector< elem< T, S > > * elems
(local) elements in partition
Definition: kdpart.hpp:80
kdpart::bbox< S > box
(global) partition bounding box
Definition: kdpart.hpp:79
kdpart::vec3< S > bounds[2]
the bounds. bounds[0] = lower left (min), bounds[1] = upper right (max)
Definition: kdpart.hpp:59
bool is_power_of_two(double val)
Check if a given value is an integer power of two.
Definition: kdpart.hpp:91
T cnt
(global) partition size
Definition: kdpart.hpp:78
the struct holding all partition data
Definition: kdpart.hpp:76
void sort_copy(std::vector< V > &v1, std::vector< W > &v2)
Definition: kdpart.hpp:156
void dsp_from_cnt(const std::vector< T > &cnt, std::vector< T > &dsp)
Compute displacements from counts.
Definition: kdpart.hpp:140
Bounding box struct.
Definition: kdpart.hpp:57
vec3< S > ctr
element center location
Definition: kdpart.hpp:70
#define KD_MIN_SIZE
Definition: kdpart.hpp:187
axis
split axis
Definition: kdpart.hpp:63
V clamp(const V val, const W start, const W end)
Clamp a value into an interval [start, end].
Definition: kdpart.hpp:132