This subroutine sorts out the least element of t, and puts the remaining elements of t in a heap.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
integer, | intent(in) | :: | n |
the dimension of the arrays t and iorder. |
||
real(kind=wp), | intent(inout) | :: | t(n) |
On entry t stores the elements to be sorted, On exit t(n) stores the least elements of t, and t(1) to t(n-1) stores the remaining elements in the form of a heap. |
||
integer, | intent(inout) | :: | Iorder(n) |
On entry iorder(i) is the index of t(i). On exit iorder(i) is still the index of t(i), but iorder may be permuted in accordance with t. |
||
integer, | intent(in) | :: | Iheap |
|
subroutine hpsolb(n,t,Iorder,Iheap) implicit none integer,intent(in) :: n !! the dimension of the arrays t and iorder. real(wp),intent(inout) :: t(n) !! On entry t stores the elements to be sorted, !! On exit t(n) stores the least elements of t, and t(1) to t(n-1) !! stores the remaining elements in the form of a heap. integer,intent(inout) :: Iorder(n) !! On entry iorder(i) is the index of t(i). !! On exit iorder(i) is still the index of t(i), but iorder may be !! permuted in accordance with t. integer,intent(in) :: Iheap !! `iheap == 0` if t(1) to t(n) is not in the form of a heap, !! `iheap /= 0` if otherwise. integer :: i , j , k , indxin , indxou real(wp) :: ddum , out if ( Iheap==0 ) then ! Rearrange the elements t(1) to t(n) to form a heap. do k = 2 , n ddum = t(k) indxin = Iorder(k) ! Add ddum to the heap. i = k do if ( i>1 ) then j = i/2 if ( ddum<t(j) ) then t(i) = t(j) Iorder(i) = Iorder(j) i = j cycle endif endif exit end do t(i) = ddum Iorder(i) = indxin enddo endif ! Assign to 'out' the value of t(1), the least member of the heap, ! and rearrange the remaining members to form a heap as ! elements 1 to n-1 of t. if ( n>1 ) then i = 1 out = t(1) indxou = Iorder(1) ddum = t(n) indxin = Iorder(n) ! Restore the heap do j = i + i if ( j<=n-1 ) then if ( t(j+1)<t(j) ) j = j + 1 if ( t(j)<ddum ) then t(i) = t(j) Iorder(i) = Iorder(j) i = j cycle endif endif exit end do t(i) = ddum Iorder(i) = indxin ! Put the least member in t(n). t(n) = out Iorder(n) = indxou endif end subroutine hpsolb