YAC 3.9.1
Yet Another Coupler
All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
sphere_part.c
Go to the documentation of this file.
1// Copyright (c) 2024 The YAC Authors
2//
3// SPDX-License-Identifier: BSD-3-Clause
4
5#ifdef HAVE_CONFIG_H
6// Get the definition of the 'restrict' keyword.
7#include "config.h"
8#endif
9
10#include <assert.h>
11#include <stdlib.h>
12#include <math.h>
13#include <string.h>
14#include <stdint.h>
15#include <float.h>
16
17#include "sphere_part.h"
18#include "geometry.h"
19#include "interval_tree.h"
20#include "utils_core.h"
21#include "ensure_array_size.h"
22
23union I_list {
24 struct {
26 size_t num_nodes;
27 } ivt;
28 size_t *list;
29};
30
31enum {
34};
35enum {
36 U_FLAG = 1,
37 T_FLAG = 2,
38};
39
45
47
48 int flags;
49
50 union I_list I;
51 void * U, * T;
52
54
56
57 double gc_norm_vector[3];
58};
59
64
67 I_NODE = 1,
68 U_NODE = 2,
69 T_NODE = 3,
70};
71
73 struct bounding_circle bnd_circle; // has to be first entry
74 size_t local_id;
76};
77
79 double coordinates_xyz[3]; // has to be the first entry
80 size_t idx; // index of the point in the coordinates array passed to
81 // the constructor
82};
83
88
90
91 int flags;
92
93 void * U, * T;
94
95 size_t U_size, T_size;
96
97 double gc_norm_vector[3];
98};
99
106
108 size_t idx;
109 int32_t coordinate_xyz[3];
110};
111
112static void init_sphere_part_node(struct sphere_part_node * node) {
113
114 node->flags = 0;
115 node->I_size = 0;
116 node->U_size = 0;
117 node->T_size = 0;
118 node->I_angle = SIN_COS_ZERO;
119 node->gc_norm_vector[0] = 0;
120 node->gc_norm_vector[1] = 0;
121 node->gc_norm_vector[2] = 1;
122}
123
125
126 struct sphere_part_node * node = xmalloc(1 * sizeof(*node));
127
129
130 return node;
131}
132
133// computes the norm vector for the plane the partitions the data in half
134// (more or less)
135static inline void compute_gc_norm_vector(
136 void * coords_data, size_t coords_size, size_t coords_count,
137 double prev_gc_norm_vector[], double gc_norm_vector[]) {
138
139
140 double balance_point[3] = {0.0,0.0,0.0};
141
142 // compute balance point
143 for (size_t i = 0; i < coords_count; ++i) {
144
145 double * coords =
146 (double*)(((unsigned char*)coords_data) + i * coords_size);
147 balance_point[0] += coords[0];
148 balance_point[1] += coords[1];
149 balance_point[2] += coords[2];
150 }
151
152 if ((fabs(balance_point[0]) > 1e-9) ||
153 (fabs(balance_point[1]) > 1e-9) ||
154 (fabs(balance_point[2]) > 1e-9)) {
155 normalise_vector(balance_point);
156 } else {
157 balance_point[0] = prev_gc_norm_vector[2];
158 balance_point[1] = prev_gc_norm_vector[0];
159 balance_point[2] = prev_gc_norm_vector[1];
160 }
161
163 balance_point, prev_gc_norm_vector, gc_norm_vector);
164
165 if ((fabs(gc_norm_vector[0]) > 1e-9) ||
166 (fabs(gc_norm_vector[1]) > 1e-9) ||
167 (fabs(gc_norm_vector[2]) > 1e-9)) {
169 } else {
170 gc_norm_vector[0] = prev_gc_norm_vector[2];
171 gc_norm_vector[1] = prev_gc_norm_vector[0];
172 gc_norm_vector[2] = prev_gc_norm_vector[1];
173 }
174}
175
177{
178 struct temp_partition_data temp = *a;
179 *a = *b;
180 *b = temp;
181}
182
183static size_t swap_node_type(struct temp_partition_data *part_data, size_t i, int node_type, size_t begin, size_t end)
184{
185 for (size_t j = begin; j < end; ++j)
186 {
187 if (part_data[j].node_type == node_type)
188 {
189 swap_partition_data(&part_data[i], &part_data[j]);
190 return j + 1;
191 }
192 }
193
194 return end + 1;
195}
196
197static void sort_partition_data(struct temp_partition_data *part_data, size_t I_FULL_size, size_t I_size, size_t U_size, size_t T_size)
198{
199 size_t I_begin = I_FULL_size;
200 size_t U_begin = I_FULL_size + I_size;
201 size_t T_begin = I_FULL_size + I_size + U_size;
202
203 size_t I_end = I_begin + I_size;
204 size_t U_end = U_begin + U_size;
205 size_t T_end = T_begin + T_size;
206
207 while (I_end > I_begin && part_data[I_end-1].node_type == I_NODE) I_end--;
208 while (U_end > U_begin && part_data[U_end-1].node_type == U_NODE) U_end--;
209 while (T_end > T_begin && part_data[T_end-1].node_type == T_NODE) T_end--;
210
211 for (size_t i = 0; i < I_FULL_size; ++i)
212 {
213 if (part_data[i].node_type != I_NODE_FULL)
214 {
215 I_begin = swap_node_type(part_data, i, I_NODE_FULL, I_begin, T_end);
216 }
217 }
218
219 I_begin = I_FULL_size;
220 for (size_t i = I_begin; i < I_end; ++i)
221 {
222 if (part_data[i].node_type != I_NODE)
223 {
224 U_begin = swap_node_type(part_data, i, I_NODE, U_begin, U_end);
225 if (U_begin > U_end)
226 {
227 T_begin = swap_node_type(part_data, i, I_NODE, T_begin, T_end);
228 }
229 }
230 }
231
232 U_begin = I_FULL_size + I_size;
233 T_begin = I_FULL_size + I_size + U_size;
234 for (size_t i = U_begin; i < U_end; ++i)
235 {
236 if (part_data[i].node_type != U_NODE)
237 {
238 T_begin = swap_node_type(part_data, i, U_NODE, T_begin, T_end);
239 }
240 }
241}
242
243static void partition_data (
244 size_t * local_cell_ids, struct temp_partition_data * part_data,
245 size_t num_cell_ids, size_t threshold, struct sphere_part_node * parent_node,
246 double prev_gc_norm_vector[]) {
247
248 if (num_cell_ids == 0) {
249 parent_node->flags = U_IS_LEAF + T_IS_LEAF;
250 return;
251 }
252
254 part_data, sizeof(*part_data), num_cell_ids,
255 prev_gc_norm_vector, parent_node->gc_norm_vector);
256
257 // partition data into cells that overlap with the great circle and cells
258 // that are one side of the circle
259
260 size_t I_FULL_size = 0;
261 size_t I_size = 0;
262 size_t U_size = 0;
263 size_t T_size = 0;
264
265 struct sin_cos_angle max_inc_angle = SIN_COS_ZERO;
266
267 for (size_t i = 0; i < num_cell_ids; ++i) {
268
269 struct bounding_circle curr_bnd_circle = part_data[i].bnd_circle;
270
271 // get angle between the norm vector of the great circle and the base
272 // point of the bounding circle
273 struct sin_cos_angle angle =
275 curr_bnd_circle.base_vector, parent_node->gc_norm_vector);
276
277 // get the angle between between the plane of the great circle and base
278 // point of the bounding circle
279 struct sin_cos_angle diff_angle_gc =
280 sin_cos_angle_new(fabs(angle.cos), angle.sin);
281
282 // if the point intersects with the great circle
283 if (compare_angles(diff_angle_gc, curr_bnd_circle.inc_angle) <= 0) {
284
285 // if gc_norm_vector or -gc_norm_vector is in the bounding circle
286 if ((angle.sin < curr_bnd_circle.inc_angle.sin) ||
287 (0.0 >= curr_bnd_circle.inc_angle.cos)) {
288
289 // set node type for current cell
290 part_data[i].node_type = I_NODE_FULL;
291 I_FULL_size++;
292
293 max_inc_angle = SIN_COS_M_PI_2;
294
295 } else {
296
297 // set node type for current cell
298 part_data[i].node_type = I_NODE;
299 I_size++;
300
301 struct sin_cos_angle inc_angle =
302 sum_angles_no_check(diff_angle_gc, curr_bnd_circle.inc_angle);
303
304 if (compare_angles(max_inc_angle, inc_angle) < 0)
305 max_inc_angle = inc_angle;
306 }
307
308 // angle > M_PI_2
309 } else if (angle.cos < 0.0) {
310
311 // set node type for current cell
312 part_data[i].node_type = U_NODE;
313 U_size++;
314
315 } else {
316
317 // set node type for current cell
318 part_data[i].node_type = T_NODE;
319 T_size++;
320 }
321 }
322
323 sort_partition_data(part_data, I_FULL_size, I_size, U_size, T_size);
324
325 I_size += I_FULL_size;
326 parent_node->I_size = I_size;
327 parent_node->U_size = U_size;
328 parent_node->T_size = T_size;
329
330 // if max_inc_angle > PI/2
331 if (compare_angles(max_inc_angle, SIN_COS_M_PI_2) >= 0) {
332 parent_node->I_angle = SIN_COS_M_PI_2;
333 } else {
334 parent_node->I_angle = max_inc_angle;
335 }
336
337 if (I_size > 0) {
338
339 if (I_size > I_list_tree_min_size) {
340
341 assert(sizeof(struct interval_node) > sizeof(size_t));
342 struct interval_node * head_node = xmalloc(I_size * sizeof(*head_node));
343 parent_node->I.ivt.head_node = head_node;
344 parent_node->I.ivt.num_nodes = I_size;
345
346 // for all bounding circles that include gc_norm_vector or
347 // -gc_norm_vector
348 for (size_t i = 0; i < I_FULL_size; ++i) {
349
350 head_node[i].range.left = -M_PI;
351 head_node[i].range.right = M_PI;
352 head_node[i].value = part_data[i].local_id;
353 }
354
355 for (size_t i = I_FULL_size; i < I_size; ++i) {
356
357 double GCp[3], bVp[3];
358 struct bounding_circle curr_bnd_circle = part_data[i].bnd_circle;
359
360 // project the base vector of the current bounding circle onto the
361 // current partitioning great circle
363 curr_bnd_circle.base_vector, GCp);
364 crossproduct_kahan(GCp, parent_node->gc_norm_vector, bVp);
366 (fabs(bVp[0]) > 1e-9) || (fabs(bVp[1]) > 1e-9) ||
367 (fabs(bVp[2]) > 1e-9),
368 "ERROR(partition_data): "
369 "projected vector is nearly identical to gc_norm_vector")
370 normalise_vector(bVp);
371
372 // the inc_angle of the bounding circle also needs to be projected
373 // onto the great circle
374 // (the formula for this is:
375 // new_inc_angle=acos(sin(inc_angle)/cos(base_angle))
376 // To find this formula, you have to determine the longitude of
377 // of longitude circle that has exactly one intersection with the
378 // bounding circle.
379 // (see https://www.dropbox.com/s/82bxzno1mjogbtf/function_2nd_try.pdf)
380 double bnd_circle_lat_cos =
381 get_vector_angle_2(bVp, curr_bnd_circle.base_vector).cos;
382 double inc_angle =
383 M_PI_2 - acos(curr_bnd_circle.inc_angle.sin/bnd_circle_lat_cos);
384
385 // By definition the previous norm vector is on the current great
386 // great circle. We use this as a reference point for I
387
388 // compute "distance" of the projected base vector to the reference
389 // point
390 double base_angle = get_vector_angle(bVp, prev_gc_norm_vector);
391
392 head_node[i].range.left = base_angle - inc_angle;
393 head_node[i].range.right = base_angle + inc_angle;
394 head_node[i].value = part_data[i].local_id;
395 }
396
397 yac_generate_interval_tree(head_node, I_size);
398 parent_node->flags |= I_IS_INTERVAL_TREE;
399 } else {
400 for (size_t i = 0; i < I_size; ++i)
401 local_cell_ids[i] = part_data[i].local_id;
402 parent_node->I.list = (void*)local_cell_ids;
403 }
404 } else
405 parent_node->I.list = NULL;
406
407 part_data += I_size;
408 local_cell_ids += I_size;
409
410 // check whether the lists are small enough (if not -> partition again)
411 if (U_size <= threshold) {
412
413 for (size_t i = 0; i < U_size; ++i)
414 local_cell_ids[i] = part_data[i].local_id;
415 parent_node->U = (void*)local_cell_ids;
416 parent_node->flags |= U_IS_LEAF;
417
418 } else {
419
420 parent_node->U = get_sphere_part_node();
421 partition_data(local_cell_ids, part_data, U_size, threshold,
422 parent_node->U, parent_node->gc_norm_vector);
423 }
424 local_cell_ids += U_size;
425 part_data += U_size;
426
427 if (T_size <= threshold) {
428
429 for (size_t i = 0; i < T_size; ++i)
430 local_cell_ids[i] = part_data[i].local_id;
431 parent_node->T = (void*)local_cell_ids;
432 parent_node->flags |= T_IS_LEAF;
433 local_cell_ids += T_size;
434
435 } else {
436
437 parent_node->T = get_sphere_part_node();
438 partition_data(local_cell_ids, part_data, T_size, threshold,
439 parent_node->T, parent_node->gc_norm_vector);
440 }
441}
442
443static int compare_point_idx_xyz(void const * a, void const * b) {
444 return (((struct point_id_xyz *)a)->idx > ((struct point_id_xyz *)b)->idx) -
445 (((struct point_id_xyz *)a)->idx < ((struct point_id_xyz *)b)->idx);
446}
447
449 struct point_id_xyz * points, size_t num_points, size_t threshold,
450 double prev_gc_norm_vector[], size_t curr_tree_depth,
451 size_t * max_tree_depth, int * list_flag) {
452
453 if (curr_tree_depth > *max_tree_depth) *max_tree_depth = curr_tree_depth;
454
455 struct point_sphere_part_node * node = xmalloc(1 * sizeof(*node));
456 double * gc_norm_vector = &(node->gc_norm_vector[0]);
457
459 points, sizeof(*points), num_points, prev_gc_norm_vector, gc_norm_vector);
460
461 // angle between a point and the great circle plane
462 // acos(dot(gc_norm_vector, point_xyz)) = angle(gc_norm_vector, point_xyz)
463 // acos(dot(gc_norm_vector, point_xyz)) - PI/2 = angle(gc_plane, point_xyz)
464 // dot <= 0.0 -> U list
465 // dot > 0.0 -> T list
466
468 {
470 for (size_t i = 0; i < num_points; ++i) {
471 double * curr_coordinates_xyz = &(points[i].coordinates_xyz[0]);
472 double dot = curr_coordinates_xyz[0] * gc_norm_vector[0] +
473 curr_coordinates_xyz[1] * gc_norm_vector[1] +
474 curr_coordinates_xyz[2] * gc_norm_vector[2];
475
476 // if (angle >= M_PI_2)
477 list_flag[i] = (dot <= 0.0);
478 }
479 }
480 size_t U_size = 0, T_size = 0;
481 for (size_t i = 0; i < num_points; ++i) {
482 if (list_flag[i]) ++U_size;
483 else ++T_size;
484 }
485
486 // The number of T-points among the first U-size number of points in the
487 // array is equal to the number of U-points in the remaining array. These
488 // have to be exchanged in order to get a sorted vertex array
489
490 // search for all T-points in the U-part of the array and swap them with a
491 // U-point in the T-part
492 for (size_t i = 0, j = U_size; i < U_size; ++i) {
493 // if the current point belongs to the T-list
494 if (!list_flag[i]) {
495 // search for a matching U-point
496 for (;!list_flag[j];++j);
497 struct point_id_xyz temp_point = points[i];
498 points[i] = points[j];
499 points[j] = temp_point;
500 ++j;
501 }
502 }
503
504 node->U_size = U_size;
505 node->T_size = T_size;
506 node->flags = 0;
507
508 // check whether the lists are small enough (if not -> partition again)
509 if ((U_size <= threshold) || (U_size == num_points)) {
510
511 node->U = points;
512 node->flags |= U_IS_LEAF;
513 qsort(points, U_size, sizeof(*points), compare_point_idx_xyz);
514
515 } else {
516
517 node->U =
519 points, U_size, threshold, gc_norm_vector,
520 curr_tree_depth + 1, max_tree_depth, list_flag);
521 }
522
523 if ((T_size <= threshold) || (T_size == num_points)) {
524
525 node->T = points + U_size;
526 node->flags |= T_IS_LEAF;
527 qsort(points + U_size, T_size, sizeof(*points), compare_point_idx_xyz);
528
529 } else {
530
531 node->T =
533 points + U_size, T_size, threshold, gc_norm_vector,
534 curr_tree_depth + 1, max_tree_depth, list_flag);
535 }
536
537 return node;
538}
539
541 struct bounding_circle * circles, size_t num_circles) {
542
543 struct bnd_sphere_part_search * search = xmalloc(1 * sizeof(*search));
544
545 double gc_norm_vector[3] = {0.0,0.0,1.0};
546
548
549 size_t * ids = xmalloc(num_circles * sizeof(*ids));
550 search->ids = ids;
551
552 struct temp_partition_data * part_data =
553 xmalloc(num_circles * sizeof(*part_data));
554 for (size_t i = 0; i < num_circles; ++i) {
555 part_data[i].bnd_circle = circles[i];
556 part_data[i].local_id = i;
557 }
558
559 partition_data(ids, part_data, num_circles, I_list_tree_min_size,
560 &(search->base_node), gc_norm_vector);
561
562 free(part_data);
563
564 return search;
565}
566
568 const void * a,const void * b) {
569
570 struct point_id_xyz_int32 * point_a = (struct point_id_xyz_int32 *)a;
571 struct point_id_xyz_int32 * point_b = (struct point_id_xyz_int32 *)b;
572
573 int ret;
574
575 for (int i = 0; i < 3; ++i)
576 if ((ret = (point_a->coordinate_xyz[i] > point_b->coordinate_xyz[i]) -
577 (point_a->coordinate_xyz[i] < point_b->coordinate_xyz[i])))
578 return ret;
579 return 0;
580}
581
582static struct point_id_xyz *
585 yac_int const * ids, int const * mask) {
586
587 struct point_id_xyz_int32 * points_int32 =
588 xmalloc(*num_points * sizeof(*points_int32));
589
590 double const scale = (double)(2 << 21);
591
592 size_t num_unmasked_points;
593
594 if (mask == NULL) {
595 num_unmasked_points = *num_points;
596 for (size_t i = 0; i < num_unmasked_points; ++i) {
597
598 points_int32[i].idx = i;
599 for (size_t j = 0; j < 3; ++j)
600 points_int32[i].coordinate_xyz[j] =
601 (int32_t)round(coordinates_xyz[i][j] * scale);
602 }
603 } else {
604 num_unmasked_points = 0;
605 for (size_t i = 0; i < *num_points; ++i) {
606
607 if (!mask[i]) continue;
608 points_int32[num_unmasked_points].idx = i;
609 for (size_t j = 0; j < 3; ++j)
610 points_int32[num_unmasked_points].coordinate_xyz[j] =
611 (int32_t)round(coordinates_xyz[i][j] * scale);
612 num_unmasked_points++;
613 }
614 }
615
616 // sort points
617 qsort(points_int32, num_unmasked_points,
618 sizeof(*points_int32), compare_points_int32_coord);
619
620 struct point_id_xyz_int32 dummy;
621 dummy.idx = SIZE_MAX;
622 dummy.coordinate_xyz[0] = INT32_MAX;
623 dummy.coordinate_xyz[1] = INT32_MAX;
624 dummy.coordinate_xyz[2] = INT32_MAX;
625 struct point_id_xyz_int32 * prev = &dummy, * curr = points_int32;
626 yac_int prev_id = XT_INT_MAX;
627 size_t new_num_points = 0;
628 for (size_t i = 0; i < num_unmasked_points; ++i, ++curr) {
629
630 size_t curr_idx = curr->idx;
631 if (compare_points_int32_coord(prev, curr)) {
632 prev = points_int32 + new_num_points++;
633 if (prev != curr) *prev = *curr;
634 prev_id = ids[curr_idx];
635 } else {
636 yac_int curr_id = ids[curr_idx];
637 if (curr_id > prev_id) {
638 prev_id = curr_id;
639 prev->idx = curr_idx;
640 }
641 }
642 }
643
644 struct point_id_xyz * points = xmalloc(new_num_points * sizeof(*points));
645 for (size_t i = 0; i < new_num_points; ++i) {
646 size_t curr_idx = points_int32[i].idx;
647 points[i].idx = curr_idx;
648 memcpy(
649 points[i].coordinates_xyz, coordinates_xyz[curr_idx], 3 * sizeof(double));
650 }
651 *num_points = new_num_points;
652 free(points_int32);
653 return points;
654}
655
657 size_t num_points, yac_const_coordinate_pointer coordinates_xyz,
658 yac_int const * ids) {
659
660 if (num_points == 0) return NULL;
661
662 struct point_sphere_part_search * search = xmalloc(1 * sizeof(*search));
663 struct point_id_xyz * points =
665 search->points = points;
666
667 size_t max_tree_depth = 0;
668
669 int * list_flag = xmalloc(num_points * sizeof(*list_flag));
670
671 // emperical measurements have given a threshold for the leaf size of 2
672 struct point_sphere_part_node * tmp_node =
674 points, num_points, I_list_tree_min_size, (double[3]){0.0,0.0,1.0},
675 1, &max_tree_depth, list_flag);
676
677 free(list_flag);
678
679 search->base_node = *tmp_node;
680 search->max_tree_depth = max_tree_depth;
681 free(tmp_node);
682
683 return search;
684}
685
687 size_t num_points, yac_const_coordinate_pointer coordinates_xyz,
688 yac_int const * ids, int const * mask) {
689
690 if (num_points == 0) return NULL;
691
692 struct point_id_xyz * points =
694
695 struct point_sphere_part_search * search = xmalloc(1 * sizeof(*search));
696 search->points = xrealloc(points, num_points * sizeof(*points));
697
698 size_t max_tree_depth = 0;
699
700 int * list_flag = xmalloc(num_points * sizeof(*list_flag));
701
702 // emperical measurements have given a threshold for the leaf size of 2
703 struct point_sphere_part_node * tmp_node =
706 (double[3]){0.0,0.0,1.0}, 1, &max_tree_depth, list_flag);
707
708 free(list_flag);
709
710 search->base_node = *tmp_node;
711 search->max_tree_depth = max_tree_depth;
712 free(tmp_node);
713
714 return search;
715}
716
718 struct sphere_part_node * node, struct bounding_circle bnd_circle,
719 size_t ** restrict overlap_cells, size_t * overlap_cells_array_size,
720 size_t * restrict num_overlap_cells,
721 struct overlaps * search_interval_tree_buffer, double prev_gc_norm_vector[]) {
722
723 if (node->flags & I_IS_INTERVAL_TREE) {
724
725 struct sin_cos_angle angle =
727 bnd_circle.base_vector, node->gc_norm_vector);
728 angle.cos = fabs(angle.cos);
729
730 // if gc_norm_vector or -gc_norm_vector is in the bounding circle
731 if (compare_angles(angle, bnd_circle.inc_angle) <= 0) {
732 ENSURE_ARRAY_SIZE(*overlap_cells, *overlap_cells_array_size,
733 *num_overlap_cells + node->I.ivt.num_nodes);
734 for (size_t i = 0; i < node->I.ivt.num_nodes; ++i)
735 (*overlap_cells)[(*num_overlap_cells)+i] =
736 node->I.ivt.head_node[i].value;
737 *num_overlap_cells += node->I.ivt.num_nodes;
738 return;
739 }
740
741 double GCp[3], bVp[3];
742
743 // project the base vector of the current bounding circle onto the
744 // current great circle
745 crossproduct_kahan(node->gc_norm_vector, bnd_circle.base_vector, GCp);
746 crossproduct_kahan(GCp, node->gc_norm_vector, bVp);
748 (fabs(bVp[0]) > 1e-9) || (fabs(bVp[1]) > 1e-9) || (fabs(bVp[2]) > 1e-9),
749 "ERROR(search_bnd_circle_I_node): "
750 "projected vector is nearly identical to gc_norm_vector")
751 normalise_vector(bVp);
752
753 // the inc_angle of the bounding circle also needs to be projected
754 // onto the great circle (see routine partition data for a detailed
755 // explanation)
756 double bnd_circle_lat_cos =
757 get_vector_angle_2(bVp, bnd_circle.base_vector).cos;
758 double inc_angle =
759 M_PI_2 - acos(bnd_circle.inc_angle.sin/bnd_circle_lat_cos);
760
761 // compute "distance" of the projected base vector to the reference point
762 double base_angle = get_vector_angle(bVp, prev_gc_norm_vector);
763
764 search_interval_tree_buffer->num_overlaps = 0;
765
767 node->I.ivt.head_node, node->I.ivt.num_nodes,
768 (struct interval){
769 .left = base_angle - inc_angle,
770 .right = base_angle + inc_angle},
771 search_interval_tree_buffer);
772
773 ENSURE_ARRAY_SIZE(*overlap_cells, *overlap_cells_array_size,
774 *num_overlap_cells +
775 search_interval_tree_buffer->num_overlaps);
776
777 for (size_t i = 0; i < search_interval_tree_buffer->num_overlaps; ++i)
778 (*overlap_cells)[(*num_overlap_cells)+i] =
779 node->I.ivt.head_node[
780 search_interval_tree_buffer->overlap_iv[i]].value;
781
782 *num_overlap_cells += search_interval_tree_buffer->num_overlaps;
783
784 } else {
785
786 ENSURE_ARRAY_SIZE(*overlap_cells, *overlap_cells_array_size,
787 *num_overlap_cells + node->I_size);
788 memcpy(*overlap_cells + *num_overlap_cells, node->I.list,
789 node->I_size * sizeof(**overlap_cells));
790 *num_overlap_cells += node->I_size;
791 }
792}
793
796 struct sphere_part_node * node, struct bounding_circle bnd_circle,
797 size_t ** restrict overlap_cells, size_t * overlap_cells_array_size,
798 size_t * restrict num_overlap_cells,
799 struct overlaps * search_interval_tree_buffer, double prev_gc_norm_vector[]) {
800
801 if (node->flags & T_IS_LEAF) {
802
803 ENSURE_ARRAY_SIZE(*overlap_cells, *overlap_cells_array_size,
804 *num_overlap_cells + node->T_size);
805 memcpy(*overlap_cells + *num_overlap_cells, node->T,
806 node->T_size * sizeof(**overlap_cells));
807 *num_overlap_cells += node->T_size;
808
809 } else {
811 node->T, bnd_circle, overlap_cells, overlap_cells_array_size,
812 num_overlap_cells, search_interval_tree_buffer, node->gc_norm_vector);
813 }
814
815 if (node->flags & U_IS_LEAF) {
816
817 ENSURE_ARRAY_SIZE(*overlap_cells, *overlap_cells_array_size,
818 *num_overlap_cells + node->U_size);
819 memcpy(*overlap_cells + *num_overlap_cells, node->U,
820 node->U_size * sizeof(**overlap_cells));
821 *num_overlap_cells += node->U_size;
822
823 } else {
825 node->U, bnd_circle, overlap_cells, overlap_cells_array_size,
826 num_overlap_cells, search_interval_tree_buffer, node->gc_norm_vector);
827 }
828
830 node, bnd_circle, overlap_cells, overlap_cells_array_size,
831 num_overlap_cells, search_interval_tree_buffer, prev_gc_norm_vector);
832}
833
835 struct sphere_part_node * node, struct bounding_circle bnd_circle,
836 size_t ** restrict overlap_cells, size_t * overlap_cells_array_size,
837 size_t * restrict num_overlap_cells,
838 struct overlaps * search_interval_tree_buffer, double prev_gc_norm_vector[]) {
839
840 double dot = bnd_circle.base_vector[0] * node->gc_norm_vector[0] +
841 bnd_circle.base_vector[1] * node->gc_norm_vector[1] +
842 bnd_circle.base_vector[2] * node->gc_norm_vector[2];
843
844 // angle < M_PI_2 + bnd_circle.inc_angle
845 if (dot > - bnd_circle.inc_angle.sin) {
846
847 if (node->flags & T_IS_LEAF) {
848
849 ENSURE_ARRAY_SIZE(*overlap_cells, *overlap_cells_array_size,
850 *num_overlap_cells + node->T_size);
851 memcpy(*overlap_cells + *num_overlap_cells, node->T,
852 node->T_size * sizeof(**overlap_cells));
853 *num_overlap_cells += node->T_size;
854
855 } else {
857 node->T, bnd_circle, overlap_cells, overlap_cells_array_size,
858 num_overlap_cells, search_interval_tree_buffer,
859 node->gc_norm_vector);
860 }
861 }
862
863 // angle > M_PI_2 - bnd_circle.inc_angle
864 if (dot < bnd_circle.inc_angle.sin) {
865
866 if (node->flags & U_IS_LEAF) {
867
868 ENSURE_ARRAY_SIZE(*overlap_cells, *overlap_cells_array_size,
869 *num_overlap_cells + node->U_size);
870 memcpy(*overlap_cells + *num_overlap_cells, node->U,
871 node->U_size * sizeof(**overlap_cells));
872 *num_overlap_cells += node->U_size;
873
874 } else {
876 node->U, bnd_circle, overlap_cells, overlap_cells_array_size,
877 num_overlap_cells, search_interval_tree_buffer,
878 node->gc_norm_vector);
879 }
880 }
881
882 struct sin_cos_angle angle_sum =
883 sum_angles_no_check(node->I_angle, bnd_circle.inc_angle);
884
885 // if (I_angle + inc_angle > PI/2) ||
886 // (fabs(angle - M_PI_2) <= (I_angle + inc_angle))
887 //
888 // assumtion:
889 // I_angle >= 0 && I_angle <= PI/2
890 // inc_angle >= 0 && inc_angle <= PI/2
891 // angle >= 0 && angle <= PI
892 //
893 // => I_angle + inc_angle >= 0 && I_angle + inc_angle <= PI
894 //
895 // I_angle + inc_angle >= PI/2
896 //
897 // fabs(angle - M_PI_2) <= (I_angle + inc_angle)
898 // => sin(fabs(angle - M_PI_2)) <= sin(I_angle + inc_angle)
899 // this is wrong for (I_angle + inc_angle) > PI/2, however that case is
900 // already covered by the first condition
901 // => fabs(cos(angle)) <= sin(I_angle + inc_angle)
902 if (((angle_sum.sin < 0.0) || (angle_sum.cos <= 0.0)) ||
903 (fabs(dot) <= angle_sum.sin)) {
905 node, bnd_circle, overlap_cells, overlap_cells_array_size,
906 num_overlap_cells, search_interval_tree_buffer, prev_gc_norm_vector);
907 }
908}
909
910static void search_bnd_circle(struct sphere_part_node * node,
911 struct bounding_circle bnd_circle,
912 size_t ** restrict overlap_cells,
913 size_t * overlap_cells_array_size,
914 size_t * restrict num_overlap_cells,
915 struct overlaps * search_interval_tree_buffer,
916 double prev_gc_norm_vector[]) {
917
918 // if the bounding circle has an angle in the range of [0;PI/2[
919 if (bnd_circle.inc_angle.cos > 0.0)
921 node, bnd_circle, overlap_cells, overlap_cells_array_size,
922 num_overlap_cells, search_interval_tree_buffer, prev_gc_norm_vector);
923 else
925 node, bnd_circle, overlap_cells, overlap_cells_array_size,
926 num_overlap_cells, search_interval_tree_buffer, prev_gc_norm_vector);
927}
928
929static inline void check_leaf_NN(
930 struct point_id_xyz * points, size_t num_points,
931 double * point_coordinates_xyz, struct sin_cos_angle * best_angle,
932 double (**result_coordinates_xyz)[3],
933 size_t * result_coordinates_xyz_array_size, size_t ** local_point_ids,
934 size_t * local_point_ids_array_size, size_t total_num_local_point_ids,
935 size_t * num_local_point_ids) {
936
937 size_t * local_point_ids_ = *local_point_ids;
938 size_t local_point_ids_array_size_ = *local_point_ids_array_size;
939 double (*result_coordinates_xyz_)[3];
940 size_t result_coordinates_xyz_array_size_;
941 size_t num_local_point_ids_ = *num_local_point_ids;
942
943 if (result_coordinates_xyz != NULL) {
944 result_coordinates_xyz_ = *result_coordinates_xyz;
945 result_coordinates_xyz_array_size_ = *result_coordinates_xyz_array_size;
947 result_coordinates_xyz_, result_coordinates_xyz_array_size_,
948 total_num_local_point_ids + num_local_point_ids_ + num_points);
949 *result_coordinates_xyz = result_coordinates_xyz_;
950 *result_coordinates_xyz_array_size = result_coordinates_xyz_array_size_;
951 result_coordinates_xyz_ += total_num_local_point_ids;
952 }
954 local_point_ids_, local_point_ids_array_size_,
955 total_num_local_point_ids + num_local_point_ids_ + num_points);
956 *local_point_ids = local_point_ids_;
957 *local_point_ids_array_size = local_point_ids_array_size_;
958 local_point_ids_ += total_num_local_point_ids;
959
960
961 // check leaf for results
962 for (size_t i = 0; i < num_points; ++i) {
963
964 struct sin_cos_angle curr_angle =
966 points[i].coordinates_xyz, point_coordinates_xyz);
967 int compare = compare_angles(curr_angle, *best_angle);
968
969 // if the point is worse than the currently best point
970 if (compare > 0) continue;
971
972 // if we found a better point
973 if (compare < 0) {
974
975 *best_angle = curr_angle;
976 num_local_point_ids_ = 1;
977 if (result_coordinates_xyz != NULL) {
978 result_coordinates_xyz_[0][0] = points[i].coordinates_xyz[0];
979 result_coordinates_xyz_[0][1] = points[i].coordinates_xyz[1];
980 result_coordinates_xyz_[0][2] = points[i].coordinates_xyz[2];
981 }
982 local_point_ids_[0] = points[i].idx;
983
984 } else {
985
986 if (result_coordinates_xyz != NULL) {
987 result_coordinates_xyz_[num_local_point_ids_][0] =
988 points[i].coordinates_xyz[0];
989 result_coordinates_xyz_[num_local_point_ids_][1] =
990 points[i].coordinates_xyz[1];
991 result_coordinates_xyz_[num_local_point_ids_][2] =
992 points[i].coordinates_xyz[2];
993 }
994 local_point_ids_[num_local_point_ids_] = points[i].idx;
995 num_local_point_ids_++;
996 }
997 }
998
999 *num_local_point_ids = num_local_point_ids_;
1000}
1001
1003 struct bounding_circle * bnd_circle, double (**result_coordinates_xyz)[3],
1004 size_t * result_coordinates_xyz_array_size, size_t ** local_point_ids,
1005 size_t * local_point_ids_array_size, size_t total_num_local_point_ids,
1006 size_t * num_local_point_ids, double * dot_stack,
1007 struct point_sphere_part_node ** node_stack,
1008 int * flags, size_t curr_tree_depth) {
1009
1010 double * point_coordinates_xyz = bnd_circle->base_vector;
1011 struct sin_cos_angle best_angle = bnd_circle->inc_angle;
1012
1013 double dot = dot_stack[curr_tree_depth];
1014 struct point_sphere_part_node * node = node_stack[curr_tree_depth];
1015 int skip_U = flags[curr_tree_depth] & U_FLAG;
1016 int skip_T = flags[curr_tree_depth] & T_FLAG;
1017
1018 do {
1019
1020 if (!skip_U) {
1021
1022 flags[curr_tree_depth] |= U_FLAG;
1023
1024 // angle + inc_angle >= M_PI_2
1025 if ((dot < best_angle.sin) | (best_angle.cos <= 0.0)) {
1026
1027 if (node->flags & U_IS_LEAF) {
1028
1030 (struct point_id_xyz *)(node->U), node->U_size,
1031 point_coordinates_xyz, &best_angle, result_coordinates_xyz,
1032 result_coordinates_xyz_array_size, local_point_ids,
1033 local_point_ids_array_size, total_num_local_point_ids,
1034 num_local_point_ids);
1035
1036 } else {
1037
1038 // traverse down one level
1039 ++curr_tree_depth;
1040 node = (struct point_sphere_part_node *)(node->U);
1041 dot = node->gc_norm_vector[0] * point_coordinates_xyz[0] +
1042 node->gc_norm_vector[1] * point_coordinates_xyz[1] +
1043 node->gc_norm_vector[2] * point_coordinates_xyz[2];
1044 dot_stack[curr_tree_depth] = dot;
1045 node_stack[curr_tree_depth] = node;
1046 flags[curr_tree_depth] = 0;
1047 // skip_U = 0;
1048 skip_T = 0;
1049 continue;
1050 }
1051 }
1052 }
1053
1054 if (!skip_T) {
1055
1056 flags[curr_tree_depth] = U_FLAG + T_FLAG;
1057
1058 // angle - inc_angle < M_PI_2
1059 if ((dot > - best_angle.sin) || (best_angle.cos <= 0.0)) {
1060
1061 if (node->flags & T_IS_LEAF) {
1062
1064 (struct point_id_xyz *)(node->T), node->T_size,
1065 point_coordinates_xyz, &best_angle, result_coordinates_xyz,
1066 result_coordinates_xyz_array_size, local_point_ids,
1067 local_point_ids_array_size, total_num_local_point_ids,
1068 num_local_point_ids);
1069
1070 } else {
1071
1072 // traverse down one level
1073 ++curr_tree_depth;
1074 node = (struct point_sphere_part_node *)(node->T);
1075 dot = node->gc_norm_vector[0] * point_coordinates_xyz[0] +
1076 node->gc_norm_vector[1] * point_coordinates_xyz[1] +
1077 node->gc_norm_vector[2] * point_coordinates_xyz[2];
1078 dot_stack[curr_tree_depth] = dot;
1079 node_stack[curr_tree_depth] = node;
1080 flags[curr_tree_depth] = 0;
1081 skip_U = 0;
1082 // skip_T = 0;
1083 continue;
1084 }
1085 }
1086 }
1087
1088 if (curr_tree_depth == 0) break;
1089
1090 // go up one level in the tree
1091
1092 curr_tree_depth--;
1093 dot = dot_stack[curr_tree_depth];
1094 node = node_stack[curr_tree_depth];
1095 skip_U = flags[curr_tree_depth] & U_FLAG;
1096 skip_T = flags[curr_tree_depth] & T_FLAG;
1097
1098 } while (1);
1099
1100 bnd_circle->inc_angle = best_angle;
1101}
1102
1103
1105 struct point_id_xyz * points, size_t num_points, double coordinate_xyz[3],
1106 size_t ** local_point_ids, size_t * local_point_ids_array_size,
1107 double (**result_coordinates_xyz)[3],
1108 size_t * result_coordinates_xyz_array_size,
1109 size_t total_num_local_point_ids, size_t * num_local_point_ids) {
1110
1111 for (size_t i = 0; i < num_points; ++i) {
1112
1113 // if the points are nearly identical
1114 if (points_are_identically(points[i].coordinates_xyz, coordinate_xyz)) {
1115
1116 ENSURE_ARRAY_SIZE(*local_point_ids, *local_point_ids_array_size,
1117 total_num_local_point_ids + 1);
1118 (*local_point_ids)[total_num_local_point_ids] = points[i].idx;
1119 if (result_coordinates_xyz != NULL) {
1120 ENSURE_ARRAY_SIZE(*result_coordinates_xyz,
1121 *result_coordinates_xyz_array_size,
1122 total_num_local_point_ids + 1);
1123 memcpy((*result_coordinates_xyz) + total_num_local_point_ids,
1124 points[i].coordinates_xyz, 3 * sizeof(double));
1125 }
1126 *num_local_point_ids = 1;
1127 return 1;
1128 }
1129 }
1130
1131 return 0;
1132}
1133
1135 size_t num_points,
1136 double (*coordinates_xyz)[3],
1137 double * cos_angles,
1138 double (**result_coordinates_xyz)[3],
1139 size_t * result_coordinates_xyz_array_size,
1140 size_t ** local_point_ids,
1141 size_t * local_point_ids_array_size,
1142 size_t * num_local_point_ids) {
1143
1144 memset(num_local_point_ids, 0, num_points * sizeof(*num_local_point_ids));
1145
1146 if (search == NULL) return;
1147
1148 struct point_sphere_part_node * base_node = &(search->base_node);
1149
1150 size_t total_num_local_point_ids = 0;
1151
1152 double * dot_stack = xmalloc(search->max_tree_depth * sizeof(*dot_stack));
1153 struct point_sphere_part_node ** node_stack =
1154 xmalloc(search->max_tree_depth * sizeof(*node_stack));
1155 int * flags = xmalloc(search->max_tree_depth * sizeof(*flags));
1156
1157 for (size_t i = 0; i < num_points; ++i) {
1158
1159 struct point_sphere_part_node * curr_node = base_node;
1160
1161 double * curr_coordinates_xyz = coordinates_xyz[i];
1162
1163 size_t curr_tree_depth = 0;
1164 struct point_id_xyz * points = NULL;
1165 size_t curr_num_points = 0;
1166
1167 // get the matching leaf for the current point
1168 do {
1169
1170 double dot = curr_node->gc_norm_vector[0]*curr_coordinates_xyz[0] +
1171 curr_node->gc_norm_vector[1]*curr_coordinates_xyz[1] +
1172 curr_node->gc_norm_vector[2]*curr_coordinates_xyz[2];
1173
1174 dot_stack[curr_tree_depth] = dot;
1175 node_stack[curr_tree_depth] = curr_node;
1176 flags[curr_tree_depth] = 0;
1177
1178 // angle > M_PI_2
1179 if (dot < yac_angle_tol) {
1180
1181 flags[curr_tree_depth] = U_FLAG;
1182
1183 if (curr_node->flags & U_IS_LEAF) {
1184 if (curr_node->U_size > 0) {
1185 points = (struct point_id_xyz*)(curr_node->U);
1186 curr_num_points = curr_node->U_size;
1187 break;
1188 } else {
1189 flags[curr_tree_depth] = T_FLAG;
1190 YAC_ASSERT(
1191 curr_node->flags & T_IS_LEAF,
1192 "ERROR(yac_point_sphere_part_search_NN): "
1193 "if one branch is empty, the other has to be a leaf");
1194 points = (struct point_id_xyz*)(curr_node->T);
1195 curr_num_points = curr_node->T_size;
1196 break;
1197 }
1198 } else curr_node = curr_node->U;
1199
1200 // angle < M_PI_2
1201 } else if (dot > -yac_angle_tol) {
1202
1203 flags[curr_tree_depth] = T_FLAG;
1204
1205 if (curr_node->flags & T_IS_LEAF) {
1206 if (curr_node->T_size > 0) {
1207 points = (struct point_id_xyz*)(curr_node->T);
1208 curr_num_points = curr_node->T_size;
1209 break;
1210 } else {
1211 flags[curr_tree_depth] = U_FLAG;
1212 YAC_ASSERT(
1213 curr_node->flags & U_IS_LEAF,
1214 "ERROR(yac_point_sphere_part_search_NN): "
1215 "if one branch is empty, the other has to be a leaf");
1216 points = (struct point_id_xyz*)(curr_node->U);
1217 curr_num_points = curr_node->U_size;
1218 break;
1219 }
1220 } else curr_node = curr_node->T;
1221 }
1222
1223 curr_tree_depth++;
1224 } while (1);
1225
1226 // if we do not have to do a finer search
1228 points, curr_num_points, curr_coordinates_xyz, local_point_ids,
1229 local_point_ids_array_size, result_coordinates_xyz,
1230 result_coordinates_xyz_array_size, total_num_local_point_ids,
1231 num_local_point_ids + i)) {
1232
1233 if (cos_angles != NULL) cos_angles[i] = 1.0;
1234
1235 } else {
1236
1237 struct bounding_circle bnd_circle;
1238 bnd_circle.base_vector[0] = curr_coordinates_xyz[0];
1239 bnd_circle.base_vector[1] = curr_coordinates_xyz[1];
1240 bnd_circle.base_vector[2] = curr_coordinates_xyz[2];
1241 bnd_circle.inc_angle = SIN_COS_M_PI;
1242 bnd_circle.sq_crd = DBL_MAX;
1243
1245 points, curr_num_points, curr_coordinates_xyz, &bnd_circle.inc_angle,
1246 result_coordinates_xyz, result_coordinates_xyz_array_size,
1247 local_point_ids, local_point_ids_array_size, total_num_local_point_ids,
1248 num_local_point_ids + i);
1249
1250 // get best result points
1252 &bnd_circle, result_coordinates_xyz, result_coordinates_xyz_array_size,
1253 local_point_ids, local_point_ids_array_size, total_num_local_point_ids,
1254 num_local_point_ids + i, dot_stack, node_stack, flags, curr_tree_depth);
1255
1256 if (cos_angles != NULL) cos_angles[i] = bnd_circle.inc_angle.cos;
1257 }
1258
1259 total_num_local_point_ids += num_local_point_ids[i];
1260 }
1261
1262 free(flags);
1263 free(node_stack);
1264 free(dot_stack);
1265}
1266
1267static inline int compare_point_id_xyz_angle(const void * a, const void * b) {
1268
1269 const struct point_id_xyz_angle * p_a = (const struct point_id_xyz_angle *)a;
1270 const struct point_id_xyz_angle * p_b = (const struct point_id_xyz_angle *)b;
1271
1272 int ret = (p_a->cos_angle < p_b->cos_angle) -
1273 (p_a->cos_angle > p_b->cos_angle);
1274
1275 if (ret != 0) return ret;
1276
1277 return (p_a->point.idx > p_b->point.idx) - (p_a->point.idx < p_b->point.idx);
1278}
1279
1281 size_t n, struct point_id_xyz * points, size_t num_points,
1282 double * point_coordinates_xyz, struct point_id_xyz_angle ** results,
1283 size_t * results_array_size) {
1284
1285 assert(num_points > 0);
1286
1287 ENSURE_ARRAY_SIZE(*results, *results_array_size, num_points);
1288 struct point_id_xyz_angle * results_ = *results;
1289
1290#ifdef __NEC__
1291// vectorization of the following loop leads to failues of
1292//
1293// test_couple_config_parallel1
1294// test_interp_method_nnn_parallel
1295// test_interp_method_hcsbb_parallel
1296//
1297// with NEC compiler when CFLAGS='-O2'
1298#pragma _NEC novector
1299#endif
1300 for (size_t i = 0; i < num_points; ++i) {
1301
1302 results_[i].point = points[i];
1303 results_[i].cos_angle = clamp_abs_one(
1304 points[i].coordinates_xyz[0] * point_coordinates_xyz[0] +
1305 points[i].coordinates_xyz[1] * point_coordinates_xyz[1] +
1306 points[i].coordinates_xyz[2] * point_coordinates_xyz[2]);
1307 }
1308
1309 qsort(results_, num_points, sizeof(*results_), compare_point_id_xyz_angle);
1310
1311 if (num_points <= n) return num_points;
1312
1313 size_t num_results;
1314 double min_cos_angle = results_[n - 1].cos_angle;
1315
1316 for (num_results = n;
1317 (num_results < num_points) &&
1318 !(fabs(min_cos_angle - results_[num_results].cos_angle) > 0.0);
1319 ++num_results);
1320
1321 return num_results;
1322}
1323
1324static inline struct sin_cos_angle check_leaf_NNN(
1325 size_t n, double * point_coordinates_xyz,
1326 struct point_id_xyz * points, size_t num_points,
1327 struct point_id_xyz_angle ** results, size_t * results_array_size,
1328 size_t * num_results, struct sin_cos_angle curr_angle) {
1329
1330 size_t num_results_ = *num_results;
1331 ENSURE_ARRAY_SIZE(*results, *results_array_size, num_results_ + num_points);
1332 struct point_id_xyz_angle * results_ = *results;
1333
1334 int flag = 0;
1335
1336 double min_cos_angle = results_[num_results_-1].cos_angle;
1337
1338 // check leaf for results
1339 for (size_t i = 0; i < num_points; ++i) {
1340
1341 double curr_cos_angle =
1342 points[i].coordinates_xyz[0] * point_coordinates_xyz[0] +
1343 points[i].coordinates_xyz[1] * point_coordinates_xyz[1] +
1344 points[i].coordinates_xyz[2] * point_coordinates_xyz[2];
1345
1346 // if the point is worse than the currently best point
1347 if (curr_cos_angle < min_cos_angle) continue;
1348
1349 struct point_id_xyz_angle point =
1350 {.point = points[i], .cos_angle = curr_cos_angle};
1351
1352 // insert point
1353 size_t j;
1354 for (j = 0; j < num_results_; ++j) {
1355
1357 &point, results_ + num_results_ - j - 1) < 0) {
1358 results_[num_results_ - j] = results_[num_results_ - j - 1];
1359 } else {
1360 break;
1361 }
1362 }
1363 results_[num_results_ - j] = point;
1364
1365 ++num_results_;
1366 flag = 1;
1367 }
1368
1369 if (flag) {
1370
1371 if (num_results_ > n) {
1372
1373 size_t new_num_results;
1374 min_cos_angle = results_[n - 1].cos_angle;
1375
1376 for (new_num_results = n;
1377 (new_num_results < num_results_) &&
1378 !(fabs(min_cos_angle - results_[new_num_results].cos_angle) > 0.0);
1379 ++new_num_results);
1380 num_results_ = new_num_results;
1381 }
1382 *num_results = num_results_;
1383
1384 return
1386 results_[num_results_-1].point.coordinates_xyz, point_coordinates_xyz);
1387 } else return curr_angle;
1388}
1389
1391 size_t n, double * point_coordinates_xyz,
1392 struct point_id_xyz_angle ** results, size_t * results_array_size,
1393 size_t * num_results, double * dot_stack,
1394 struct point_sphere_part_node ** node_stack, int * flags,
1395 size_t curr_tree_depth) {
1396
1397 struct sin_cos_angle angle =
1399 (*results)[(*num_results)-1].point.coordinates_xyz,
1400 point_coordinates_xyz);
1401
1402 // if we have already found at least n exactly matching points
1403 if ((*num_results >= n) && (angle.sin <= yac_angle_tol)) return;
1404
1405 double dot = dot_stack[curr_tree_depth];
1406 struct point_sphere_part_node * node = node_stack[curr_tree_depth];
1407 int skip_U = flags[curr_tree_depth] & U_FLAG;
1408 int skip_T = flags[curr_tree_depth] & T_FLAG;
1409
1410 do {
1411
1412 if (!skip_U) {
1413
1414 flags[curr_tree_depth] |= U_FLAG;
1415
1416 // angle + inc_angle >= M_PI_2
1417 if ((dot < angle.sin) | (angle.cos <= 0.0)) {
1418
1419 if (node->flags & U_IS_LEAF) {
1420
1421 angle = check_leaf_NNN(
1422 n, point_coordinates_xyz, (struct point_id_xyz *)(node->U),
1423 node->U_size, results, results_array_size, num_results, angle);
1424
1425 } else {
1426
1427 // traverse down one level
1428 ++curr_tree_depth;
1429 node = (struct point_sphere_part_node *)(node->U);
1430 dot = node->gc_norm_vector[0] * point_coordinates_xyz[0] +
1431 node->gc_norm_vector[1] * point_coordinates_xyz[1] +
1432 node->gc_norm_vector[2] * point_coordinates_xyz[2];
1433 dot_stack[curr_tree_depth] = dot;
1434 node_stack[curr_tree_depth] = node;
1435 flags[curr_tree_depth] = 0;
1436 // skip_U = 0;
1437 skip_T = 0;
1438 continue;
1439 }
1440 }
1441 }
1442
1443 if (!skip_T) {
1444
1445 flags[curr_tree_depth] = U_FLAG + T_FLAG;
1446
1447 // angle - inc_angle < M_PI_2
1448 if ((dot > - angle.sin) || (angle.cos <= 0.0)) {
1449
1450 if (node->flags & T_IS_LEAF) {
1451
1452 angle = check_leaf_NNN(
1453 n, point_coordinates_xyz, (struct point_id_xyz *)(node->T),
1454 node->T_size, results, results_array_size, num_results, angle);
1455
1456 } else {
1457
1458 // traverse down one level
1459 ++curr_tree_depth;
1460 node = (struct point_sphere_part_node *)(node->T);
1461 dot = node->gc_norm_vector[0] * point_coordinates_xyz[0] +
1462 node->gc_norm_vector[1] * point_coordinates_xyz[1] +
1463 node->gc_norm_vector[2] * point_coordinates_xyz[2];
1464 dot_stack[curr_tree_depth] = dot;
1465 node_stack[curr_tree_depth] = node;
1466 flags[curr_tree_depth] = 0;
1467 skip_U = 0;
1468 // skip_T = 0;
1469 continue;
1470 }
1471 }
1472 }
1473
1474 if (curr_tree_depth == 0) break;
1475
1476 // go up one level in the tree
1477
1478 curr_tree_depth--;
1479 dot = dot_stack[curr_tree_depth];
1480 node = node_stack[curr_tree_depth];
1481 skip_U = flags[curr_tree_depth] & U_FLAG;
1482 skip_T = flags[curr_tree_depth] & T_FLAG;
1483
1484 } while (1);
1485}
1486
1488 size_t num_points,
1489 double (*coordinates_xyz)[3], size_t n,
1490 double ** cos_angles,
1491 size_t * cos_angles_array_size,
1492 double (**result_coordinates_xyz)[3],
1493 size_t * result_coordinates_xyz_array_size,
1494 size_t ** local_point_ids,
1495 size_t * local_point_ids_array_size,
1496 size_t * num_local_point_ids) {
1497
1498 if (num_points == 0) return;
1499
1500 if (cos_angles != NULL)
1501 ENSURE_ARRAY_SIZE(*cos_angles, *cos_angles_array_size, num_points * n);
1502
1503 if (n == 1) {
1505 search, num_points, coordinates_xyz, (cos_angles!=NULL)?*cos_angles:NULL,
1506 result_coordinates_xyz, result_coordinates_xyz_array_size,
1507 local_point_ids, local_point_ids_array_size, num_local_point_ids);
1508
1509 size_t total_num_local_points = 0;
1510 for (size_t i = 0; i < num_points; ++i)
1511 total_num_local_points += num_local_point_ids[i];
1512
1513 if ((cos_angles != NULL) && (total_num_local_points > num_points)) {
1514
1515 ENSURE_ARRAY_SIZE(*cos_angles, *cos_angles_array_size,
1516 total_num_local_points);
1517
1518 for (size_t i = num_points - 1, offset = total_num_local_points - 1;
1519 i != (size_t)-1; i--) {
1520
1521 for (size_t j = 0; j < num_local_point_ids[i]; ++j, --offset)
1522 (*cos_angles)[offset] = (*cos_angles)[i];
1523 }
1524 }
1525 return;
1526 }
1527
1528
1529 if (search == NULL) {
1530 memset(num_local_point_ids, 0, num_points * sizeof(*num_local_point_ids));
1531 return;
1532 }
1533
1534 struct point_sphere_part_node * base_node = &(search->base_node);
1535
1536 size_t total_num_local_point_ids = 0;
1537
1538 double * dot_stack = xmalloc(search->max_tree_depth * sizeof(*dot_stack));
1539 struct point_sphere_part_node ** node_stack =
1540 xmalloc(search->max_tree_depth * sizeof(*node_stack));
1541 int * flags = xmalloc(search->max_tree_depth * sizeof(*flags));
1542
1543 struct point_id_xyz_angle * results = NULL;
1544 size_t results_array_size = 0;
1545
1546 for (size_t i = 0; i < num_points; ++i) {
1547
1548 struct point_sphere_part_node * curr_node = base_node;
1549
1550 double * curr_coordinates_xyz = coordinates_xyz[i];
1551
1552 size_t curr_tree_depth = 0;
1553 struct point_id_xyz * points = search->points;
1554 size_t curr_num_points = 0;
1555
1556 // get the matching leaf for the current point
1557 do {
1558
1559 double dot = curr_node->gc_norm_vector[0]*curr_coordinates_xyz[0] +
1560 curr_node->gc_norm_vector[1]*curr_coordinates_xyz[1] +
1561 curr_node->gc_norm_vector[2]*curr_coordinates_xyz[2];
1562
1563 dot_stack[curr_tree_depth] = dot;
1564 node_stack[curr_tree_depth] = curr_node;
1565 flags[curr_tree_depth] = 0;
1566
1567 // angle >= M_PI_2
1568 if (dot <= 0.0) {
1569
1570 if (curr_node->U_size < n) {
1571
1572 flags[curr_tree_depth] = U_FLAG + T_FLAG;
1573 curr_num_points = curr_node->U_size + curr_node->T_size;
1574 break;
1575 } else if (curr_node->flags & U_IS_LEAF) {
1576
1577 flags[curr_tree_depth] = U_FLAG;
1578 curr_num_points = curr_node->U_size;
1579 break;
1580 } else {
1581
1582 flags[curr_tree_depth] = U_FLAG;
1583 curr_node = curr_node->U;
1584 }
1585
1586 } else {
1587
1588 if (curr_node->T_size < n) {
1589
1590 flags[curr_tree_depth] = U_FLAG + T_FLAG;
1591 curr_num_points = curr_node->U_size + curr_node->T_size;
1592 break;
1593 } else if (curr_node->flags & T_IS_LEAF) {
1594
1595 points += curr_node->U_size;
1596 flags[curr_tree_depth] = T_FLAG;
1597 curr_num_points = curr_node->T_size;
1598 break;
1599 } else {
1600
1601 points += curr_node->U_size;
1602 flags[curr_tree_depth] = T_FLAG;
1603 curr_node = curr_node->T;
1604 }
1605 }
1606
1607 curr_tree_depth++;
1608 } while (1);
1609
1610 assert(curr_num_points > 0);
1611
1612 size_t num_results =
1614 n, points, curr_num_points, curr_coordinates_xyz,
1615 &results, &results_array_size);
1616
1617 // do a detailed search
1619 n, curr_coordinates_xyz, &results, &results_array_size, &num_results,
1620 dot_stack, node_stack, flags, curr_tree_depth);
1621
1622 // extract the results
1623 ENSURE_ARRAY_SIZE(*local_point_ids, *local_point_ids_array_size,
1624 total_num_local_point_ids + num_results);
1625 size_t * local_point_ids_ =
1626 (*local_point_ids) + total_num_local_point_ids;
1627 double * cos_angles_;
1628 if (cos_angles != NULL) {
1629 ENSURE_ARRAY_SIZE(*cos_angles, *cos_angles_array_size,
1630 total_num_local_point_ids + num_results);
1631 cos_angles_ = (*cos_angles) + total_num_local_point_ids;
1632 } else {
1633 cos_angles_ = NULL;
1634 }
1635 double (*result_coordinates_xyz_)[3];
1636 if (result_coordinates_xyz != NULL) {
1637 ENSURE_ARRAY_SIZE(*result_coordinates_xyz,
1638 *result_coordinates_xyz_array_size,
1639 total_num_local_point_ids + num_results);
1640 result_coordinates_xyz_ =
1641 (*result_coordinates_xyz) + total_num_local_point_ids;
1642 } else {
1643 result_coordinates_xyz_ = NULL;
1644 }
1645
1646 for (size_t j = 0; j < num_results; ++j) {
1647
1648 local_point_ids_[j] = results[j].point.idx;
1649 if (cos_angles_ != NULL) cos_angles_[j] = results[j].cos_angle;
1650 if (result_coordinates_xyz_ != NULL) {
1651 result_coordinates_xyz_[j][0] = results[j].point.coordinates_xyz[0];
1652 result_coordinates_xyz_[j][1] = results[j].point.coordinates_xyz[1];
1653 result_coordinates_xyz_[j][2] = results[j].point.coordinates_xyz[2];
1654 }
1655 }
1656
1657 num_local_point_ids[i] = num_results;
1658 total_num_local_point_ids += num_results;
1659 }
1660
1661 free(results);
1662 free(flags);
1663 free(node_stack);
1664 free(dot_stack);
1665}
1666
1668 struct point_sphere_part_search * search,
1669 size_t num_bnd_circles, struct bounding_circle * bnd_circles,
1670 size_t n, size_t ** local_point_ids, size_t * local_point_ids_array_size,
1671 size_t * num_local_point_ids) {
1672
1673 if (num_bnd_circles == 0) return;
1674
1675 if (search == NULL) {
1676 memset(
1677 num_local_point_ids, 0, num_bnd_circles * sizeof(*num_local_point_ids));
1678 return;
1679 }
1680
1681 struct point_sphere_part_node * base_node = &(search->base_node);
1682
1683 size_t total_num_local_point_ids = 0;
1684
1685 double * dot_stack = xmalloc(search->max_tree_depth * sizeof(*dot_stack));
1686 struct point_sphere_part_node ** node_stack =
1687 xmalloc(search->max_tree_depth * sizeof(*node_stack));
1688 int * flags = xmalloc(search->max_tree_depth * sizeof(*flags));
1689
1690 struct point_id_xyz_angle * results = NULL;
1691 size_t results_array_size = 0;
1692
1693 for (size_t i = 0; i < num_bnd_circles; ++i) {
1694
1695 struct point_sphere_part_node * curr_node = base_node;
1696
1697 double * curr_coordinates_xyz = bnd_circles[i].base_vector;
1698
1699 size_t curr_tree_depth = 0;
1700 struct point_id_xyz * points = search->points;
1701 size_t curr_num_points = 0;
1702
1703 // get the matching leaf for the current point
1704 do {
1705
1706 double dot = curr_node->gc_norm_vector[0]*curr_coordinates_xyz[0] +
1707 curr_node->gc_norm_vector[1]*curr_coordinates_xyz[1] +
1708 curr_node->gc_norm_vector[2]*curr_coordinates_xyz[2];
1709
1710 dot_stack[curr_tree_depth] = dot;
1711 node_stack[curr_tree_depth] = curr_node;
1712 flags[curr_tree_depth] = 0;
1713
1714 // angle >= M_PI_2
1715 if (dot <= 0.0) {
1716
1717 if (curr_node->U_size < n) {
1718
1719 flags[curr_tree_depth] = U_FLAG + T_FLAG;
1720 curr_num_points = curr_node->U_size + curr_node->T_size;
1721 break;
1722 } else if (curr_node->flags & U_IS_LEAF) {
1723
1724 flags[curr_tree_depth] = U_FLAG;
1725 curr_num_points = curr_node->U_size;
1726 break;
1727 } else {
1728
1729 flags[curr_tree_depth] = U_FLAG;
1730 curr_node = curr_node->U;
1731 }
1732
1733 } else {
1734
1735 if (curr_node->T_size < n) {
1736
1737 flags[curr_tree_depth] = U_FLAG + T_FLAG;
1738 curr_num_points = curr_node->U_size + curr_node->T_size;
1739 break;
1740 } else if (curr_node->flags & T_IS_LEAF) {
1741
1742 points += curr_node->U_size;
1743 flags[curr_tree_depth] = T_FLAG;
1744 curr_num_points = curr_node->T_size;
1745 break;
1746 } else {
1747
1748 points += curr_node->U_size;
1749 flags[curr_tree_depth] = T_FLAG;
1750 curr_node = curr_node->T;
1751 }
1752 }
1753
1754 curr_tree_depth++;
1755 } while (1);
1756
1757 YAC_ASSERT(
1758 curr_num_points > 0,
1759 "ERROR(yac_point_sphere_part_search_NNN_bnd_circle): "
1760 "insufficient number of points");
1761
1762 size_t num_results =
1764 n, points, curr_num_points, curr_coordinates_xyz,
1765 &results, &results_array_size);
1766
1767 // do a detailed search
1769 n, curr_coordinates_xyz, &results, &results_array_size, &num_results,
1770 dot_stack, node_stack, flags, curr_tree_depth);
1771
1772 for (; num_results > 0; --num_results)
1773 if (results[num_results-1].cos_angle >= bnd_circles[i].inc_angle.cos)
1774 break;
1775
1776 // extract the results
1777 ENSURE_ARRAY_SIZE(*local_point_ids, *local_point_ids_array_size,
1778 total_num_local_point_ids + num_results);
1779 size_t * local_point_ids_ =
1780 (*local_point_ids) + total_num_local_point_ids;
1781
1782 for (size_t j = 0; j < num_results; ++j)
1783 local_point_ids_[j] = results[j].point.idx;
1784
1785 num_local_point_ids[i] = num_results;
1786 total_num_local_point_ids += num_results;
1787 }
1788
1789 free(results);
1790 free(flags);
1791 free(node_stack);
1792 free(dot_stack);
1793}
1794
1795static int compare_angles_(void const * a, void const * b) {
1796 struct sin_cos_angle * a_ = (struct sin_cos_angle *)a;
1797 struct sin_cos_angle * b_ = (struct sin_cos_angle *)b;
1798 return compare_angles(*a_, *b_);
1799}
1800
1802 struct point_sphere_part_search * search,
1803 size_t num_points, yac_coordinate_pointer coordinates_xyz,
1804 size_t n, struct sin_cos_angle * angles) {
1805
1806 YAC_ASSERT(
1807 search != NULL,
1808 "ERRROR(yac_point_sphere_part_search_NNN_ubound): "
1809 "invalid point sphere part search (has to be != NULL)");
1810 YAC_ASSERT(
1811 n > 0, "ERROR(yac_point_sphere_part_search_NNN_ubound): "
1812 "invalid n (has to be > 0)")
1813
1814 struct sin_cos_angle * temp_angles = NULL;
1815 size_t temp_angles_array_size = 0;
1816
1817 struct point_sphere_part_node * base_node = &(search->base_node);
1818
1819 for (size_t i = 0; i < num_points; ++i) {
1820
1821 struct point_sphere_part_node * curr_node = base_node;
1822
1823 double * curr_coordinates_xyz = coordinates_xyz[i];
1824
1825 struct point_id_xyz * points = search->points;
1826 size_t curr_num_points = 0;
1827
1828 // get the matching leaf for the current point
1829 do {
1830
1831 double dot = curr_node->gc_norm_vector[0]*curr_coordinates_xyz[0] +
1832 curr_node->gc_norm_vector[1]*curr_coordinates_xyz[1] +
1833 curr_node->gc_norm_vector[2]*curr_coordinates_xyz[2];
1834
1835 // angle >= M_PI_2
1836 if (dot <= 0.0) {
1837
1838 if (curr_node->U_size < n) {
1839
1840 curr_num_points = curr_node->U_size + curr_node->T_size;
1841 break;
1842 } else if (curr_node->flags & U_IS_LEAF) {
1843
1844 curr_num_points = curr_node->U_size;
1845 break;
1846 } else {
1847
1848 curr_node = curr_node->U;
1849 }
1850
1851 } else {
1852
1853 if (curr_node->T_size < n) {
1854
1855 curr_num_points = curr_node->U_size + curr_node->T_size;
1856 break;
1857 } else if (curr_node->flags & T_IS_LEAF) {
1858
1859 points += curr_node->U_size;
1860 curr_num_points = curr_node->T_size;
1861 break;
1862 } else {
1863
1864 points += curr_node->U_size;
1865 curr_node = curr_node->T;
1866 }
1867 }
1868 } while (1);
1869
1870 YAC_ASSERT(
1871 curr_num_points >= n,
1872 "ERROR(yac_point_sphere_part_search_NNN_ubound): "
1873 "failed to find a sufficient number of points");
1874
1875 // search of the closest "n" points in the current
1876 // list of points
1877
1878 if (n == 1) {
1879
1880 struct sin_cos_angle best_angle =
1882 curr_coordinates_xyz, points[0].coordinates_xyz);
1883 for (size_t j = 1; j < curr_num_points; ++j) {
1884 struct sin_cos_angle curr_angle =
1886 curr_coordinates_xyz, points[j].coordinates_xyz);
1887 if (compare_angles(best_angle, curr_angle) > 0)
1888 best_angle = curr_angle;
1889 }
1890 angles[i] = best_angle;
1891
1892 } else {
1893
1894 // compute the angles for all current points
1896 temp_angles, temp_angles_array_size, curr_num_points);
1897 for (size_t j = 0; j < curr_num_points; ++j)
1898 temp_angles[j] =
1900 curr_coordinates_xyz, points[j].coordinates_xyz);
1901
1902 qsort(
1903 temp_angles, curr_num_points, sizeof(*temp_angles),
1905
1906 angles[i] = temp_angles[n-1];
1907 }
1908 }
1909
1910 free(temp_angles);
1911}
1912
1913static void search_point(struct sphere_part_node * node,
1914 double point[],
1915 size_t ** overlap_cells,
1916 size_t * overlap_cells_array_size,
1917 size_t * num_overlap_cells,
1918 struct overlaps * search_interval_tree_buffer,
1919 double prev_gc_norm_vector[]) {
1920
1921 double dot = point[0] * node->gc_norm_vector[0] +
1922 point[1] * node->gc_norm_vector[1] +
1923 point[2] * node->gc_norm_vector[2];
1924
1925 // angle < M_PI_2
1926 if (dot > -yac_angle_tol) {
1927
1928 if (node->flags & T_IS_LEAF) {
1929
1930 ENSURE_ARRAY_SIZE(*overlap_cells, *overlap_cells_array_size,
1931 *num_overlap_cells + node->T_size);
1932 memcpy(*overlap_cells + *num_overlap_cells, node->T,
1933 node->T_size * sizeof(**overlap_cells));
1934 *num_overlap_cells += node->T_size;
1935
1936 } else {
1937 search_point(node->T, point, overlap_cells,
1938 overlap_cells_array_size, num_overlap_cells,
1939 search_interval_tree_buffer, node->gc_norm_vector);
1940 }
1941 }
1942
1943 // angle > M_PI_2
1944 if (dot < yac_angle_tol) {
1945
1946 if (node->flags & U_IS_LEAF) {
1947
1948 ENSURE_ARRAY_SIZE(*overlap_cells, *overlap_cells_array_size,
1949 *num_overlap_cells + node->U_size);
1950 memcpy(*overlap_cells + *num_overlap_cells, node->U,
1951 node->U_size * sizeof(**overlap_cells));
1952 *num_overlap_cells += node->U_size;
1953
1954 } else {
1955 search_point(node->U, point, overlap_cells,
1956 overlap_cells_array_size, num_overlap_cells,
1957 search_interval_tree_buffer, node->gc_norm_vector);
1958 }
1959 }
1960
1961 // fabs(angle - M_PI_2) <= (node->I_angle)
1962 // fabs(cos(angle)) <= sin(node->I_angle)
1963 if (fabs(dot) <= node->I_angle.sin) {
1964
1965 if (node->flags & I_IS_INTERVAL_TREE) {
1966
1967 // project the point onto the current great circle
1968 double GCp[3], bVp[3];
1969 crossproduct_kahan(node->gc_norm_vector, point, GCp);
1970 crossproduct_kahan(GCp, node->gc_norm_vector, bVp);
1971
1972 struct interval search_interval;
1973 // if the projected point does not coincide with the norm vector of
1974 // the plane through the great circle
1975 if ((fabs(bVp[0]) > 1e-9) || (fabs(bVp[1]) > 1e-9) ||
1976 (fabs(bVp[2]) > 1e-9)) {
1977 normalise_vector(bVp);
1978 double base_angle = get_vector_angle(bVp, prev_gc_norm_vector);
1979 search_interval.left=base_angle;
1980 search_interval.right=base_angle;
1981 } else {
1982 search_interval.left = -M_PI;
1983 search_interval.right = M_PI;
1984 }
1985
1986 search_interval_tree_buffer->num_overlaps = 0;
1987
1989 search_interval, search_interval_tree_buffer);
1990
1991 ENSURE_ARRAY_SIZE(*overlap_cells, *overlap_cells_array_size,
1992 *num_overlap_cells +
1993 search_interval_tree_buffer->num_overlaps);
1994
1995 for (size_t i = 0; i < search_interval_tree_buffer->num_overlaps;
1996 ++i) {
1997 (*overlap_cells)[(*num_overlap_cells)+i] =
1998 node->I.ivt.head_node[
1999 search_interval_tree_buffer->overlap_iv[i]].value;
2000 }
2001
2002 *num_overlap_cells += search_interval_tree_buffer->num_overlaps;
2003
2004 } else {
2005
2006 ENSURE_ARRAY_SIZE(*overlap_cells, *overlap_cells_array_size,
2007 *num_overlap_cells + node->I_size);
2008 memcpy(*overlap_cells + *num_overlap_cells, node->I.list,
2009 node->I_size * sizeof(**overlap_cells));
2010 *num_overlap_cells += node->I_size;
2011 }
2012 }
2013}
2014
2016 struct bnd_sphere_part_search * search, yac_coordinate_pointer coordinates_xyz,
2017 size_t count, size_t ** cells, size_t * num_cells_per_coordinate) {
2018
2019 struct sphere_part_node * base_node = &(search->base_node);
2020
2021 size_t * temp_search_results = NULL;
2022 size_t temp_search_results_array_size = 0;
2023 size_t num_temp_search_results = 0;
2024
2025 struct overlaps search_interval_tree_buffer = {0, 0, NULL};
2026
2027 memset(
2028 num_cells_per_coordinate, 0, count * sizeof(*num_cells_per_coordinate));
2029 size_t * cells_ = NULL;
2030 size_t cells_array_size = 0;
2031 size_t num_cells = 0;
2032 ENSURE_ARRAY_SIZE(cells_, cells_array_size, count);
2033
2034 for (size_t i = 0; i < count; ++i) {
2035
2036 double * curr_coordinates_xyz = &(coordinates_xyz[i][0]);
2037
2038 num_temp_search_results = 0;
2039
2040 double gc_norm_vector[3] = {0.0,0.0,1.0};
2041
2042 search_point(base_node, curr_coordinates_xyz, &temp_search_results,
2043 &temp_search_results_array_size, &num_temp_search_results,
2044 &search_interval_tree_buffer, gc_norm_vector);
2045
2047 cells_, cells_array_size, num_cells + num_temp_search_results);
2048
2049 memcpy(cells_ + num_cells, temp_search_results,
2050 num_temp_search_results * sizeof(*temp_search_results));
2051 num_cells_per_coordinate[i] = num_temp_search_results;
2052 num_cells += num_temp_search_results;
2053 }
2054
2055 free(temp_search_results);
2056 free(search_interval_tree_buffer.overlap_iv);
2057
2058 *cells = xrealloc(cells_, num_cells * sizeof(*cells_));
2059}
2060
2062 struct bnd_sphere_part_search * search, struct bounding_circle * bnd_circles,
2063 size_t count, size_t ** cells, size_t * num_cells_per_bnd_circle) {
2064
2065 struct sphere_part_node * base_node = &(search->base_node);
2066
2067 size_t * temp_search_results = NULL;
2068 size_t temp_search_results_array_size = 0;
2069 size_t num_temp_search_results = 0;
2070
2071 struct overlaps search_interval_tree_buffer = {0, 0, NULL};
2072
2073 memset(
2074 num_cells_per_bnd_circle, 0, count * sizeof(*num_cells_per_bnd_circle));
2075 size_t * cells_ = NULL;
2076 size_t cells_array_size = 0;
2077 size_t num_cells = 0;
2078 ENSURE_ARRAY_SIZE(cells_, cells_array_size, count);
2079
2080 for (size_t i = 0; i < count; ++i) {
2081
2082 num_temp_search_results = 0;
2083
2084 double gc_norm_vector[3] = {0.0,0.0,1.0};
2085
2086 search_bnd_circle(base_node, bnd_circles[i], &temp_search_results,
2087 &temp_search_results_array_size, &num_temp_search_results,
2088 &search_interval_tree_buffer, gc_norm_vector);
2089
2091 cells_, cells_array_size, num_cells + num_temp_search_results);
2092
2093 memcpy(cells_ + num_cells, temp_search_results,
2094 num_temp_search_results * sizeof(*temp_search_results));
2095 num_cells_per_bnd_circle[i] = num_temp_search_results;
2096 num_cells += num_temp_search_results;
2097 }
2098
2099 free(temp_search_results);
2100 free(search_interval_tree_buffer.overlap_iv);
2101
2102 *cells = xrealloc(cells_, num_cells * sizeof(*cells_));
2103}
2104
2105static void free_sphere_part_tree (struct sphere_part_node tree) {
2106
2107 // free I_list
2108 if (tree.flags & I_IS_INTERVAL_TREE)
2109 free(tree.I.ivt.head_node);
2110
2111 if ((tree.flags & U_IS_LEAF) == 0) {
2112 free_sphere_part_tree(*(struct sphere_part_node*)(tree.U));
2113 free(tree.U);
2114 }
2115
2116 if ((tree.flags & T_IS_LEAF) == 0) {
2117 free_sphere_part_tree(*(struct sphere_part_node*)(tree.T));
2118 free(tree.T);
2119 }
2120}
2121
2123
2124 if ((tree->flags & U_IS_LEAF) == 0) {
2126 free(tree->U);
2127 }
2128
2129 if ((tree->flags & T_IS_LEAF) == 0) {
2131 free(tree->T);
2132 }
2133}
2134
2136 struct point_sphere_part_search * search) {
2137
2138 if (search == NULL) return;
2139
2141 free(search->points);
2142 free(search);
2143}
2144
2146
2147 if (search == NULL) return;
2148
2150 free(search->ids);
2151 free(search);
2152}
#define YAC_ASSERT(exp, msg)
#define ENSURE_ARRAY_SIZE(arrayp, curr_array_size, req_size)
static struct sin_cos_angle get_vector_angle_2(double const a[3], double const b[3])
Definition geometry.h:376
static double clamp_abs_one(double val)
Definition geometry.h:208
static const struct sin_cos_angle SIN_COS_ZERO
Definition geometry.h:36
static int points_are_identically(double const *a, double const *b)
Definition geometry.h:609
static const struct sin_cos_angle SIN_COS_M_PI
Definition geometry.h:41
static struct sin_cos_angle sum_angles_no_check(struct sin_cos_angle a, struct sin_cos_angle b)
Definition geometry.h:485
static void crossproduct_kahan(double const a[], double const b[], double cross[])
Definition geometry.h:297
static const struct sin_cos_angle SIN_COS_M_PI_2
Definition geometry.h:40
static int compare_angles(struct sin_cos_angle a, struct sin_cos_angle b)
Definition geometry.h:389
#define yac_angle_tol
Definition geometry.h:26
static void normalise_vector(double v[])
Definition geometry.h:650
static double get_vector_angle(double const a[3], double const b[3])
Definition geometry.h:329
static struct sin_cos_angle sin_cos_angle_new(double sin, double cos)
Definition geometry.h:366
void yac_search_interval_tree(struct interval_node tree[], size_t num_nodes, struct interval query, struct overlaps *overlaps)
void yac_generate_interval_tree(struct interval_node intervals[], size_t num_nodes)
#define xrealloc(ptr, size)
Definition ppm_xfuncs.h:67
#define xmalloc(size)
Definition ppm_xfuncs.h:66
static int leaf_contains_matching_point(struct point_id_xyz *points, size_t num_points, double coordinate_xyz[3], size_t **local_point_ids, size_t *local_point_ids_array_size, double(**result_coordinates_xyz)[3], size_t *result_coordinates_xyz_array_size, size_t total_num_local_point_ids, size_t *num_local_point_ids)
static size_t initial_point_bnd_search_NNN(size_t n, struct point_id_xyz *points, size_t num_points, double *point_coordinates_xyz, struct point_id_xyz_angle **results, size_t *results_array_size)
static void search_point(struct sphere_part_node *node, double point[], size_t **overlap_cells, size_t *overlap_cells_array_size, size_t *num_overlap_cells, struct overlaps *search_interval_tree_buffer, double prev_gc_norm_vector[])
static int compare_point_idx_xyz(void const *a, void const *b)
void yac_point_sphere_part_search_NNN_ubound(struct point_sphere_part_search *search, size_t num_points, yac_coordinate_pointer coordinates_xyz, size_t n, struct sin_cos_angle *angles)
static size_t swap_node_type(struct temp_partition_data *part_data, size_t i, int node_type, size_t begin, size_t end)
static void check_leaf_NN(struct point_id_xyz *points, size_t num_points, double *point_coordinates_xyz, struct sin_cos_angle *best_angle, double(**result_coordinates_xyz)[3], size_t *result_coordinates_xyz_array_size, size_t **local_point_ids, size_t *local_point_ids_array_size, size_t total_num_local_point_ids, size_t *num_local_point_ids)
static struct point_id_xyz * get_unique_points(size_t *num_points, yac_const_coordinate_pointer coordinates_xyz, yac_int const *ids, int const *mask)
void yac_bnd_sphere_part_search_do_bnd_circle_search(struct bnd_sphere_part_search *search, struct bounding_circle *bnd_circles, size_t count, size_t **cells, size_t *num_cells_per_bnd_circle)
static void compute_gc_norm_vector(void *coords_data, size_t coords_size, size_t coords_count, double prev_gc_norm_vector[], double gc_norm_vector[])
static struct point_sphere_part_node * partition_point_data(struct point_id_xyz *points, size_t num_points, size_t threshold, double prev_gc_norm_vector[], size_t curr_tree_depth, size_t *max_tree_depth, int *list_flag)
static void point_search_NN(struct bounding_circle *bnd_circle, double(**result_coordinates_xyz)[3], size_t *result_coordinates_xyz_array_size, size_t **local_point_ids, size_t *local_point_ids_array_size, size_t total_num_local_point_ids, size_t *num_local_point_ids, double *dot_stack, struct point_sphere_part_node **node_stack, int *flags, size_t curr_tree_depth)
void yac_point_sphere_part_search_NNN_bnd_circle(struct point_sphere_part_search *search, size_t num_bnd_circles, struct bounding_circle *bnd_circles, size_t n, size_t **local_point_ids, size_t *local_point_ids_array_size, size_t *num_local_point_ids)
struct bnd_sphere_part_search * yac_bnd_sphere_part_search_new(struct bounding_circle *circles, size_t num_circles)
static void search_bnd_circle_I_node(struct sphere_part_node *node, struct bounding_circle bnd_circle, size_t **restrict overlap_cells, size_t *overlap_cells_array_size, size_t *restrict num_overlap_cells, struct overlaps *search_interval_tree_buffer, double prev_gc_norm_vector[])
node_type
Definition sphere_part.c:65
@ I_NODE
Definition sphere_part.c:67
@ T_NODE
Definition sphere_part.c:69
@ U_NODE
Definition sphere_part.c:68
@ I_NODE_FULL
Definition sphere_part.c:66
static void search_big_bnd_circle(struct sphere_part_node *node, struct bounding_circle bnd_circle, size_t **restrict overlap_cells, size_t *overlap_cells_array_size, size_t *restrict num_overlap_cells, struct overlaps *search_interval_tree_buffer, double prev_gc_norm_vector[])
TODO change to iterative implementation and allocate overlap_cells first.
static void sort_partition_data(struct temp_partition_data *part_data, size_t I_FULL_size, size_t I_size, size_t U_size, size_t T_size)
static void partition_data(size_t *local_cell_ids, struct temp_partition_data *part_data, size_t num_cell_ids, size_t threshold, struct sphere_part_node *parent_node, double prev_gc_norm_vector[])
void yac_delete_point_sphere_part_search(struct point_sphere_part_search *search)
struct point_sphere_part_search * yac_point_sphere_part_search_mask_new(size_t num_points, yac_const_coordinate_pointer coordinates_xyz, yac_int const *ids, int const *mask)
void yac_bnd_sphere_part_search_delete(struct bnd_sphere_part_search *search)
@ I_list_tree_min_size
Definition sphere_part.c:32
struct point_sphere_part_search * yac_point_sphere_part_search_new(size_t num_points, yac_const_coordinate_pointer coordinates_xyz, yac_int const *ids)
static struct sin_cos_angle check_leaf_NNN(size_t n, double *point_coordinates_xyz, struct point_id_xyz *points, size_t num_points, struct point_id_xyz_angle **results, size_t *results_array_size, size_t *num_results, struct sin_cos_angle curr_angle)
static int compare_angles_(void const *a, void const *b)
static struct sphere_part_node * get_sphere_part_node()
static void free_sphere_part_tree(struct sphere_part_node tree)
static int compare_point_id_xyz_angle(const void *a, const void *b)
void yac_bnd_sphere_part_search_do_point_search(struct bnd_sphere_part_search *search, yac_coordinate_pointer coordinates_xyz, size_t count, size_t **cells, size_t *num_cells_per_coordinate)
static void search_small_bnd_circle(struct sphere_part_node *node, struct bounding_circle bnd_circle, size_t **restrict overlap_cells, size_t *overlap_cells_array_size, size_t *restrict num_overlap_cells, struct overlaps *search_interval_tree_buffer, double prev_gc_norm_vector[])
static void swap_partition_data(struct temp_partition_data *a, struct temp_partition_data *b)
static void init_sphere_part_node(struct sphere_part_node *node)
yac_node_flags
Definition sphere_part.c:40
@ I_IS_INTERVAL_TREE
Definition sphere_part.c:43
@ T_IS_LEAF
Definition sphere_part.c:42
@ U_IS_LEAF
Definition sphere_part.c:41
static void search_bnd_circle(struct sphere_part_node *node, struct bounding_circle bnd_circle, size_t **restrict overlap_cells, size_t *overlap_cells_array_size, size_t *restrict num_overlap_cells, struct overlaps *search_interval_tree_buffer, double prev_gc_norm_vector[])
void yac_point_sphere_part_search_NN(struct point_sphere_part_search *search, size_t num_points, double(*coordinates_xyz)[3], double *cos_angles, double(**result_coordinates_xyz)[3], size_t *result_coordinates_xyz_array_size, size_t **local_point_ids, size_t *local_point_ids_array_size, size_t *num_local_point_ids)
void yac_point_sphere_part_search_NNN(struct point_sphere_part_search *search, size_t num_points, double(*coordinates_xyz)[3], size_t n, double **cos_angles, size_t *cos_angles_array_size, double(**result_coordinates_xyz)[3], size_t *result_coordinates_xyz_array_size, size_t **local_point_ids, size_t *local_point_ids_array_size, size_t *num_local_point_ids)
static void free_point_sphere_part_tree(struct point_sphere_part_node *tree)
@ U_FLAG
Definition sphere_part.c:36
@ T_FLAG
Definition sphere_part.c:37
static int compare_points_int32_coord(const void *a, const void *b)
static void point_search_NNN(size_t n, double *point_coordinates_xyz, struct point_id_xyz_angle **results, size_t *results_array_size, size_t *num_results, double *dot_stack, struct point_sphere_part_node **node_stack, int *flags, size_t curr_tree_depth)
algorithm for searching cells and points on a grid
struct sphere_part_node base_node
Definition sphere_part.c:61
struct sin_cos_angle inc_angle
angle between the middle point and the boundary of the spherical cap
Definition geometry.h:53
double base_vector[3]
Definition geometry.h:51
double sq_crd
Definition geometry.h:56
struct interval range
double right
double left
size_t * overlap_iv
size_t num_overlaps
struct point_id_xyz point
Definition sphere_part.c:85
int32_t coordinate_xyz[3]
double coordinates_xyz[3]
Definition sphere_part.c:79
struct point_sphere_part_node base_node
struct point_id_xyz * points
double sin
Definition geometry.h:33
double cos
Definition geometry.h:33
double gc_norm_vector[3]
Definition sphere_part.c:57
struct sin_cos_angle I_angle
Definition sphere_part.c:55
union I_list I
Definition sphere_part.c:50
struct bounding_circle bnd_circle
Definition sphere_part.c:73
size_t num_cells[2]
static int mask[16]
struct I_list::@62 ivt
struct interval_node * head_node
Definition sphere_part.c:25
size_t * list
Definition sphere_part.c:28
size_t num_nodes
Definition sphere_part.c:26
#define YAC_OMP_PARALLEL
#define YAC_OMP_FOR
static struct user_input_data_points ** points
Definition yac.c:157
static size_t num_points
Definition yac.c:158
Xt_int yac_int
Definition yac_types.h:15
double const (* yac_const_coordinate_pointer)[3]
Definition yac_types.h:20
double(* yac_coordinate_pointer)[3]
Definition yac_types.h:19