MODFLOW 6  version 6.6.0.dev0
USGS Modular Hydrologic Model
sortmodule Module Reference

Data Types

interface  qsort
 
interface  unique_values
 

Functions/Subroutines

subroutine qsort_int1d (indx, v, reverse)
 
subroutine qsort_dbl1d (indx, v, reverse)
 
subroutine unique_values_int1d (a, b)
 
subroutine unique_values_dbl1d (a, b)
 
subroutine, public selectn (indx, v, reverse)
 
subroutine rswap (a, b)
 
subroutine iswap (ia, ib)
 

Function/Subroutine Documentation

◆ iswap()

subroutine sortmodule::iswap ( integer(i4b), intent(inout)  ia,
integer(i4b), intent(inout)  ib 
)
private

Definition at line 476 of file sort.f90.

477  ! -- dummy arguments
478  integer(I4B), intent(inout) :: ia
479  integer(I4B), intent(inout) :: ib
480  ! -- local variables
481  integer(I4B) :: id
482  ! -- functions
483  ! -- code
484  id = ia
485  ia = ib
486  ib = id
Here is the caller graph for this function:

◆ qsort_dbl1d()

subroutine sortmodule::qsort_dbl1d ( integer(i4b), dimension(:), intent(inout)  indx,
real(dp), dimension(:), intent(inout)  v,
logical, intent(in), optional  reverse 
)
private

Definition at line 151 of file sort.f90.

152  ! -- dummy arguments
153  integer(I4B), dimension(:), intent(inout) :: indx
154  real(DP), dimension(:), intent(inout) :: v
155  logical, intent(in), optional :: reverse
156  ! -- local variables
157  logical :: lrev
158  integer(I4B), parameter :: nn = 15
159  integer(I4B), parameter :: nstack = 50
160  integer(I4B) :: nsize
161  integer(I4B) :: k
162  integer(I4B) :: i
163  integer(I4B) :: j
164  integer(I4B) :: jstack
165  integer(I4B) :: ileft
166  integer(I4B) :: iright
167  integer(I4B), dimension(nstack) :: istack
168  integer(I4B) :: iidx
169  integer(I4B) :: ia
170  real(DP) :: a
171  ! -- functions
172  ! -- code
173  !
174  ! -- process optional dummy variables
175  if (present(reverse)) then
176  lrev = reverse
177  else
178  lrev = .false.
179  end if
180  !
181  ! -- initialize variables
182  nsize = size(v)
183  jstack = 0
184  ileft = 1
185  iright = nsize
186  !
187  ! -- perform quicksort
188  do
189  if (iright - ileft < nn) then
190  do j = (ileft + 1), iright
191  a = v(j)
192  iidx = indx(j)
193  do i = (j - 1), ileft, -1
194  if (v(i) <= a) exit
195  v(i + 1) = v(i)
196  indx(i + 1) = indx(i)
197  end do
198  v(i + 1) = a
199  indx(i + 1) = iidx
200  end do
201  if (jstack == 0) return
202  iright = istack(jstack)
203  ileft = istack(jstack - 1)
204  jstack = jstack - 2
205  else
206  k = (ileft + iright) / 2
207  call rswap(v(k), v(ileft + 1))
208  call iswap(indx(k), indx(ileft + 1))
209  if (v(ileft) > v(iright)) then
210  call rswap(v(ileft), v(iright))
211  call iswap(indx(ileft), indx(iright))
212  end if
213  if (v(ileft + 1) > v(iright)) then
214  call rswap(v(ileft + 1), v(iright))
215  call iswap(indx(ileft + 1), indx(iright))
216  end if
217  if (v(ileft) > v(ileft + 1)) then
218  call rswap(v(ileft), v(ileft + 1))
219  call iswap(indx(ileft), indx(ileft + 1))
220  end if
221  i = ileft + 1
222  j = iright
223  a = v(ileft + 1)
224  ia = indx(ileft + 1)
225  do
226  do
227  i = i + 1
228  if (v(i) >= a) then
229  exit
230  end if
231  end do
232  do
233  j = j - 1
234  if (v(j) <= a) then
235  exit
236  end if
237  end do
238  if (j < i) then
239  exit
240  end if
241  call rswap(v(i), v(j))
242  call iswap(indx(i), indx(j))
243  end do
244  v(ileft + 1) = v(j)
245  indx(ileft + 1) = indx(j)
246  v(j) = a
247  indx(j) = ia
248  jstack = jstack + 2
249  if (jstack > nstack) then
250  write (errmsg, '(a,3(1x,a))') &
251  'JSTACK > NSTACK IN SortModule::qsort'
252  call store_error(errmsg, terminate=.true.)
253  end if
254  if ((iright - i + 1) >= (j - 1)) then
255  istack(jstack) = iright
256  istack(jstack - 1) = i
257  iright = j - 1
258  else
259  istack(jstack) = j - 1
260  istack(jstack - 1) = ileft
261  ileft = i
262  end if
263  end if
264  end do
265  !
266  ! -- reverse order of the heap index
267  if (lrev) then
268  j = nsize
269  do i = 1, nsize / 2
270  call rswap(v(i), v(j))
271  call iswap(indx(i), indx(j))
272  j = j - 1
273  end do
274  end if

◆ qsort_int1d()

subroutine sortmodule::qsort_int1d ( integer(i4b), dimension(:), intent(inout)  indx,
integer(i4b), dimension(:), intent(inout)  v,
logical, intent(in), optional  reverse 
)
private

Definition at line 23 of file sort.f90.

24  ! -- dummy arguments
25  integer(I4B), dimension(:), intent(inout) :: indx
26  integer(I4B), dimension(:), intent(inout) :: v
27  logical, intent(in), optional :: reverse
28  ! -- local variables
29  logical :: lrev
30  integer(I4B), parameter :: nn = 15
31  integer(I4B), parameter :: nstack = 50
32  integer(I4B) :: nsize
33  integer(I4B) :: k
34  integer(I4B) :: i
35  integer(I4B) :: j
36  integer(I4B) :: jstack
37  integer(I4B) :: ileft
38  integer(I4B) :: iright
39  integer(I4B), dimension(nstack) :: istack
40  integer(I4B) :: iidx
41  integer(I4B) :: ia
42  integer(I4B) :: a
43  ! -- functions
44  ! -- code
45  !
46  ! -- process optional dummy variables
47  if (present(reverse)) then
48  lrev = reverse
49  else
50  lrev = .false.
51  end if
52  !
53  ! -- initialize variables
54  nsize = size(v)
55  jstack = 0
56  ileft = 1
57  iright = nsize
58  !
59  ! -- perform quicksort
60  do
61  if (iright - ileft < nn) then
62  do j = (ileft + 1), iright
63  a = v(j)
64  iidx = indx(j)
65  do i = (j - 1), ileft, -1
66  if (v(i) <= a) exit
67  v(i + 1) = v(i)
68  indx(i + 1) = indx(i)
69  end do
70  v(i + 1) = a
71  indx(i + 1) = iidx
72  end do
73  if (jstack == 0) return
74  iright = istack(jstack)
75  ileft = istack(jstack - 1)
76  jstack = jstack - 2
77  else
78  k = (ileft + iright) / 2
79  call iswap(v(k), v(ileft + 1))
80  call iswap(indx(k), indx(ileft + 1))
81  if (v(ileft) > v(iright)) then
82  call iswap(v(ileft), v(iright))
83  call iswap(indx(ileft), indx(iright))
84  end if
85  if (v(ileft + 1) > v(iright)) then
86  call iswap(v(ileft + 1), v(iright))
87  call iswap(indx(ileft + 1), indx(iright))
88  end if
89  if (v(ileft) > v(ileft + 1)) then
90  call iswap(v(ileft), v(ileft + 1))
91  call iswap(indx(ileft), indx(ileft + 1))
92  end if
93  i = ileft + 1
94  j = iright
95  a = v(ileft + 1)
96  ia = indx(ileft + 1)
97  do
98  do
99  i = i + 1
100  if (v(i) >= a) then
101  exit
102  end if
103  end do
104  do
105  j = j - 1
106  if (v(j) <= a) then
107  exit
108  end if
109  end do
110  if (j < i) then
111  exit
112  end if
113  call iswap(v(i), v(j))
114  call iswap(indx(i), indx(j))
115  end do
116  v(ileft + 1) = v(j)
117  indx(ileft + 1) = indx(j)
118  v(j) = a
119  indx(j) = ia
120  jstack = jstack + 2
121  if (jstack > nstack) then
122  write (errmsg, '(a,3(1x,a))') &
123  'JSTACK > NSTACK IN SortModule::qsort'
124  call store_error(errmsg, terminate=.true.)
125  end if
126  if ((iright - i + 1) >= (j - 1)) then
127  istack(jstack) = iright
128  istack(jstack - 1) = i
129  iright = j - 1
130  else
131  istack(jstack) = j - 1
132  istack(jstack - 1) = ileft
133  ileft = i
134  end if
135  end if
136  end do
137  !
138  ! -- reverse order of the heap index
139  if (lrev) then
140  j = nsize
141  do i = 1, nsize / 2
142  call iswap(v(i), v(j))
143  call iswap(indx(i), indx(j))
144  j = j - 1
145  end do
146  end if

◆ rswap()

subroutine sortmodule::rswap ( real(dp), intent(inout)  a,
real(dp), intent(inout)  b 
)
private

Definition at line 463 of file sort.f90.

464  ! -- dummy arguments
465  real(DP), intent(inout) :: a
466  real(DP), intent(inout) :: b
467  ! -- local variables
468  real(DP) :: d
469  ! -- functions
470  ! -- code
471  d = a
472  a = b
473  b = d
Here is the caller graph for this function:

◆ selectn()

subroutine, public sortmodule::selectn ( integer(i4b), dimension(:), intent(inout)  indx,
real(dp), dimension(:), intent(inout)  v,
logical, intent(in), optional  reverse 
)

Definition at line 383 of file sort.f90.

384  ! -- dummy arguments
385  integer(I4B), dimension(:), intent(inout) :: indx
386  real(DP), dimension(:), intent(inout) :: v
387  logical, intent(in), optional :: reverse
388  ! -- local variables
389  logical :: lrev
390  integer(I4B) :: nsizei
391  integer(I4B) :: nsizev
392  integer(I4B) :: i
393  integer(I4B) :: j
394  integer(I4B) :: k
395  integer(I4B) :: n
396  !integer(I4B) :: iidx
397  real(DP), dimension(:), allocatable :: vv
398  ! -- functions
399  ! -- code
400  !
401  ! -- process optional dummy variables
402  if (present(reverse)) then
403  lrev = reverse
404  else
405  lrev = .false.
406  end if
407  !
408  ! -- initialize heap
409  nsizev = size(v)
410  nsizei = min(nsizev, size(indx))
411  allocate (vv(nsizei))
412  !
413  ! -- initialize heap index (indx) and heap (vv)
414  do n = 1, nsizei
415  vv(n) = v(n)
416  indx(n) = n
417  end do
418  !
419  ! -- initial sort
420  call qsort(indx, vv)
421  !
422  ! -- evaluate the remaining elements in v
423  do i = nsizei + 1, nsizev
424  !
425  ! -- put the current value on the heap
426  if (v(i) > vv(1)) then
427  vv(1) = v(i)
428  indx(1) = i
429  j = 1
430  do
431  k = 2 * j
432  if (k > nsizei) then
433  exit
434  end if
435  if (k /= nsizei) then
436  if (vv(k) > vv(k + 1)) then
437  k = k + 1
438  end if
439  end if
440  if (vv(j) <= vv(k)) then
441  exit
442  end if
443  call rswap(vv(k), vv(j))
444  call iswap(indx(k), indx(j))
445  j = k
446  end do
447  end if
448  end do
449  !
450  ! -- final sort
451  call qsort(indx, vv)
452  !
453  ! -- reverse order of the heap index
454  if (lrev) then
455  j = nsizei
456  do i = 1, nsizei / 2
457  call iswap(indx(i), indx(j))
458  j = j - 1
459  end do
460  end if
Here is the call graph for this function:
Here is the caller graph for this function:

◆ unique_values_dbl1d()

subroutine sortmodule::unique_values_dbl1d ( real(dp), dimension(:), intent(in), allocatable  a,
real(dp), dimension(:), intent(inout), allocatable  b 
)
private

Definition at line 329 of file sort.f90.

330  ! - dummy arguments
331  real(DP), dimension(:), allocatable, intent(in) :: a
332  real(DP), dimension(:), allocatable, intent(inout) :: b
333  ! -- local variables
334  integer(I4B) :: count
335  integer(I4B) :: n
336  integer(I4B), dimension(:), allocatable :: indxarr
337  real(DP), dimension(:), allocatable :: tarr
338  ! -- functions
339  ! -- code
340  !
341  ! -- allocate tarr and create idxarr
342  allocate (tarr(size(a)))
343  allocate (indxarr(size(a)))
344  !
345  ! -- fill tarr with a and create index
346  do n = 1, size(a)
347  tarr(n) = a(n)
348  indxarr(n) = n
349  end do
350  !
351  ! -- sort a in increasing order
352  call qsort(indxarr, tarr, reverse=.true.)
353  !
354  ! -- determine the number of unique values
355  count = 1
356  do n = 2, size(tarr)
357  if (tarr(n) > tarr(n - 1)) count = count + 1
358  end do
359  !
360  ! -- allocate b for unique values
361  if (allocated(b)) then
362  deallocate (b)
363  end if
364  allocate (b(count))
365  !
366  ! -- fill b with unique values
367  b(1) = tarr(1)
368  count = 1
369  do n = 2, size(a)
370  if (tarr(n) > b(count)) then
371  count = count + 1
372  b(count) = tarr(n)
373  end if
374  end do
375  !
376  ! -- allocate tarr and create idxarr
377  deallocate (tarr)
378  deallocate (indxarr)

◆ unique_values_int1d()

subroutine sortmodule::unique_values_int1d ( integer(i4b), dimension(:), intent(in), allocatable  a,
integer(i4b), dimension(:), intent(inout), allocatable  b 
)
private

Definition at line 277 of file sort.f90.

278  ! - dummy arguments
279  integer(I4B), dimension(:), allocatable, intent(in) :: a
280  integer(I4B), dimension(:), allocatable, intent(inout) :: b
281  ! -- local variables
282  integer(I4B) :: count
283  integer(I4B) :: n
284  integer(I4B), dimension(:), allocatable :: indxarr
285  integer(I4B), dimension(:), allocatable :: tarr
286  ! -- functions
287  ! -- code
288  !
289  ! -- allocate tarr and create idxarr
290  allocate (tarr(size(a)))
291  allocate (indxarr(size(a)))
292  !
293  ! -- fill tarr with a and create index
294  do n = 1, size(a)
295  tarr(n) = a(n)
296  indxarr(n) = n
297  end do
298  !
299  ! -- sort a in increasing order
300  call qsort(indxarr, tarr, reverse=.true.)
301  !
302  ! -- determine the number of unique values
303  count = 1
304  do n = 2, size(tarr)
305  if (tarr(n) > tarr(n - 1)) count = count + 1
306  end do
307  !
308  ! -- allocate b for unique values
309  if (allocated(b)) then
310  deallocate (b)
311  end if
312  allocate (b(count))
313  !
314  ! -- fill b with unique values
315  b(1) = tarr(1)
316  count = 1
317  do n = 2, size(a)
318  if (tarr(n) > b(count)) then
319  count = count + 1
320  b(count) = tarr(n)
321  end if
322  end do
323  !
324  ! -- allocate tarr and create idxarr
325  deallocate (tarr)
326  deallocate (indxarr)