9 integer,
parameter :: sp = kind(0.0)
10 integer,
parameter :: dp = kind(0.0d0)
85 integer :: heap_size = 0
120 nalloc =
size(results_in,1)
121 if (nalloc .lt. 1)
then
122 write (*,*)
'PQ_CREATE: error, input arrays must be allocated.'
124 res%elems => results_in
162 subroutine heapify(a,i_in)
168 type(
pq),
pointer :: a
169 integer,
intent(in) :: i_in
171 integer :: i, l, r, largest
173 real(
kdkind) :: pri_i, pri_l, pri_r, pri_largest
192 if (l .gt. a%heap_size)
then
196 pri_i = a%elems(i)%dis
197 pri_l = a%elems(l)%dis
198 if (pri_l .gt. pri_i)
then
210 if (r .le. a%heap_size)
then
211 pri_r = a%elems(r)%dis
212 if (pri_r .gt. pri_largest)
then
218 if (largest .ne. i)
then
222 a%elems(i) = a%elems(largest)
223 a%elems(largest) = temp
237 end subroutine heapify
245 type(
pq),
pointer :: a
248 if (a%heap_size .gt. 0)
then
251 write (*,*)
'PQ_MAX: ERROR, heap_size < 1'
258 type(
pq),
pointer :: a
260 if (a%heap_size .gt. 0)
then
263 write (*,*)
'PQ_MAX_PRI: ERROR, heapsize < 1'
275 type(
pq),
pointer :: a
278 if (a%heap_size .ge. 1)
then
287 a%elems(1) = a%elems(a%heap_size)
288 a%heap_size = a%heap_size-1
292 write (*,*)
'PQ_EXTRACT_MAX: error, attempted to pop non-positive PQ'
304 type(
pq),
pointer :: a
305 real(kdkind),
intent(in) :: dis
306 integer,
intent(in) :: idx
309 integer :: i, isparent
310 real(kdkind) :: parentdis
317 a%heap_size = a%heap_size + 1
322 parentdis = a%elems(isparent)%dis
323 if (dis .gt. parentdis)
then
325 a%elems(i)%dis = parentdis
326 a%elems(i)%idx = a%elems(isparent)%idx
343 subroutine pq_adjust_heap(a,i)
344 type(
pq),
pointer :: a
345 integer,
intent(in) :: i
352 real(kdkind) :: prichild
353 integer :: parent, child, n
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
369 prichild = a%elems(child)%dis
370 if (e%dis .ge. prichild)
then
374 a%elems(parent) = a%elems(child)
381 end subroutine pq_adjust_heap
391 type(
pq),
pointer :: a
392 real(kdkind),
intent(in) :: dis
393 integer,
intent(in) :: idx
397 integer :: parent, child, n
398 real(kdkind) :: prichild, prichildp1
408 loop:
do while (child .le. n)
409 prichild = a%elems(child)%dis
416 if (child .lt. n)
then
417 prichildp1 = a%elems(child+1)%dis
418 if (prichild .lt. prichildp1)
then
420 prichild = prichildp1
424 if (dis .ge. prichild)
then
430 a%elems(parent) = a%elems(child)
435 a%elems(parent)%dis = dis
436 a%elems(parent)%idx = idx
459 type(
pq),
pointer :: a
462 if ((i .lt. 1) .or. (i .gt. a%heap_size))
then
463 write (*,*)
'PQ_DELETE: error, attempt to remove out of bounds element.'
469 a%elems(i) = a%elems(a%heap_size)
470 a%heap_size = a%heap_size - 1
533 real(
kdkind) :: cut_val_left, cut_val_right
546 integer :: dimen=0, n=0
548 real(
kdkind),
pointer :: the_data(:,:) => null()
562 integer,
pointer :: ind(:) => null()
566 logical :: sort = .false.
568 logical :: rearrange = .false.
569 real(
kdkind),
pointer :: rearranged_data(:,:) => null()
579 type :: tree_search_record
589 integer :: nn, nfound
591 integer :: centeridx=999, correltime=9999
597 real(
kdkind),
pointer :: qv(:)
600 real(
kdkind),
pointer :: data(:,:)
601 integer,
pointer :: ind(:)
602 end type tree_search_record
607 type(tree_search_record),
save,
target :: sr
636 integer,
intent(in),
optional :: dim
637 logical,
intent(in),
optional :: sort
638 logical,
intent(in),
optional :: rearrange
641 real(
kdkind),
target :: input_data(:,:)
646 mr%the_data => input_data
649 if (
present(dim))
then
652 mr%dimen =
size(input_data,1)
654 mr%n =
size(input_data,2)
656 if (mr%dimen > mr%n)
then
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.'
669 if (
present(sort))
then
675 if (
present(rearrange))
then
676 mr%rearrange = rearrange
678 mr%rearrange = .true.
681 if (mr%rearrange)
then
682 allocate(mr%rearranged_data(mr%dimen,mr%n))
684 mr%rearranged_data(:,i) = mr%the_data(:, &
688 nullify(mr%rearranged_data)
693 subroutine build_tree(tp)
697 type(
tree_node),
pointer :: dummy => null()
699 allocate (tp%ind(tp%n))
703 tp%root => build_tree_for_range(tp,1,tp%n, dummy)
704 end subroutine build_tree
706 recursive function build_tree_for_range(tp,l,u,parent)
result (res)
715 integer,
intent (In) :: l, u
718 integer :: i, c, m, dimen
737 allocate(res%box(dimen))
751 call spread_in_coordinate(tp,i,l,u,res%box(i))
772 if (
associated(parent))
then
773 if (i .ne. parent%cut_dim)
then
778 call spread_in_coordinate(tp,i,l,u,res%box(i))
780 res%box(i) = parent%box(i)
785 c = maxloc(res%box(1:dimen)%upper-res%box(1:dimen)%lower,1)
793 call select_on_coordinate(tp%the_data,tp%ind,c,m,l,u)
802 average = sum(tp%the_data(c,tp%ind(l:u))) / real(u-l+1,
kdkind)
804 average = (res%box(c)%upper + res%box(c)%lower)/2.0
807 res%cut_val = average
808 m = select_on_coordinate_value(tp%the_data,tp%ind,c,average,l,u)
817 res%left => build_tree_for_range(tp,l,m,res)
818 res%right => build_tree_for_range(tp,m+1,u,res)
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
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
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)
842 end function build_tree_for_range
844 integer function select_on_coordinate_value(v,ind,c,alpha,li,ui) &
863 integer,
intent (In) :: c, li, ui
864 real(
kdkind),
intent(in) :: alpha
886 if ( v(c,ind(lb)) <= alpha )
then
891 tmp = ind(lb); ind(lb) = ind(rb); ind(rb) = tmp
897 if (v(c,ind(lb)) <= alpha)
then
903 end function select_on_coordinate_value
905 subroutine select_on_coordinate(v,ind,c,k,li,ui)
910 integer,
intent (In) :: c, k, li, ui
912 integer :: i, l, m, s, t, u
923 if (v(c,ind(i))<v(c,t))
then
936 end subroutine select_on_coordinate
938 subroutine spread_in_coordinate(tp,c,l,u,interv)
945 type(
interval),
intent(out) :: interv
948 integer,
intent (In) :: c, l, u
951 real(
kdkind) :: last, lmax, lmin, t, smin,smax
955 real(
kdkind),
pointer :: v(:,:)
956 integer,
pointer :: ind(:)
958 v => tp%the_data(1:,1:)
965 do i = l + 2, ulocal, 2
973 if (smin>lmin) smin = lmin
974 if (smax<lmax) smax = lmax
976 if (i==ulocal+1)
then
977 last = v(c,ind(ulocal))
978 if (smin>last) smin = last
979 if (smax<last) smax = last
985 end subroutine spread_in_coordinate
993 call destroy_node(tp%root)
998 if (tp%rearrange)
then
999 deallocate(tp%rearranged_data)
1000 nullify(tp%rearranged_data)
1007 recursive subroutine destroy_node(np)
1012 intrinsic ASSOCIATED
1014 if (
associated(np%left))
then
1015 call destroy_node(np%left)
1018 if (
associated(np%right))
then
1019 call destroy_node(np%right)
1022 if (
associated(np%box))
deallocate(np%box)
1026 end subroutine destroy_node
1035 real(
kdkind),
target,
intent (In) :: qv(:)
1036 integer,
intent (In) :: nn
1040 sr%ballsize = huge(1.0)
1046 sr%overflow = .false.
1048 sr%results => results
1053 sr%rearrange = tp%rearrange
1054 if (tp%rearrange)
then
1055 sr%Data => tp%rearranged_data
1057 sr%Data => tp%the_data
1061 call validate_query_storage(nn)
1064 call search(tp%root)
1078 integer,
intent (In) :: idxin, correltime, nn
1081 allocate (sr%qv(tp%dimen))
1082 sr%qv = tp%the_data(:,idxin)
1083 sr%ballsize = huge(1.0)
1084 sr%centeridx = idxin
1085 sr%correltime = correltime
1093 sr%results => results
1096 sr%rearrange = tp%rearrange
1098 if (sr%rearrange)
then
1099 sr%Data => tp%rearranged_data
1101 sr%Data => tp%the_data
1104 call validate_query_storage(nn)
1107 call search(tp%root)
1128 real(
kdkind),
target,
intent (In) :: qv(:)
1129 real(
kdkind),
intent(in) :: r2
1130 integer,
intent(out) :: nfound
1131 integer,
intent (In) :: nalloc
1142 sr%results => results
1144 call validate_query_storage(nalloc)
1146 sr%overflow = .false.
1148 sr%rearrange= tp%rearrange
1150 if (tp%rearrange)
then
1151 sr%Data => tp%rearranged_data
1153 sr%Data => tp%the_data
1162 call search(tp%root)
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.'
1178 nfound,nalloc,results)
1186 integer,
intent (In) :: idxin, correltime, nalloc
1187 real(
kdkind),
intent(in) :: r2
1188 integer,
intent(out) :: nfound
1194 allocate (sr%qv(tp%dimen))
1195 sr%qv = tp%the_data(:,idxin)
1199 sr%centeridx = idxin
1200 sr%correltime = correltime
1202 sr%results => results
1205 sr%overflow = .false.
1207 call validate_query_storage(nalloc)
1213 sr%rearrange = tp%rearrange
1215 if (tp%rearrange)
then
1216 sr%Data => tp%rearranged_data
1218 sr%Data => tp%the_data
1220 sr%rearrange = tp%rearrange
1228 call search(tp%root)
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.'
1247 real(
kdkind),
target,
intent (In) :: qv(:)
1248 real(
kdkind),
intent(in) :: r2
1267 sr%rearrange = tp%rearrange
1268 if (tp%rearrange)
then
1269 sr%Data => tp%rearranged_data
1271 sr%Data => tp%the_data
1279 sr%overflow = .false.
1281 call search(tp%root)
1294 integer,
intent (In) :: correltime, idxin
1295 real(
kdkind),
intent(in) :: r2
1302 allocate (sr%qv(tp%dimen))
1303 sr%qv = tp%the_data(:,idxin)
1308 sr%centeridx = idxin
1309 sr%correltime = correltime
1316 sr%rearrange = tp%rearrange
1318 if (sr%rearrange)
then
1319 sr%Data => tp%rearranged_data
1321 sr%Data => tp%the_data
1329 sr%overflow = .false.
1331 call search(tp%root)
1339 subroutine validate_query_storage(n)
1343 integer,
intent(in) :: n
1345 if (
size(sr%results,1) .lt. n)
then
1346 write (*,*)
'KD_TREE_TRANS: you did not provide enough storage for results(1:n)'
1352 end subroutine validate_query_storage
1354 function square_distance(d, iv,qv)
result (res)
1365 real(
kdkind) :: iv(:),qv(:)
1368 res = sum( (iv(1:d)-qv(1:d))**2 )
1369 end function square_distance
1371 recursive subroutine search(node)
1381 type(
tree_node),
pointer :: ncloser, nfarther
1383 integer :: cut_dim, i
1385 real(
kdkind) :: qval, dis
1387 real(
kdkind),
pointer :: qv(:)
1390 if ((
associated(node%left) .and.
associated(node%right)) .eqv. .false.)
then
1392 if (sr%nn .eq. 0)
then
1393 call process_terminal_node_fixedball(node)
1395 call process_terminal_node(node)
1400 cut_dim = node%cut_dim
1403 if (qval < node%cut_val)
then
1404 ncloser => node%left
1405 nfarther => node%right
1406 dis = (node%cut_val_right - qval)**2
1409 ncloser => node%right
1410 nfarther => node%left
1411 dis = (node%cut_val_left - qval)**2
1415 if (
associated(ncloser))
call search(ncloser)
1418 if (
associated(nfarther))
then
1419 ballsize = sr%ballsize
1421 if (dis <= ballsize)
then
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
1441 call search(nfarther)
1445 end subroutine search
1448 real(
kdkind) function dis2_from_bnd(x,amin,amax) result (res)
1449 real(
kdkind),
intent(in) :: x, amin,amax
1464 end function dis2_from_bnd
1466 logical function box_in_search_range(node, sr)
result(res)
1474 type (tree_search_record),
pointer :: sr
1477 real(kdkind) :: dis, ballsize
1478 real(kdkind) :: l, u
1481 ballsize = sr%ballsize
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
1495 end function box_in_search_range
1498 subroutine process_terminal_node(node)
1505 real(kdkind),
pointer :: qv(:)
1506 integer,
pointer :: ind(:)
1507 real(kdkind),
pointer :: data(:,:)
1509 integer :: dimen, i, indexofi, k, centeridx, correltime
1510 real(kdkind) :: ballsize, sd, newpri
1511 logical :: rearrange
1512 type(pq),
pointer :: pqp
1523 ballsize = sr%ballsize
1524 rearrange = sr%rearrange
1526 data => sr%Data(1:,1:)
1527 centeridx = sr%centeridx
1528 correltime = sr%correltime
1534 mainloop:
do i = node%l, node%u
1538 sd = sd + (
data(k,i) - qv(k))**2
1539 if (sd>ballsize) cycle mainloop
1546 sd = sd + (
data(k,indexofi) - qv(k))**2
1547 if (sd>ballsize) cycle mainloop
1551 if (centeridx > 0)
then
1552 if (abs(indexofi-centeridx) < correltime) cycle mainloop
1574 if (sr%nfound .lt. sr%nn)
then
1578 sr%nfound = sr%nfound +1
1579 newpri = pq_insert(pqp,sd,indexofi)
1580 if (sr%nfound .eq. sr%nn) ballsize = newpri
1593 ballsize = pq_replace_max(pqp,sd,indexofi)
1599 sr%ballsize = ballsize
1601 end subroutine process_terminal_node
1603 subroutine process_terminal_node_fixedball(node)
1611 real(kdkind),
pointer :: qv(:)
1612 integer,
pointer :: ind(:)
1613 real(kdkind),
pointer :: data(:,:)
1616 integer :: dimen, i, indexofi, k
1617 integer :: centeridx, correltime, nn
1618 real(kdkind) :: ballsize, sd
1619 logical :: rearrange
1626 ballsize = sr%ballsize
1627 rearrange = sr%rearrange
1629 data => sr%Data(1:,1:)
1630 centeridx = sr%centeridx
1631 correltime = sr%correltime
1636 mainloop:
do i = node%l, node%u
1661 sd = sd + (
data(k,i) - qv(k))**2
1662 if (sd>ballsize) cycle mainloop
1669 sd = sd + (
data(k,indexofi) - qv(k))**2
1670 if (sd>ballsize) cycle mainloop
1674 if (centeridx > 0)
then
1675 if (abs(indexofi-centeridx)<correltime) cycle mainloop
1679 if (nfound .gt. sr%nalloc)
then
1682 sr%overflow = .true.
1684 sr%results(nfound)%dis = sd
1685 sr%results(nfound)%idx = indexofi
1692 end subroutine process_terminal_node_fixedball
1700 real(kdkind),
intent (In) :: qv(:)
1701 integer,
intent (In) :: nn
1702 type(kdtree2_result) :: results(:)
1705 real(kdkind),
allocatable :: all_distances(:)
1707 allocate (all_distances(tp%n))
1709 all_distances(i) = square_distance(tp%dimen,qv,tp%the_data(:,i))
1713 results(i)%dis = huge(1.0)
1717 if (all_distances(i)<results(nn)%dis)
then
1720 if (all_distances(i)<results(j)%dis)
exit
1723 do k = nn - 1, j, -1
1724 results(k+1) = results(k)
1726 results(j)%dis = all_distances(i)
1730 deallocate (all_distances)
1740 real(kdkind),
intent (In) :: qv(:)
1741 real(kdkind),
intent (In) :: r2
1742 integer,
intent(out) :: nfound
1743 type(kdtree2_result) :: results(:)
1745 integer :: i, nalloc
1746 real(kdkind),
allocatable :: all_distances(:)
1748 allocate (all_distances(tp%n))
1750 all_distances(i) = square_distance(tp%dimen,qv,tp%the_data(:,i))
1754 nalloc =
size(results,1)
1757 if (all_distances(i)< r2)
then
1759 if (nfound .lt. nalloc)
then
1761 results(nfound)%dis = all_distances(i)
1762 results(nfound)%idx = i
1766 deallocate (all_distances)
1776 integer,
intent(in) :: nfound
1777 type(kdtree2_result),
target :: results(:)
1784 if (nfound .gt. 1)
call heapsort_struct(results,nfound)
1789 subroutine heapsort(a,ind,n)
1795 integer,
intent(in) :: n
1796 real(kdkind),
intent(inout) :: a(:)
1797 integer,
intent(inout) :: ind(:)
1801 real(kdkind) ::
value
1805 integer :: ileft,iright
1820 value=a(ileft); ivalue=ind(ileft)
1822 value=a(iright); ivalue=ind(iright)
1823 a(iright)=a(1); ind(iright)=ind(1)
1825 if (iright == 1)
then
1826 a(1)=
value;ind(1)=ivalue
1832 do while (j <= iright)
1834 if(a(j) < a(j+1)) j=j+1
1836 if(
value < a(j))
then
1837 a(i)=a(j); ind(i)=ind(j)
1844 a(i)=
value; ind(i)=ivalue
1846 end subroutine heapsort
1848 subroutine heapsort_struct(a,n)
1853 integer,
intent(in) :: n
1854 type(kdtree2_result),
intent(inout) :: a(:)
1858 type(kdtree2_result) ::
value
1861 integer :: ileft,iright
1881 if (iright == 1)
then
1888 do while (j <= iright)
1890 if(a(j)%dis < a(j+1)%dis) j=j+1
1892 if(
value%dis < a(j)%dis)
then
1902 end subroutine heapsort_struct
type(kdtree2) function, pointer, public kdtree2_create(input_data, dim, sort, rearrange)
subroutine, public kdtree2_r_nearest_brute_force(tp, qv, r2, nfound, results)
subroutine, public kdtree2_sort_results(nfound, results)
subroutine, public kdtree2_n_nearest_brute_force(tp, qv, nn, results)
integer function, public kdtree2_r_count_around_point(tp, idxin, correltime, r2)
subroutine, public kdtree2_n_nearest_around_point(tp, idxin, correltime, nn, results)
integer, parameter bucket_size
subroutine, public kdtree2_destroy(tp)
subroutine, public kdtree2_n_nearest(tp, qv, nn, results)
subroutine, public kdtree2_r_nearest(tp, qv, r2, nfound, nalloc, results)
subroutine, public kdtree2_r_nearest_around_point(tp, idxin, correltime, r2, nfound, nalloc, results)
integer function, public kdtree2_r_count(tp, qv, r2)
integer, parameter, public kdkind
subroutine, public pq_extract_max(a, e)
subroutine, public pq_delete(a, i)
real(kdkind) function, public pq_replace_max(a, dis, idx)
real(kdkind) function, public pq_maxpri(a)
subroutine, public pq_max(a, e)
type(pq) function, public pq_create(results_in)
real(kdkind) function, public pq_insert(a, dis, idx)