SPEED
kdtree2.f90
Go to the documentation of this file.
1!
2!(c) Matthew Kennel, Institute for Nonlinear Science (2004)
3!
4! Licensed under the Academic Free License version 1.1 found in file LICENSE
5! with additional provisions found in that same file.
6!
8
9 integer, parameter :: sp = kind(0.0)
10 integer, parameter :: dp = kind(0.0d0)
11
12 private :: sp, dp
13
14 !
15 ! You must comment out exactly one
16 ! of the two lines. If you comment
17 ! out kdkind = sp then you get single precision
18 ! and if you comment out kdkind = dp
19 ! you get double precision.
20 !
21
22 ! Single precision
23 !integer, parameter :: kdkind = sp
24
25 ! Double Precision
26 integer, parameter :: kdkind = dp
27
28 public :: kdkind
29
31
34 !
35 ! maintain a priority queue (PQ) of data, pairs of 'priority/payload',
36 ! implemented with a binary heap. This is the type, and the 'dis' field
37 ! is the priority.
38 !
40 ! a pair of distances, indexes
41 real(kdkind) :: dis!=0.0
42 integer :: idx!=-1 Initializers cause some bugs in compilers.
43 end type kdtree2_result
44 !
45 ! A heap-based priority queue lets one efficiently implement the following
46 ! operations, each in log(N) time, as opposed to linear time.
47 !
48 ! 1) add a datum (push a datum onto the queue, increasing its length)
49 ! 2) return the priority value of the maximum priority element
50 ! 3) pop-off (and delete) the element with the maximum priority, decreasing
51 ! the size of the queue.
52 ! 4) replace the datum with the maximum priority with a supplied datum
53 ! (of either higher or lower priority), maintaining the size of the
54 ! queue.
55 !
56 !
57 ! In the k-d tree case, the 'priority' is the square distance of a point in
58 ! the data set to a reference point. The goal is to keep the smallest M
59 ! distances to a reference point. The tree algorithm searches terminal
60 ! nodes to decide whether to add points under consideration.
61 !
62 ! A priority queue is useful here because it lets one quickly return the
63 ! largest distance currently existing in the list. If a new candidate
64 ! distance is smaller than this, then the new candidate ought to replace
65 ! the old candidate. In priority queue terms, this means removing the
66 ! highest priority element, and inserting the new one.
67 !
68 ! Algorithms based on Cormen, Leiserson, Rivest, _Introduction
69 ! to Algorithms_, 1990, with further optimization by the author.
70 !
71 ! Originally informed by a C implementation by Sriranga Veeraraghavan.
72 !
73 ! This module is not written in the most clear way, but is implemented such
74 ! for speed, as it its operations will be called many times during searches
75 ! of large numbers of neighbors.
76 !
77 type pq
78 !
79 ! The priority queue consists of elements
80 ! priority(1:heap_size), with associated payload(:).
81 !
82 ! There are heap_size active elements.
83 ! Assumes the allocation is always sufficient. Will NOT increase it
84 ! to match.
85 integer :: heap_size = 0
86 type(kdtree2_result), pointer :: elems(:)
87 end type pq
88
89 public :: kdtree2_result
90
91 public :: pq
92 public :: pq_create
93 public :: pq_delete, pq_insert
95 private
96
97contains
98
99
100 function pq_create(results_in) result(res)
101 !
102 ! Create a priority queue from ALREADY allocated
103 ! array pointers for storage. NOTE! It will NOT
104 ! add any alements to the heap, i.e. any existing
105 ! data in the input arrays will NOT be used and may
106 ! be overwritten.
107 !
108 ! usage:
109 ! real(kdkind), pointer :: x(:)
110 ! integer, pointer :: k(:)
111 ! allocate(x(1000),k(1000))
112 ! pq => pq_create(x,k)
113 !
114 type(kdtree2_result), target:: results_in(:)
115 type(pq) :: res
116 !
117 !
118 integer :: nalloc
119
120 nalloc = size(results_in,1)
121 if (nalloc .lt. 1) then
122 write (*,*) 'PQ_CREATE: error, input arrays must be allocated.'
123 end if
124 res%elems => results_in
125 res%heap_size = 0
126 return
127 end function pq_create
128
129 !
130 ! operations for getting parents and left + right children
131 ! of elements in a binary heap.
132 !
133
134!
135! These are written inline for speed.
136!
137! integer function parent(i)
138! integer, intent(in) :: i
139! parent = (i/2)
140! return
141! end function parent
142
143! integer function left(i)
144! integer, intent(in) ::i
145! left = (2*i)
146! return
147! end function left
148
149! integer function right(i)
150! integer, intent(in) :: i
151! right = (2*i)+1
152! return
153! end function right
154
155! logical function compare_priority(p1,p2)
156! real(kdkind), intent(in) :: p1, p2
157!
158! compare_priority = (p1 .gt. p2)
159! return
160! end function compare_priority
161
162 subroutine heapify(a,i_in)
163 !
164 ! take a heap rooted at 'i' and force it to be in the
165 ! heap canonical form. This is performance critical
166 ! and has been tweaked a little to reflect this.
167 !
168 type(pq),pointer :: a
169 integer, intent(in) :: i_in
170 !
171 integer :: i, l, r, largest
172
173 real(kdkind) :: pri_i, pri_l, pri_r, pri_largest
174
175
176 type(kdtree2_result) :: temp
177
178 i = i_in
179
180bigloop: do
181 l = 2*i ! left(i)
182 r = l+1 ! right(i)
183 !
184 ! set 'largest' to the index of either i, l, r
185 ! depending on whose priority is largest.
186 !
187 ! note that l or r can be larger than the heap size
188 ! in which case they do not count.
189
190
191 ! does left child have higher priority?
192 if (l .gt. a%heap_size) then
193 ! we know that i is the largest as both l and r are invalid.
194 exit
195 else
196 pri_i = a%elems(i)%dis
197 pri_l = a%elems(l)%dis
198 if (pri_l .gt. pri_i) then
199 largest = l
200 pri_largest = pri_l
201 else
202 largest = i
203 pri_largest = pri_i
204 endif
205
206 !
207 ! between i and l we have a winner
208 ! now choose between that and r.
209 !
210 if (r .le. a%heap_size) then
211 pri_r = a%elems(r)%dis
212 if (pri_r .gt. pri_largest) then
213 largest = r
214 endif
215 endif
216 endif
217
218 if (largest .ne. i) then
219 ! swap data in nodes largest and i, then heapify
220
221 temp = a%elems(i)
222 a%elems(i) = a%elems(largest)
223 a%elems(largest) = temp
224 !
225 ! Canonical heapify() algorithm has tail-ecursive call:
226 !
227 ! call heapify(a,largest)
228 ! we will simulate with cycle
229 !
230 i = largest
231 cycle bigloop ! continue the loop
232 else
233 return ! break from the loop
234 end if
235 enddo bigloop
236 return
237 end subroutine heapify
238
239 subroutine pq_max(a,e)
240 !
241 ! return the priority and its payload of the maximum priority element
242 ! on the queue, which should be the first one, if it is
243 ! in heapified form.
244 !
245 type(pq),pointer :: a
246 type(kdtree2_result),intent(out) :: e
247
248 if (a%heap_size .gt. 0) then
249 e = a%elems(1)
250 else
251 write (*,*) 'PQ_MAX: ERROR, heap_size < 1'
252 stop
253 endif
254 return
255 end subroutine pq_max
256
257 real(kdkind) function pq_maxpri(a)
258 type(pq), pointer :: a
259
260 if (a%heap_size .gt. 0) then
261 pq_maxpri = a%elems(1)%dis
262 else
263 write (*,*) 'PQ_MAX_PRI: ERROR, heapsize < 1'
264 stop
265 endif
266 return
267 end function pq_maxpri
268
269 subroutine pq_extract_max(a,e)
270 !
271 ! return the priority and payload of maximum priority
272 ! element, and remove it from the queue.
273 ! (equivalent to 'pop()' on a stack)
274 !
275 type(pq),pointer :: a
276 type(kdtree2_result), intent(out) :: e
277
278 if (a%heap_size .ge. 1) then
279 !
280 ! return max as first element
281 !
282 e = a%elems(1)
283
284 !
285 ! move last element to first
286 !
287 a%elems(1) = a%elems(a%heap_size)
288 a%heap_size = a%heap_size-1
289 call heapify(a,1)
290 return
291 else
292 write (*,*) 'PQ_EXTRACT_MAX: error, attempted to pop non-positive PQ'
293 stop
294 end if
295
296 end subroutine pq_extract_max
297
298
299 real(kdkind) function pq_insert(a,dis,idx)
300 !
301 ! Insert a new element and return the new maximum priority,
302 ! which may or may not be the same as the old maximum priority.
303 !
304 type(pq),pointer :: a
305 real(kdkind), intent(in) :: dis
306 integer, intent(in) :: idx
307 ! type(kdtree2_result), intent(in) :: e
308 !
309 integer :: i, isparent
310 real(kdkind) :: parentdis
311 !
312
313 ! if (a%heap_size .ge. a%max_elems) then
314 ! write (*,*) 'PQ_INSERT: error, attempt made to insert element on full PQ'
315 ! stop
316 ! else
317 a%heap_size = a%heap_size + 1
318 i = a%heap_size
319
320 do while (i .gt. 1)
321 isparent = int(i/2)
322 parentdis = a%elems(isparent)%dis
323 if (dis .gt. parentdis) then
324 ! move what was in i's parent into i.
325 a%elems(i)%dis = parentdis
326 a%elems(i)%idx = a%elems(isparent)%idx
327 i = isparent
328 else
329 exit
330 endif
331 end do
332
333 ! insert the element at the determined position
334 a%elems(i)%dis = dis
335 a%elems(i)%idx = idx
336
337 pq_insert = a%elems(1)%dis
338 return
339 ! end if
340
341 end function pq_insert
342
343 subroutine pq_adjust_heap(a,i)
344 type(pq),pointer :: a
345 integer, intent(in) :: i
346 !
347 ! nominally arguments (a,i), but specialize for a=1
348 !
349 ! This routine assumes that the trees with roots 2 and 3 are already heaps, i.e.
350 ! the children of '1' are heaps. When the procedure is completed, the
351 ! tree rooted at 1 is a heap.
352 real(kdkind) :: prichild
353 integer :: parent, child, n
354
355 type(kdtree2_result) :: e
356
357 e = a%elems(i)
358
359 parent = i
360 child = 2*i
361 n = a%heap_size
362
363 do while (child .le. n)
364 if (child .lt. n) then
365 if (a%elems(child)%dis .lt. a%elems(child+1)%dis) then
366 child = child+1
367 endif
368 endif
369 prichild = a%elems(child)%dis
370 if (e%dis .ge. prichild) then
371 exit
372 else
373 ! move child into parent.
374 a%elems(parent) = a%elems(child)
375 parent = child
376 child = 2*parent
377 end if
378 end do
379 a%elems(parent) = e
380 return
381 end subroutine pq_adjust_heap
382
383
384 real(kdkind) function pq_replace_max(a,dis,idx)
385 !
386 ! Replace the extant maximum priority element
387 ! in the PQ with (dis,idx). Return
388 ! the new maximum priority, which may be larger
389 ! or smaller than the old one.
390 !
391 type(pq),pointer :: a
392 real(kdkind), intent(in) :: dis
393 integer, intent(in) :: idx
394! type(kdtree2_result), intent(in) :: e
395 ! not tested as well!
396
397 integer :: parent, child, n
398 real(kdkind) :: prichild, prichildp1
399
400 type(kdtree2_result) :: etmp
401
402 if (.true.) then
403 n=a%heap_size
404 if (n .ge. 1) then
405 parent =1
406 child=2
407
408 loop: do while (child .le. n)
409 prichild = a%elems(child)%dis
410
411 !
412 ! posibly child+1 has higher priority, and if
413 ! so, get it, and increment child.
414 !
415
416 if (child .lt. n) then
417 prichildp1 = a%elems(child+1)%dis
418 if (prichild .lt. prichildp1) then
419 child = child+1
420 prichild = prichildp1
421 endif
422 endif
423
424 if (dis .ge. prichild) then
425 exit loop
426 ! we have a proper place for our new element,
427 ! bigger than either children's priority.
428 else
429 ! move child into parent.
430 a%elems(parent) = a%elems(child)
431 parent = child
432 child = 2*parent
433 end if
434 end do loop
435 a%elems(parent)%dis = dis
436 a%elems(parent)%idx = idx
437 pq_replace_max = a%elems(1)%dis
438 else
439 a%elems(1)%dis = dis
440 a%elems(1)%idx = idx
441 pq_replace_max = dis
442 endif
443 else
444 !
445 ! slower version using elementary pop and push operations.
446 !
447 call pq_extract_max(a,etmp)
448 etmp%dis = dis
449 etmp%idx = idx
450 pq_replace_max = pq_insert(a,dis,idx)
451 endif
452 return
453 end function pq_replace_max
454
455 subroutine pq_delete(a,i)
456 !
457 ! delete item with index 'i'
458 !
459 type(pq),pointer :: a
460 integer :: i
461
462 if ((i .lt. 1) .or. (i .gt. a%heap_size)) then
463 write (*,*) 'PQ_DELETE: error, attempt to remove out of bounds element.'
464 stop
465 endif
466
467 ! swap the item to be deleted with the last element
468 ! and shorten heap by one.
469 a%elems(i) = a%elems(a%heap_size)
470 a%heap_size = a%heap_size - 1
471
472 call heapify(a,i)
473
474 end subroutine pq_delete
475
477
478
482 ! K-D tree routines in Fortran 90 by Matt Kennel.
483 ! Original program was written in Sather by Steve Omohundro and
484 ! Matt Kennel. Only the Euclidean metric is supported.
485 !
486 !
487 ! This module is identical to 'kd_tree', except that the order
488 ! of subscripts is reversed in the data file.
489 ! In otherwords for an embedding of N D-dimensional vectors, the
490 ! data file is here, in natural Fortran order data(1:D, 1:N)
491 ! because Fortran lays out columns first,
492 !
493 ! whereas conventionally (C-style) it is data(1:N,1:D)
494 ! as in the original kd_tree module.
495 !
496 !-------------DATA TYPE, CREATION, DELETION---------------------
497 public :: kdkind
499 !---------------------------------------------------------------
500 !-------------------SEARCH ROUTINES-----------------------------
502 ! Return fixed number of nearest neighbors around arbitrary vector,
503 ! or extant point in dataset, with decorrelation window.
504 !
506 ! Return points within a fixed ball of arb vector/extant point
507 !
508 public :: kdtree2_sort_results
509 ! Sort, in order of increasing distance, rseults from above.
510 !
512 ! Count points within a fixed ball of arb vector/extant point
513 !
515 ! brute force of kdtree2_[n|r]_nearest
516 !----------------------------------------------------------------
517
518
519 integer, parameter :: bucket_size = 12
520 ! The maximum number of points to keep in a terminal node.
521
523 real(kdkind) :: lower,upper
524 end type interval
525
526 type :: tree_node
527 ! an internal tree node
528 private
529 integer :: cut_dim
530 ! the dimension to cut
531 real(kdkind) :: cut_val
532 ! where to cut the dimension
533 real(kdkind) :: cut_val_left, cut_val_right
534 ! improved cutoffs knowing the spread in child boxes.
535 integer :: l, u
536 type (tree_node), pointer :: left, right
537 type(interval), pointer :: box(:) => null()
538 ! child pointers
539 ! Points included in this node are indexes[k] with k \in [l,u]
540
541
542 end type tree_node
543
544 type :: kdtree2
545 ! Global information about the tree, one per tree
546 integer :: dimen=0, n=0
547 ! dimensionality and total # of points
548 real(kdkind), pointer :: the_data(:,:) => null()
549 ! pointer to the actual data array
550 !
551 ! IMPORTANT NOTE: IT IS DIMENSIONED the_data(1:d,1:N)
552 ! which may be opposite of what may be conventional.
553 ! This is, because in Fortran, the memory layout is such that
554 ! the first dimension is in sequential order. Hence, with
555 ! (1:d,1:N), all components of the vector will be in consecutive
556 ! memory locations. The search time is dominated by the
557 ! evaluation of distances in the terminal nodes. Putting all
558 ! vector components in consecutive memory location improves
559 ! memory cache locality, and hence search speed, and may enable
560 ! vectorization on some processors and compilers.
561
562 integer, pointer :: ind(:) => null()
563 ! permuted index into the data, so that indexes[l..u] of some
564 ! bucket represent the indexes of the actual points in that
565 ! bucket.
566 logical :: sort = .false.
567 ! do we always sort output results?
568 logical :: rearrange = .false.
569 real(kdkind), pointer :: rearranged_data(:,:) => null()
570 ! if (rearrange .eqv. .true.) then rearranged_data has been
571 ! created so that rearranged_data(:,i) = the_data(:,ind(i)),
572 ! permitting search to use more cache-friendly rearranged_data, at
573 ! some initial computation and storage cost.
574 type (tree_node), pointer :: root => null()
575 ! root pointer of the tree
576 end type kdtree2
577
578
579 type :: tree_search_record
580 !
581 ! One of these is created for each search.
582 !
583 private
584 !
585 ! Many fields are copied from the tree structure, in order to
586 ! speed up the search.
587 !
588 integer :: dimen
589 integer :: nn, nfound
590 real(kdkind) :: ballsize
591 integer :: centeridx=999, correltime=9999
592 ! exclude points within 'correltime' of 'centeridx', iff centeridx >= 0
593 integer :: nalloc ! how much allocated for results(:)?
594 logical :: rearrange ! are the data rearranged or original?
595 ! did the # of points found overflow the storage provided?
596 logical :: overflow
597 real(kdkind), pointer :: qv(:) ! query vector
598 type(kdtree2_result), pointer :: results(:) ! results
599 type(pq) :: pq
600 real(kdkind), pointer :: data(:,:) ! temp pointer to data
601 integer, pointer :: ind(:) ! temp pointer to indexes
602 end type tree_search_record
603
604 private
605 ! everything else is private.
606
607 type(tree_search_record), save, target :: sr ! A GLOBAL VARIABLE for search
608
609contains
610
611 function kdtree2_create(input_data,dim,sort,rearrange) result (mr)
612 !
613 ! create the actual tree structure, given an input array of data.
614 !
615 ! Note, input data is input_data(1:d,1:N), NOT the other way around.
616 ! THIS IS THE REVERSE OF THE PREVIOUS VERSION OF THIS MODULE.
617 ! The reason for it is cache friendliness, improving performance.
618 !
619 ! Optional arguments: If 'dim' is specified, then the tree
620 ! will only search the first 'dim' components
621 ! of input_data, otherwise, dim is inferred
622 ! from SIZE(input_data,1).
623 !
624 ! if sort .eqv. .true. then output results
625 ! will be sorted by increasing distance.
626 ! default=.false., as it is faster to not sort.
627 !
628 ! if rearrange .eqv. .true. then an internal
629 ! copy of the data, rearranged by terminal node,
630 ! will be made for cache friendliness.
631 ! default=.true., as it speeds searches, but
632 ! building takes longer, and extra memory is used.
633 !
634 ! .. Function Return Cut_value ..
635 type (kdtree2), pointer :: mr
636 integer, intent(in), optional :: dim
637 logical, intent(in), optional :: sort
638 logical, intent(in), optional :: rearrange
639 ! ..
640 ! .. Array Arguments ..
641 real(kdkind), target :: input_data(:,:)
642 !
643 integer :: i
644 ! ..
645 allocate (mr)
646 mr%the_data => input_data
647 ! pointer assignment
648
649 if (present(dim)) then
650 mr%dimen = dim
651 else
652 mr%dimen = size(input_data,1)
653 end if
654 mr%n = size(input_data,2)
655
656 if (mr%dimen > mr%n) then
657 ! unlikely to be correct
658 write (*,*) 'KD_TREE_TRANS: likely user error.'
659 write (*,*) 'KD_TREE_TRANS: You passed in matrix with D=',mr%dimen
660 write (*,*) 'KD_TREE_TRANS: and N=',mr%n
661 write (*,*) 'KD_TREE_TRANS: note, that new format is data(1:D,1:N)'
662 write (*,*) 'KD_TREE_TRANS: with usually N >> D. If N =approx= D, then a k-d tree'
663 write (*,*) 'KD_TREE_TRANS: is not an appropriate data structure.'
664 stop
665 end if
666
667 call build_tree(mr)
668
669 if (present(sort)) then
670 mr%sort = sort
671 else
672 mr%sort = .false.
673 endif
674
675 if (present(rearrange)) then
676 mr%rearrange = rearrange
677 else
678 mr%rearrange = .true.
679 endif
680
681 if (mr%rearrange) then
682 allocate(mr%rearranged_data(mr%dimen,mr%n))
683 do i=1,mr%n
684 mr%rearranged_data(:,i) = mr%the_data(:, &
685 mr%ind(i))
686 enddo
687 else
688 nullify(mr%rearranged_data)
689 endif
690
691 end function kdtree2_create
692
693 subroutine build_tree(tp)
694 type (kdtree2), pointer :: tp
695 ! ..
696 integer :: j
697 type(tree_node), pointer :: dummy => null()
698 ! ..
699 allocate (tp%ind(tp%n))
700 forall (j=1:tp%n)
701 tp%ind(j) = j
702 end forall
703 tp%root => build_tree_for_range(tp,1,tp%n, dummy)
704 end subroutine build_tree
705
706 recursive function build_tree_for_range(tp,l,u,parent) result (res)
707 ! .. Function Return Cut_value ..
708 type (tree_node), pointer :: res
709 ! ..
710 ! .. Structure Arguments ..
711 type (kdtree2), pointer :: tp
712 type (tree_node),pointer :: parent
713 ! ..
714 ! .. Scalar Arguments ..
715 integer, intent (In) :: l, u
716 ! ..
717 ! .. Local Scalars ..
718 integer :: i, c, m, dimen
719 logical :: recompute
720 real(kdkind) :: average
721
722!!$ If (.False.) Then
723!!$ If ((l .Lt. 1) .Or. (l .Gt. tp%n)) Then
724!!$ Stop 'illegal L value in build_tree_for_range'
725!!$ End If
726!!$ If ((u .Lt. 1) .Or. (u .Gt. tp%n)) Then
727!!$ Stop 'illegal u value in build_tree_for_range'
728!!$ End If
729!!$ If (u .Lt. l) Then
730!!$ Stop 'U is less than L, thats illegal.'
731!!$ End If
732!!$ Endif
733!!$
734 ! first compute min and max
735 dimen = tp%dimen
736 allocate (res)
737 allocate(res%box(dimen))
738
739 ! First, compute an APPROXIMATE bounding box of all points associated with this node.
740 if ( u < l ) then
741 ! no points in this box
742 nullify(res)
743 return
744 end if
745
746 if ((u-l)<=bucket_size) then
747 !
748 ! always compute true bounding box for terminal nodes.
749 !
750 do i=1,dimen
751 call spread_in_coordinate(tp,i,l,u,res%box(i))
752 end do
753 res%cut_dim = 0
754 res%cut_val = 0.0
755 res%l = l
756 res%u = u
757 res%left =>null()
758 res%right => null()
759 else
760 !
761 ! modify approximate bounding box. This will be an
762 ! overestimate of the true bounding box, as we are only recomputing
763 ! the bounding box for the dimension that the parent split on.
764 !
765 ! Going to a true bounding box computation would significantly
766 ! increase the time necessary to build the tree, and usually
767 ! has only a very small difference. This box is not used
768 ! for searching but only for deciding which coordinate to split on.
769 !
770 do i=1,dimen
771 recompute=.true.
772 if (associated(parent)) then
773 if (i .ne. parent%cut_dim) then
774 recompute=.false.
775 end if
776 endif
777 if (recompute) then
778 call spread_in_coordinate(tp,i,l,u,res%box(i))
779 else
780 res%box(i) = parent%box(i)
781 endif
782 end do
783
784
785 c = maxloc(res%box(1:dimen)%upper-res%box(1:dimen)%lower,1)
786 !
787 ! c is the identity of which coordinate has the greatest spread.
788 !
789
790 if (.false.) then
791 ! select exact median to have fully balanced tree.
792 m = (l+u)/2
793 call select_on_coordinate(tp%the_data,tp%ind,c,m,l,u)
794 else
795 !
796 ! select point halfway between min and max, as per A. Moore,
797 ! who says this helps in some degenerate cases, or
798 ! actual arithmetic average.
799 !
800 if (.true.) then
801 ! actually compute average
802 average = sum(tp%the_data(c,tp%ind(l:u))) / real(u-l+1,kdkind)
803 else
804 average = (res%box(c)%upper + res%box(c)%lower)/2.0
805 endif
806
807 res%cut_val = average
808 m = select_on_coordinate_value(tp%the_data,tp%ind,c,average,l,u)
809 endif
810
811 ! moves indexes around
812 res%cut_dim = c
813 res%l = l
814 res%u = u
815! res%cut_val = tp%the_data(c,tp%ind(m))
816
817 res%left => build_tree_for_range(tp,l,m,res)
818 res%right => build_tree_for_range(tp,m+1,u,res)
819
820 if (associated(res%right) .eqv. .false.) then
821 res%box = res%left%box
822 res%cut_val_left = res%left%box(c)%upper
823 res%cut_val = res%cut_val_left
824 elseif (associated(res%left) .eqv. .false.) then
825 res%box = res%right%box
826 res%cut_val_right = res%right%box(c)%lower
827 res%cut_val = res%cut_val_right
828 else
829 res%cut_val_right = res%right%box(c)%lower
830 res%cut_val_left = res%left%box(c)%upper
831 res%cut_val = (res%cut_val_left + res%cut_val_right)/2
832
833
834 ! now remake the true bounding box for self.
835 ! Since we are taking unions (in effect) of a tree structure,
836 ! this is much faster than doing an exhaustive
837 ! search over all points
838 res%box%upper = max(res%left%box%upper,res%right%box%upper)
839 res%box%lower = min(res%left%box%lower,res%right%box%lower)
840 endif
841 end if
842 end function build_tree_for_range
843
844 integer function select_on_coordinate_value(v,ind,c,alpha,li,ui) &
845 result(res)
846 ! Move elts of ind around between l and u, so that all points
847 ! <= than alpha (in c cooordinate) are first, and then
848 ! all points > alpha are second.
849
850 !
851 ! Algorithm (matt kennel).
852 !
853 ! Consider the list as having three parts: on the left,
854 ! the points known to be <= alpha. On the right, the points
855 ! known to be > alpha, and in the middle, the currently unknown
856 ! points. The algorithm is to scan the unknown points, starting
857 ! from the left, and swapping them so that they are added to
858 ! the left stack or the right stack, as appropriate.
859 !
860 ! The algorithm finishes when the unknown stack is empty.
861 !
862 ! .. Scalar Arguments ..
863 integer, intent (In) :: c, li, ui
864 real(kdkind), intent(in) :: alpha
865 ! ..
866 real(kdkind) :: v(1:,1:)
867 integer :: ind(1:)
868 integer :: tmp
869 ! ..
870 integer :: lb, rb
871 !
872 ! The points known to be <= alpha are in
873 ! [l,lb-1]
874 !
875 ! The points known to be > alpha are in
876 ! [rb+1,u].
877 !
878 ! Therefore we add new points into lb or
879 ! rb as appropriate. When lb=rb
880 ! we are done. We return the location of the last point <= alpha.
881 !
882 !
883 lb = li; rb = ui
884
885 do while (lb < rb)
886 if ( v(c,ind(lb)) <= alpha ) then
887 ! it is good where it is.
888 lb = lb+1
889 else
890 ! swap it with rb.
891 tmp = ind(lb); ind(lb) = ind(rb); ind(rb) = tmp
892 rb = rb-1
893 endif
894 end do
895
896 ! now lb .eq. ub
897 if (v(c,ind(lb)) <= alpha) then
898 res = lb
899 else
900 res = lb-1
901 endif
902
903 end function select_on_coordinate_value
904
905 subroutine select_on_coordinate(v,ind,c,k,li,ui)
906 ! Move elts of ind around between l and u, so that the kth
907 ! element
908 ! is >= those below, <= those above, in the coordinate c.
909 ! .. Scalar Arguments ..
910 integer, intent (In) :: c, k, li, ui
911 ! ..
912 integer :: i, l, m, s, t, u
913 ! ..
914 real(kdkind) :: v(:,:)
915 integer :: ind(:)
916 ! ..
917 l = li
918 u = ui
919 do while (l<u)
920 t = ind(l)
921 m = l
922 do i = l + 1, u
923 if (v(c,ind(i))<v(c,t)) then
924 m = m + 1
925 s = ind(m)
926 ind(m) = ind(i)
927 ind(i) = s
928 end if
929 end do
930 s = ind(l)
931 ind(l) = ind(m)
932 ind(m) = s
933 if (m<=k) l = m + 1
934 if (m>=k) u = m - 1
935 end do
936 end subroutine select_on_coordinate
937
938 subroutine spread_in_coordinate(tp,c,l,u,interv)
939 ! the spread in coordinate 'c', between l and u.
940 !
941 ! Return lower bound in 'smin', and upper in 'smax',
942 ! ..
943 ! .. Structure Arguments ..
944 type (kdtree2), pointer :: tp
945 type(interval), intent(out) :: interv
946 ! ..
947 ! .. Scalar Arguments ..
948 integer, intent (In) :: c, l, u
949 ! ..
950 ! .. Local Scalars ..
951 real(kdkind) :: last, lmax, lmin, t, smin,smax
952 integer :: i, ulocal
953 ! ..
954 ! .. Local Arrays ..
955 real(kdkind), pointer :: v(:,:)
956 integer, pointer :: ind(:)
957 ! ..
958 v => tp%the_data(1:,1:)
959 ind => tp%ind(1:)
960 smin = v(c,ind(l))
961 smax = smin
962
963 ulocal = u
964
965 do i = l + 2, ulocal, 2
966 lmin = v(c,ind(i-1))
967 lmax = v(c,ind(i))
968 if (lmin>lmax) then
969 t = lmin
970 lmin = lmax
971 lmax = t
972 end if
973 if (smin>lmin) smin = lmin
974 if (smax<lmax) smax = lmax
975 end do
976 if (i==ulocal+1) then
977 last = v(c,ind(ulocal))
978 if (smin>last) smin = last
979 if (smax<last) smax = last
980 end if
981
982 interv%lower = smin
983 interv%upper = smax
984
985 end subroutine spread_in_coordinate
986
987
988 subroutine kdtree2_destroy(tp)
989 ! Deallocates all memory for the tree, except input data matrix
990 ! .. Structure Arguments ..
991 type (kdtree2), pointer :: tp
992 ! ..
993 call destroy_node(tp%root)
994
995 deallocate (tp%ind)
996 nullify (tp%ind)
997
998 if (tp%rearrange) then
999 deallocate(tp%rearranged_data)
1000 nullify(tp%rearranged_data)
1001 endif
1002
1003 deallocate(tp)
1004 return
1005
1006 contains
1007 recursive subroutine destroy_node(np)
1008 ! .. Structure Arguments ..
1009 type (tree_node), pointer :: np
1010 ! ..
1011 ! .. Intrinsic Functions ..
1012 intrinsic ASSOCIATED
1013 ! ..
1014 if (associated(np%left)) then
1015 call destroy_node(np%left)
1016 nullify (np%left)
1017 end if
1018 if (associated(np%right)) then
1019 call destroy_node(np%right)
1020 nullify (np%right)
1021 end if
1022 if (associated(np%box)) deallocate(np%box)
1023 deallocate(np)
1024 return
1025
1026 end subroutine destroy_node
1027
1028 end subroutine kdtree2_destroy
1029
1030 subroutine kdtree2_n_nearest(tp,qv,nn,results)
1031 ! Find the 'nn' vectors in the tree nearest to 'qv' in euclidean norm
1032 ! returning their indexes and distances in 'indexes' and 'distances'
1033 ! arrays already allocated passed to this subroutine.
1034 type (kdtree2), pointer :: tp
1035 real(kdkind), target, intent (In) :: qv(:)
1036 integer, intent (In) :: nn
1037 type(kdtree2_result), target :: results(:)
1038
1039
1040 sr%ballsize = huge(1.0)
1041 sr%qv => qv
1042 sr%nn = nn
1043 sr%nfound = 0
1044 sr%centeridx = -1
1045 sr%correltime = 0
1046 sr%overflow = .false.
1047
1048 sr%results => results
1049
1050 sr%nalloc = nn ! will be checked
1051
1052 sr%ind => tp%ind
1053 sr%rearrange = tp%rearrange
1054 if (tp%rearrange) then
1055 sr%Data => tp%rearranged_data
1056 else
1057 sr%Data => tp%the_data
1058 endif
1059 sr%dimen = tp%dimen
1060
1061 call validate_query_storage(nn)
1062 sr%pq = pq_create(results)
1063
1064 call search(tp%root)
1065
1066 if (tp%sort) then
1067 call kdtree2_sort_results(nn, results)
1068 endif
1069! deallocate(sr%pqp)
1070 return
1071 end subroutine kdtree2_n_nearest
1072
1073 subroutine kdtree2_n_nearest_around_point(tp,idxin,correltime,nn,results)
1074 ! Find the 'nn' vectors in the tree nearest to point 'idxin',
1075 ! with correlation window 'correltime', returing results in
1076 ! results(:), which must be pre-allocated upon entry.
1077 type (kdtree2), pointer :: tp
1078 integer, intent (In) :: idxin, correltime, nn
1079 type(kdtree2_result), target :: results(:)
1080
1081 allocate (sr%qv(tp%dimen))
1082 sr%qv = tp%the_data(:,idxin) ! copy the vector
1083 sr%ballsize = huge(1.0) ! the largest real(kdkind) number
1084 sr%centeridx = idxin
1085 sr%correltime = correltime
1086
1087 sr%nn = nn
1088 sr%nfound = 0
1089
1090 sr%dimen = tp%dimen
1091 sr%nalloc = nn
1092
1093 sr%results => results
1094
1095 sr%ind => tp%ind
1096 sr%rearrange = tp%rearrange
1097
1098 if (sr%rearrange) then
1099 sr%Data => tp%rearranged_data
1100 else
1101 sr%Data => tp%the_data
1102 endif
1103
1104 call validate_query_storage(nn)
1105 sr%pq = pq_create(results)
1106
1107 call search(tp%root)
1108
1109 if (tp%sort) then
1110 call kdtree2_sort_results(nn, results)
1111 endif
1112 deallocate (sr%qv)
1113 return
1114 end subroutine kdtree2_n_nearest_around_point
1115
1116 subroutine kdtree2_r_nearest(tp,qv,r2,nfound,nalloc,results)
1117 ! find the nearest neighbors to point 'idxin', within SQUARED
1118 ! Euclidean distance 'r2'. Upon ENTRY, nalloc must be the
1119 ! size of memory allocated for results(1:nalloc). Upon
1120 ! EXIT, nfound is the number actually found within the ball.
1121 !
1122 ! Note that if nfound .gt. nalloc then more neighbors were found
1123 ! than there were storage to store. The resulting list is NOT
1124 ! the smallest ball inside norm r^2
1125 !
1126 ! Results are NOT sorted unless tree was created with sort option.
1127 type (kdtree2), pointer :: tp
1128 real(kdkind), target, intent (In) :: qv(:)
1129 real(kdkind), intent(in) :: r2
1130 integer, intent(out) :: nfound
1131 integer, intent (In) :: nalloc
1132 type(kdtree2_result), target :: results(:)
1133
1134 !
1135 sr%qv => qv
1136 sr%ballsize = r2
1137 sr%nn = 0 ! flag for fixed ball search
1138 sr%nfound = 0
1139 sr%centeridx = -1
1140 sr%correltime = 0
1141
1142 sr%results => results
1143
1144 call validate_query_storage(nalloc)
1145 sr%nalloc = nalloc
1146 sr%overflow = .false.
1147 sr%ind => tp%ind
1148 sr%rearrange= tp%rearrange
1149
1150 if (tp%rearrange) then
1151 sr%Data => tp%rearranged_data
1152 else
1153 sr%Data => tp%the_data
1154 endif
1155 sr%dimen = tp%dimen
1156
1157 !
1158 !sr%dsl = Huge(sr%dsl) ! set to huge positive values
1159 !sr%il = -1 ! set to invalid indexes
1160 !
1161
1162 call search(tp%root)
1163 nfound = sr%nfound
1164 if (tp%sort) then
1165 call kdtree2_sort_results(nfound, results)
1166 endif
1167
1168 if (sr%overflow) then
1169 write (*,*) 'KD_TREE_TRANS: warning! return from kdtree2_r_nearest found more neighbors'
1170 write (*,*) 'KD_TREE_TRANS: than storage was provided for. Answer is NOT smallest ball'
1171 write (*,*) 'KD_TREE_TRANS: with that number of neighbors! I.e. it is wrong.'
1172 endif
1173
1174 return
1175 end subroutine kdtree2_r_nearest
1176
1177 subroutine kdtree2_r_nearest_around_point(tp,idxin,correltime,r2,&
1178 nfound,nalloc,results)
1179 !
1180 ! Like kdtree2_r_nearest, but around a point 'idxin' already existing
1181 ! in the data set.
1182 !
1183 ! Results are NOT sorted unless tree was created with sort option.
1184 !
1185 type (kdtree2), pointer :: tp
1186 integer, intent (In) :: idxin, correltime, nalloc
1187 real(kdkind), intent(in) :: r2
1188 integer, intent(out) :: nfound
1189 type(kdtree2_result), target :: results(:)
1190 ! ..
1191 ! .. Intrinsic Functions ..
1192 intrinsic huge
1193 ! ..
1194 allocate (sr%qv(tp%dimen))
1195 sr%qv = tp%the_data(:,idxin) ! copy the vector
1196 sr%ballsize = r2
1197 sr%nn = 0 ! flag for fixed r search
1198 sr%nfound = 0
1199 sr%centeridx = idxin
1200 sr%correltime = correltime
1201
1202 sr%results => results
1203
1204 sr%nalloc = nalloc
1205 sr%overflow = .false.
1206
1207 call validate_query_storage(nalloc)
1208
1209 ! sr%dsl = HUGE(sr%dsl) ! set to huge positive values
1210 ! sr%il = -1 ! set to invalid indexes
1211
1212 sr%ind => tp%ind
1213 sr%rearrange = tp%rearrange
1214
1215 if (tp%rearrange) then
1216 sr%Data => tp%rearranged_data
1217 else
1218 sr%Data => tp%the_data
1219 endif
1220 sr%rearrange = tp%rearrange
1221 sr%dimen = tp%dimen
1222
1223 !
1224 !sr%dsl = Huge(sr%dsl) ! set to huge positive values
1225 !sr%il = -1 ! set to invalid indexes
1226 !
1227
1228 call search(tp%root)
1229 nfound = sr%nfound
1230 if (tp%sort) then
1231 call kdtree2_sort_results(nfound,results)
1232 endif
1233
1234 if (sr%overflow) then
1235 write (*,*) 'KD_TREE_TRANS: warning! return from kdtree2_r_nearest found more neighbors'
1236 write (*,*) 'KD_TREE_TRANS: than storage was provided for. Answer is NOT smallest ball'
1237 write (*,*) 'KD_TREE_TRANS: with that number of neighbors! I.e. it is wrong.'
1238 endif
1239
1240 deallocate (sr%qv)
1241 return
1242 end subroutine kdtree2_r_nearest_around_point
1243
1244 function kdtree2_r_count(tp,qv,r2) result(nfound)
1245 ! Count the number of neighbors within square distance 'r2'.
1246 type (kdtree2), pointer :: tp
1247 real(kdkind), target, intent (In) :: qv(:)
1248 real(kdkind), intent(in) :: r2
1249 integer :: nfound
1250 ! ..
1251 ! .. Intrinsic Functions ..
1252 intrinsic huge
1253 ! ..
1254 sr%qv => qv
1255 sr%ballsize = r2
1256
1257 sr%nn = 0 ! flag for fixed r search
1258 sr%nfound = 0
1259 sr%centeridx = -1
1260 sr%correltime = 0
1261
1262 nullify(sr%results) ! for some reason, FTN 95 chokes on '=> null()'
1263
1264 sr%nalloc = 0 ! we do not allocate any storage but that's OK
1265 ! for counting.
1266 sr%ind => tp%ind
1267 sr%rearrange = tp%rearrange
1268 if (tp%rearrange) then
1269 sr%Data => tp%rearranged_data
1270 else
1271 sr%Data => tp%the_data
1272 endif
1273 sr%dimen = tp%dimen
1274
1275 !
1276 !sr%dsl = Huge(sr%dsl) ! set to huge positive values
1277 !sr%il = -1 ! set to invalid indexes
1278 !
1279 sr%overflow = .false.
1280
1281 call search(tp%root)
1282
1283 nfound = sr%nfound
1284
1285 return
1286 end function kdtree2_r_count
1287
1288 function kdtree2_r_count_around_point(tp,idxin,correltime,r2) &
1289 result(nfound)
1290 ! Count the number of neighbors within square distance 'r2' around
1291 ! point 'idxin' with decorrelation time 'correltime'.
1292 !
1293 type (kdtree2), pointer :: tp
1294 integer, intent (In) :: correltime, idxin
1295 real(kdkind), intent(in) :: r2
1296 integer :: nfound
1297 ! ..
1298 ! ..
1299 ! .. Intrinsic Functions ..
1300 intrinsic huge
1301 ! ..
1302 allocate (sr%qv(tp%dimen))
1303 sr%qv = tp%the_data(:,idxin)
1304 sr%ballsize = r2
1305
1306 sr%nn = 0 ! flag for fixed r search
1307 sr%nfound = 0
1308 sr%centeridx = idxin
1309 sr%correltime = correltime
1310 nullify(sr%results)
1311
1312 sr%nalloc = 0 ! we do not allocate any storage but that's OK
1313 ! for counting.
1314
1315 sr%ind => tp%ind
1316 sr%rearrange = tp%rearrange
1317
1318 if (sr%rearrange) then
1319 sr%Data => tp%rearranged_data
1320 else
1321 sr%Data => tp%the_data
1322 endif
1323 sr%dimen = tp%dimen
1324
1325 !
1326 !sr%dsl = Huge(sr%dsl) ! set to huge positive values
1327 !sr%il = -1 ! set to invalid indexes
1328 !
1329 sr%overflow = .false.
1330
1331 call search(tp%root)
1332
1333 nfound = sr%nfound
1334
1335 return
1336 end function kdtree2_r_count_around_point
1337
1338
1339 subroutine validate_query_storage(n)
1340 !
1341 ! make sure we have enough storage for n
1342 !
1343 integer, intent(in) :: n
1344
1345 if (size(sr%results,1) .lt. n) then
1346 write (*,*) 'KD_TREE_TRANS: you did not provide enough storage for results(1:n)'
1347 stop
1348 return
1349 endif
1350
1351 return
1352 end subroutine validate_query_storage
1353
1354 function square_distance(d, iv,qv) result (res)
1355 ! distance between iv[1:n] and qv[1:n]
1356 ! .. Function Return Value ..
1357 ! re-implemented to improve vectorization.
1358 real(kdkind) :: res
1359 ! ..
1360 ! ..
1361 ! .. Scalar Arguments ..
1362 integer :: d
1363 ! ..
1364 ! .. Array Arguments ..
1365 real(kdkind) :: iv(:),qv(:)
1366 ! ..
1367 ! ..
1368 res = sum( (iv(1:d)-qv(1:d))**2 )
1369 end function square_distance
1370
1371 recursive subroutine search(node)
1372 !
1373 ! This is the innermost core routine of the kd-tree search. Along
1374 ! with "process_terminal_node", it is the performance bottleneck.
1375 !
1376 ! This version uses a logically complete secondary search of
1377 ! "box in bounds", whether the sear
1378 !
1379 type (tree_node), pointer :: node
1380 ! ..
1381 type(tree_node),pointer :: ncloser, nfarther
1382 !
1383 integer :: cut_dim, i
1384 ! ..
1385 real(kdkind) :: qval, dis
1386 real(kdkind) :: ballsize
1387 real(kdkind), pointer :: qv(:)
1388 type(interval), pointer :: box(:)
1389
1390 if ((associated(node%left) .and. associated(node%right)) .eqv. .false.) then
1391 ! we are on a terminal node
1392 if (sr%nn .eq. 0) then
1393 call process_terminal_node_fixedball(node)
1394 else
1395 call process_terminal_node(node)
1396 endif
1397 else
1398 ! we are not on a terminal node
1399 qv => sr%qv(1:)
1400 cut_dim = node%cut_dim
1401 qval = qv(cut_dim)
1402
1403 if (qval < node%cut_val) then
1404 ncloser => node%left
1405 nfarther => node%right
1406 dis = (node%cut_val_right - qval)**2
1407! extra = node%cut_val - qval
1408 else
1409 ncloser => node%right
1410 nfarther => node%left
1411 dis = (node%cut_val_left - qval)**2
1412! extra = qval- node%cut_val_left
1413 endif
1414
1415 if (associated(ncloser)) call search(ncloser)
1416
1417 ! we may need to search the second node.
1418 if (associated(nfarther)) then
1419 ballsize = sr%ballsize
1420! dis=extra**2
1421 if (dis <= ballsize) then
1422 !
1423 ! we do this separately as going on the first cut dimen is often
1424 ! a good idea.
1425 ! note that if extra**2 < sr%ballsize, then the next
1426 ! check will also be false.
1427 !
1428 box => node%box(1:)
1429 do i=1,sr%dimen
1430 if (i .ne. cut_dim) then
1431 dis = dis + dis2_from_bnd(qv(i),box(i)%lower,box(i)%upper)
1432 if (dis > ballsize) then
1433 return
1434 endif
1435 endif
1436 end do
1437
1438 !
1439 ! if we are still here then we need to search mroe.
1440 !
1441 call search(nfarther)
1442 endif
1443 endif
1444 end if
1445 end subroutine search
1446
1447
1448 real(kdkind) function dis2_from_bnd(x,amin,amax) result (res)
1449 real(kdkind), intent(in) :: x, amin,amax
1450
1451 if (x > amax) then
1452 res = (x-amax)**2;
1453 return
1454 else
1455 if (x < amin) then
1456 res = (amin-x)**2;
1457 return
1458 else
1459 res = 0.0
1460 return
1461 endif
1462 endif
1463 return
1464 end function dis2_from_bnd
1465
1466 logical function box_in_search_range(node, sr) result(res)
1467 !
1468 ! Return the distance from 'qv' to the CLOSEST corner of node's
1469 ! bounding box
1470 ! for all coordinates outside the box. Coordinates inside the box
1471 ! contribute nothing to the distance.
1472 !
1473 type (tree_node), pointer :: node
1474 type (tree_search_record), pointer :: sr
1475
1476 integer :: dimen, i
1477 real(kdkind) :: dis, ballsize
1478 real(kdkind) :: l, u
1479
1480 dimen = sr%dimen
1481 ballsize = sr%ballsize
1482 dis = 0.0
1483 res = .true.
1484 do i=1,dimen
1485 l = node%box(i)%lower
1486 u = node%box(i)%upper
1487 dis = dis + (dis2_from_bnd(sr%qv(i),l,u))
1488 if (dis > ballsize) then
1489 res = .false.
1490 return
1491 endif
1492 end do
1493 res = .true.
1494 return
1495 end function box_in_search_range
1496
1497
1498 subroutine process_terminal_node(node)
1499 !
1500 ! Look for actual near neighbors in 'node', and update
1501 ! the search results on the sr data structure.
1502 !
1503 type (tree_node), pointer :: node
1504 !
1505 real(kdkind), pointer :: qv(:)
1506 integer, pointer :: ind(:)
1507 real(kdkind), pointer :: data(:,:)
1508 !
1509 integer :: dimen, i, indexofi, k, centeridx, correltime
1510 real(kdkind) :: ballsize, sd, newpri
1511 logical :: rearrange
1512 type(pq), pointer :: pqp
1513 !
1514 ! copy values from sr to local variables
1515 !
1516 !
1517 ! Notice, making local pointers with an EXPLICIT lower bound
1518 ! seems to generate faster code.
1519 ! why? I don't know.
1520 qv => sr%qv(1:)
1521 pqp => sr%pq
1522 dimen = sr%dimen
1523 ballsize = sr%ballsize
1524 rearrange = sr%rearrange
1525 ind => sr%ind(1:)
1526 data => sr%Data(1:,1:)
1527 centeridx = sr%centeridx
1528 correltime = sr%correltime
1529
1530 ! doing_correl = (centeridx >= 0) ! Do we have a decorrelation window?
1531 ! include_point = .true. ! by default include all points
1532 ! search through terminal bucket.
1533
1534 mainloop: do i = node%l, node%u
1535 if (rearrange) then
1536 sd = 0.0
1537 do k = 1,dimen
1538 sd = sd + (data(k,i) - qv(k))**2
1539 if (sd>ballsize) cycle mainloop
1540 end do
1541 indexofi = ind(i) ! only read it if we have not broken out
1542 else
1543 indexofi = ind(i)
1544 sd = 0.0
1545 do k = 1,dimen
1546 sd = sd + (data(k,indexofi) - qv(k))**2
1547 if (sd>ballsize) cycle mainloop
1548 end do
1549 endif
1550
1551 if (centeridx > 0) then ! doing correlation interval?
1552 if (abs(indexofi-centeridx) < correltime) cycle mainloop
1553 endif
1554
1555
1556 !
1557 ! two choices for any point. The list so far is either undersized,
1558 ! or it is not.
1559 !
1560 ! If it is undersized, then add the point and its distance
1561 ! unconditionally. If the point added fills up the working
1562 ! list then set the sr%ballsize, maximum distance bound (largest distance on
1563 ! list) to be that distance, instead of the initialized +infinity.
1564 !
1565 ! If the running list is full size, then compute the
1566 ! distance but break out immediately if it is larger
1567 ! than sr%ballsize, "best squared distance" (of the largest element),
1568 ! as it cannot be a good neighbor.
1569 !
1570 ! Once computed, compare to best_square distance.
1571 ! if it is smaller, then delete the previous largest
1572 ! element and add the new one.
1573
1574 if (sr%nfound .lt. sr%nn) then
1575 !
1576 ! add this point unconditionally to fill list.
1577 !
1578 sr%nfound = sr%nfound +1
1579 newpri = pq_insert(pqp,sd,indexofi)
1580 if (sr%nfound .eq. sr%nn) ballsize = newpri
1581 ! we have just filled the working list.
1582 ! put the best square distance to the maximum value
1583 ! on the list, which is extractable from the PQ.
1584
1585 else
1586 !
1587 ! now, if we get here,
1588 ! we know that the current node has a squared
1589 ! distance smaller than the largest one on the list, and
1590 ! belongs on the list.
1591 ! Hence we replace that with the current one.
1592 !
1593 ballsize = pq_replace_max(pqp,sd,indexofi)
1594 endif
1595 end do mainloop
1596 !
1597 ! Reset sr variables which may have changed during loop
1598 !
1599 sr%ballsize = ballsize
1600
1601 end subroutine process_terminal_node
1602
1603 subroutine process_terminal_node_fixedball(node)
1604 !
1605 ! Look for actual near neighbors in 'node', and update
1606 ! the search results on the sr data structure, i.e.
1607 ! save all within a fixed ball.
1608 !
1609 type (tree_node), pointer :: node
1610 !
1611 real(kdkind), pointer :: qv(:)
1612 integer, pointer :: ind(:)
1613 real(kdkind), pointer :: data(:,:)
1614 !
1615 integer :: nfound
1616 integer :: dimen, i, indexofi, k
1617 integer :: centeridx, correltime, nn
1618 real(kdkind) :: ballsize, sd
1619 logical :: rearrange
1620
1621 !
1622 ! copy values from sr to local variables
1623 !
1624 qv => sr%qv(1:)
1625 dimen = sr%dimen
1626 ballsize = sr%ballsize
1627 rearrange = sr%rearrange
1628 ind => sr%ind(1:)
1629 data => sr%Data(1:,1:)
1630 centeridx = sr%centeridx
1631 correltime = sr%correltime
1632 nn = sr%nn ! number to search for
1633 nfound = sr%nfound
1634
1635 ! search through terminal bucket.
1636 mainloop: do i = node%l, node%u
1637
1638 !
1639 ! two choices for any point. The list so far is either undersized,
1640 ! or it is not.
1641 !
1642 ! If it is undersized, then add the point and its distance
1643 ! unconditionally. If the point added fills up the working
1644 ! list then set the sr%ballsize, maximum distance bound (largest distance on
1645 ! list) to be that distance, instead of the initialized +infinity.
1646 !
1647 ! If the running list is full size, then compute the
1648 ! distance but break out immediately if it is larger
1649 ! than sr%ballsize, "best squared distance" (of the largest element),
1650 ! as it cannot be a good neighbor.
1651 !
1652 ! Once computed, compare to best_square distance.
1653 ! if it is smaller, then delete the previous largest
1654 ! element and add the new one.
1655
1656 ! which index to the point do we use?
1657
1658 if (rearrange) then
1659 sd = 0.0
1660 do k = 1,dimen
1661 sd = sd + (data(k,i) - qv(k))**2
1662 if (sd>ballsize) cycle mainloop
1663 end do
1664 indexofi = ind(i) ! only read it if we have not broken out
1665 else
1666 indexofi = ind(i)
1667 sd = 0.0
1668 do k = 1,dimen
1669 sd = sd + (data(k,indexofi) - qv(k))**2
1670 if (sd>ballsize) cycle mainloop
1671 end do
1672 endif
1673
1674 if (centeridx > 0) then ! doing correlation interval?
1675 if (abs(indexofi-centeridx)<correltime) cycle mainloop
1676 endif
1677
1678 nfound = nfound+1
1679 if (nfound .gt. sr%nalloc) then
1680 ! oh nuts, we have to add another one to the tree but
1681 ! there isn't enough room.
1682 sr%overflow = .true.
1683 else
1684 sr%results(nfound)%dis = sd
1685 sr%results(nfound)%idx = indexofi
1686 endif
1687 end do mainloop
1688 !
1689 ! Reset sr variables which may have changed during loop
1690 !
1691 sr%nfound = nfound
1692 end subroutine process_terminal_node_fixedball
1693
1694 subroutine kdtree2_n_nearest_brute_force(tp,qv,nn,results)
1695 ! find the 'n' nearest neighbors to 'qv' by exhaustive search.
1696 ! only use this subroutine for testing, as it is SLOW! The
1697 ! whole point of a k-d tree is to avoid doing what this subroutine
1698 ! does.
1699 type (kdtree2), pointer :: tp
1700 real(kdkind), intent (In) :: qv(:)
1701 integer, intent (In) :: nn
1702 type(kdtree2_result) :: results(:)
1703
1704 integer :: i, j, k
1705 real(kdkind), allocatable :: all_distances(:)
1706 ! ..
1707 allocate (all_distances(tp%n))
1708 do i = 1, tp%n
1709 all_distances(i) = square_distance(tp%dimen,qv,tp%the_data(:,i))
1710 end do
1711 ! now find 'n' smallest distances
1712 do i = 1, nn
1713 results(i)%dis = huge(1.0)
1714 results(i)%idx = -1
1715 end do
1716 do i = 1, tp%n
1717 if (all_distances(i)<results(nn)%dis) then
1718 ! insert it somewhere on the list
1719 do j = 1, nn
1720 if (all_distances(i)<results(j)%dis) exit
1721 end do
1722 ! now we know 'j'
1723 do k = nn - 1, j, -1
1724 results(k+1) = results(k)
1725 end do
1726 results(j)%dis = all_distances(i)
1727 results(j)%idx = i
1728 end if
1729 end do
1730 deallocate (all_distances)
1731 end subroutine kdtree2_n_nearest_brute_force
1732
1733
1734 subroutine kdtree2_r_nearest_brute_force(tp,qv,r2,nfound,results)
1735 ! find the nearest neighbors to 'qv' with distance**2 <= r2 by exhaustive search.
1736 ! only use this subroutine for testing, as it is SLOW! The
1737 ! whole point of a k-d tree is to avoid doing what this subroutine
1738 ! does.
1739 type (kdtree2), pointer :: tp
1740 real(kdkind), intent (In) :: qv(:)
1741 real(kdkind), intent (In) :: r2
1742 integer, intent(out) :: nfound
1743 type(kdtree2_result) :: results(:)
1744
1745 integer :: i, nalloc
1746 real(kdkind), allocatable :: all_distances(:)
1747 ! ..
1748 allocate (all_distances(tp%n))
1749 do i = 1, tp%n
1750 all_distances(i) = square_distance(tp%dimen,qv,tp%the_data(:,i))
1751 end do
1752
1753 nfound = 0
1754 nalloc = size(results,1)
1755
1756 do i = 1, tp%n
1757 if (all_distances(i)< r2) then
1758 ! insert it somewhere on the list
1759 if (nfound .lt. nalloc) then
1760 nfound = nfound+1
1761 results(nfound)%dis = all_distances(i)
1762 results(nfound)%idx = i
1763 endif
1764 end if
1765 enddo
1766 deallocate (all_distances)
1767
1768 call kdtree2_sort_results(nfound,results)
1769
1770
1771 end subroutine kdtree2_r_nearest_brute_force
1772
1773 subroutine kdtree2_sort_results(nfound,results)
1774 ! Use after search to sort results(1:nfound) in order of increasing
1775 ! distance.
1776 integer, intent(in) :: nfound
1777 type(kdtree2_result), target :: results(:)
1778 !
1779 !
1780
1781 !THIS IS BUGGY WITH INTEL FORTRAN
1782 ! If (nfound .Gt. 1) Call heapsort(results(1:nfound)%dis,results(1:nfound)%ind,nfound)
1783 !
1784 if (nfound .gt. 1) call heapsort_struct(results,nfound)
1785
1786 return
1787 end subroutine kdtree2_sort_results
1788
1789 subroutine heapsort(a,ind,n)
1790 !
1791 ! Sort a(1:n) in ascending order, permuting ind(1:n) similarly.
1792 !
1793 ! If ind(k) = k upon input, then it will give a sort index upon output.
1794 !
1795 integer,intent(in) :: n
1796 real(kdkind), intent(inout) :: a(:)
1797 integer, intent(inout) :: ind(:)
1798
1799 !
1800 !
1801 real(kdkind) :: value ! temporary for a value from a()
1802 integer :: ivalue ! temporary for a value from ind()
1803
1804 integer :: i,j
1805 integer :: ileft,iright
1806
1807 ileft=n/2+1
1808 iright=n
1809
1810 ! do i=1,n
1811 ! ind(i)=i
1812 ! Generate initial idum array
1813 ! end do
1814
1815 if(n.eq.1) return
1816
1817 do
1818 if(ileft > 1)then
1819 ileft=ileft-1
1820 value=a(ileft); ivalue=ind(ileft)
1821 else
1822 value=a(iright); ivalue=ind(iright)
1823 a(iright)=a(1); ind(iright)=ind(1)
1824 iright=iright-1
1825 if (iright == 1) then
1826 a(1)=value;ind(1)=ivalue
1827 return
1828 endif
1829 endif
1830 i=ileft
1831 j=2*ileft
1832 do while (j <= iright)
1833 if(j < iright) then
1834 if(a(j) < a(j+1)) j=j+1
1835 endif
1836 if(value < a(j)) then
1837 a(i)=a(j); ind(i)=ind(j)
1838 i=j
1839 j=j+j
1840 else
1841 j=iright+1
1842 endif
1843 end do
1844 a(i)=value; ind(i)=ivalue
1845 end do
1846 end subroutine heapsort
1847
1848 subroutine heapsort_struct(a,n)
1849 !
1850 ! Sort a(1:n) in ascending order
1851 !
1852 !
1853 integer,intent(in) :: n
1854 type(kdtree2_result),intent(inout) :: a(:)
1855
1856 !
1857 !
1858 type(kdtree2_result) :: value ! temporary value
1859
1860 integer :: i,j
1861 integer :: ileft,iright
1862
1863 ileft=n/2+1
1864 iright=n
1865
1866 ! do i=1,n
1867 ! ind(i)=i
1868 ! Generate initial idum array
1869 ! end do
1870
1871 if(n.eq.1) return
1872
1873 do
1874 if(ileft > 1)then
1875 ileft=ileft-1
1876 value=a(ileft)
1877 else
1878 value=a(iright)
1879 a(iright)=a(1)
1880 iright=iright-1
1881 if (iright == 1) then
1882 a(1) = value
1883 return
1884 endif
1885 endif
1886 i=ileft
1887 j=2*ileft
1888 do while (j <= iright)
1889 if(j < iright) then
1890 if(a(j)%dis < a(j+1)%dis) j=j+1
1891 endif
1892 if(value%dis < a(j)%dis) then
1893 a(i)=a(j);
1894 i=j
1895 j=j+j
1896 else
1897 j=iright+1
1898 endif
1899 end do
1900 a(i)=value
1901 end do
1902 end subroutine heapsort_struct
1903
1904end module kdtree2_module
1905
type(kdtree2) function, pointer, public kdtree2_create(input_data, dim, sort, rearrange)
Definition kdtree2.f90:612
subroutine, public kdtree2_r_nearest_brute_force(tp, qv, r2, nfound, results)
Definition kdtree2.f90:1735
subroutine, public kdtree2_sort_results(nfound, results)
Definition kdtree2.f90:1774
subroutine, public kdtree2_n_nearest_brute_force(tp, qv, nn, results)
Definition kdtree2.f90:1695
integer function, public kdtree2_r_count_around_point(tp, idxin, correltime, r2)
Definition kdtree2.f90:1290
subroutine, public kdtree2_n_nearest_around_point(tp, idxin, correltime, nn, results)
Definition kdtree2.f90:1074
integer, parameter bucket_size
Definition kdtree2.f90:519
subroutine, public kdtree2_destroy(tp)
Definition kdtree2.f90:989
subroutine, public kdtree2_n_nearest(tp, qv, nn, results)
Definition kdtree2.f90:1031
subroutine, public kdtree2_r_nearest(tp, qv, r2, nfound, nalloc, results)
Definition kdtree2.f90:1117
subroutine, public kdtree2_r_nearest_around_point(tp, idxin, correltime, r2, nfound, nalloc, results)
Definition kdtree2.f90:1179
integer function, public kdtree2_r_count(tp, qv, r2)
Definition kdtree2.f90:1245
integer, parameter, public kdkind
Definition kdtree2.f90:26
subroutine, public pq_extract_max(a, e)
Definition kdtree2.f90:270
subroutine, public pq_delete(a, i)
Definition kdtree2.f90:456
real(kdkind) function, public pq_replace_max(a, dis, idx)
Definition kdtree2.f90:385
real(kdkind) function, public pq_maxpri(a)
Definition kdtree2.f90:258
subroutine, public pq_max(a, e)
Definition kdtree2.f90:240
type(pq) function, public pq_create(results_in)
Definition kdtree2.f90:101
real(kdkind) function, public pq_insert(a, dis, idx)
Definition kdtree2.f90:300