MODFLOW 6  version 6.7.0.dev0
USGS Modular Hydrologic Model
qsort_inline.inc
Go to the documentation of this file.
1 !======================================================================
2 ! Fast in-line QSORT+INSERTION SORT for Fortran.
3 ! Author: Joseph M. Krahn
4 ! FILE: qsort_inline.inc
5 ! SOURCE: http://fortranwiki.org/fortran/show/qsort_inline
6 ! PURPOSE:
7 ! Generate a custom array sort procedure for a specific type,
8 ! without the comparison-callback overhead of a generic sort procedure.
9 ! This is essentially the same as an in-line optimization, which generally
10 ! is not feasible for a library-based generic sort procedure.
11 !
12 ! This implementation is as generic as possible, while avoiding the need
13 ! for a code pre-processor. The success of this approach assumes that
14 ! internal procedures are always in-lined with optimized Fortran compilation.
15 !
16 ! USAGE:
17 ! This file contains the sort subroutine body. You must supply
18 ! an integer parameter QSORT_THRESHOLD, and internal procedures:
19 ! subroutine INIT()
20 ! logical function lessThan(a,b)
21 ! subroutine SWAP(a,b)
22 ! subroutine RSHIFT(left,right)
23 !
24 ! Descriptions:
25 !
26 ! SUBROUTINE INIT()
27 ! Any user initialization code. This is needed because executable
28 ! statements cannot precede this code, which begins with declarations.
29 ! In many cases, this is just an empty procedure.
30 !
31 ! LOGICAL FUNCTION lessThan(a,b)
32 ! Return TRUE if array member 'a' is less than array member 'b'.
33 ! Only a TRUE value causes a change in sort order. This minimizes data
34 ! manipulation, and maintains the original order for values that are
35 ! equivalent by the sort comparison. It also avoids the need to
36 ! distinguish equality from greater-than.
37 !
38 ! SUBROUTINE SWAP(A,B)
39 ! Swap array members 'a' and 'b'
40 !
41 ! SUBROUTINE RSHIFT(LEFT,RIGHT)
42 ! Perform a circular shift of array members LEFT through RIGHT,
43 ! shifting the element at RIGHT back to the position at LEFT.
44 !
45 ! QSORT_THRESHOLD:
46 ! The QSORT is used down to the QSORT_THRESHOLD size sorted blocks.
47 ! Then insertion sort is used for the remainder, because it is faster
48 ! for small sort ranges. The optimal size is not critical. Most of
49 ! the benefit is in blocks of 8 or less, and values of 16 to 128
50 ! are generally about equal speed. However, the optimal value
51 ! depends a lot on the hardware and the data being sorted, so this
52 ! is left as a tunable parameter for cases where there is an
53 ! effect on performance.
54 !
55 !---------------------------------------------------------------------
56 ! NOTES:
57 ! The procedure uses a optimized combination of QSORT and INSERTION
58 ! sorting. The algorithm is based on code used in GLIBC.
59 ! A stack is used in place of recursive calls. The stack size must
60 ! be at least as big as the number of bits in the largest array index.
61 !
62 ! Sorting vectors of a multidimensional allocatable array can be
63 ! VERY slow. In this case, or with large derived types, it is better
64 ! to sort a simple derived type of key/index pairs, then reorder
65 ! the actual data using the sorted indices.
66 !
67 !---------------------------------------------------------------------
68  integer :: stack_top, right_size, left_size
69  integer :: mid, left, right, low, high
70 
71 ! A stack of 32 can handle the entire extent of a 32-bit
72 ! index, so this value is fixed. If you have 64-bit indexed
73 ! arrays, which might contain more than 2^32 elements, this
74 ! should be set to 64.
75  integer, parameter :: QSORT_STACK_SIZE = 64
76  type qsort_stack; integer :: low, high; end type
77  type(qsort_stack) :: stack(QSORT_STACK_SIZE)
78 
79  call init()
80 
81  if (arraySize > QSORT_THRESHOLD) then
82  low = 1
83  high = arraySize
84  stack_top = 0
85 
86  QSORT_LOOP: &
87  do
88  mid = (low + high)/2
89  if (lessThan (mid, low)) then
90  call SWAP(mid,low)
91  end if
92  if (lessThan (high, mid)) then
93  call SWAP(high,mid)
94  if (lessThan (mid, low)) then
95  call SWAP(mid,low)
96  end if
97  end if
98  left = low + 1
99  right = high - 1
100 
101  COLLAPSE_WALLS: &
102  do
103  do while (lessThan (left, mid))
104  left=left+1
105  end do
106  do while (lessThan (mid, right))
107  right=right-1
108  end do
109  if (left < right) then
110  call SWAP(left,right)
111  if (mid == left) then
112  mid = right
113  else if (mid == right) then
114  mid = left
115  end if
116  left=left+1
117  right=right-1
118  else
119  if (left == right) then
120  left=left+1
121  right=right-1
122  end if
123  exit COLLAPSE_WALLS
124  end if
125  end do COLLAPSE_WALLS
126 
127 ! Set up indices for the next iteration.
128 ! Determine left and right partition sizes.
129 ! Defer partitions smaller than the QSORT_THRESHOLD.
130 ! If both partitions are significant,
131 ! push the larger one onto the stack.
132  right_size = right - low
133  left_size = high - left
134  if (right_size <= QSORT_THRESHOLD) then
135  if (left_size <= QSORT_THRESHOLD) then
136  ! Ignore both small partitions: Pop a partition or exit.
137  if (stack_top<1) exit QSORT_LOOP
138  low=stack(stack_top)%low; high=stack(stack_top)%high
139  stack_top=stack_top-1
140  else
141  ! Ignore small left partition.
142  low = left
143  end if
144  else if (left_size <= QSORT_THRESHOLD) then
145  ! Ignore small right partition.
146  high = right
147  else if (right_size > left_size) then
148  ! Push larger left partition indices.
149  stack_top=stack_top+1
150  stack(stack_top)=qsort_stack(low,right)
151  low = left
152  else
153  ! Push larger right partition indices.
154  stack_top=stack_top+1
155  stack(stack_top)=qsort_stack(left,high)
156  high = right
157  end if
158  end do QSORT_LOOP
159  end if ! (arraySize > QSORT_THRESHOLD)
160 
161 ! Sort the remaining small partitions using insertion sort,
162 ! which should be faster for partitions smaller than the
163 ! appropriate QSORT_THRESHOLD.
164 
165 ! First, find smallest element in first QSORT_THRESHOLD and
166 ! place it at the array's beginning. This places a lower
167 ! bound 'guard' position, and speeds up the inner loop
168 ! below, because it will not need a lower-bound test.
169  low = 1
170  high = arraySize
171 
172 ! left is the MIN_LOC index here:
173  left=low
174  do right = low+1, min(low+QSORT_THRESHOLD,high)
175  if (lessThan(right,left)) left=right
176  end do
177  if (left/=low) call SWAP(left,low)
178 
179 ! Insertion sort, from left to right.
180 ! (assuming that the left is the lowest numbered index)
181  INSERTION_SORT: &
182  do right = low+2,high
183  left=right-1
184  if (lessThan(right,left)) then
185  do while (lessThan(right,left-1))
186  left=left-1
187  end do
188  call RSHIFT(left,right)
189  end if
190  end do INSERTION_SORT
191 !--------------------------------------------------------------
subroutine init()
Definition: GridSorting.f90:24
subroutine next(this)
Increment the iterator to the next node.
class(*) function, pointer value(this)
Get the value the iterator is pointing at.
logical function contains(this, key)
Boolean indicating if an item exists in the hashtable.
integer(i4b) function at(this, idx)
Definition: STLVecInt.f90:102
subroutine sort(this)
Sort the time selection in increasing order.
Definition: TimeSelect.f90:253