MODFLOW 6  version 6.5.0.dev2
MODFLOW 6 Code Documentation
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 505 of file sort.f90.

506  ! -- dummy arguments
507  integer(I4B), intent(inout) :: ia
508  integer(I4B), intent(inout) :: ib
509  ! -- local variables
510  integer(I4B) :: id
511  ! -- functions
512  ! -- code
513  id = ia
514  ia = ib
515  ib = id
516  !
517  ! -- return
518  return
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 155 of file sort.f90.

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

◆ 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 20 of file sort.f90.

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

◆ rswap()

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

Definition at line 489 of file sort.f90.

490  ! -- dummy arguments
491  real(DP), intent(inout) :: a
492  real(DP), intent(inout) :: b
493  ! -- local variables
494  real(DP) :: d
495  ! -- functions
496  ! -- code
497  d = a
498  a = b
499  b = d
500  !
501  ! -- return
502  return
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 400 of file sort.f90.

401 ! **************************************************************************
402 ! selectn -- heap selection
403 ! **************************************************************************
404 !
405 ! SPECIFICATIONS:
406 ! --------------------------------------------------------------------------
407  ! -- dummy arguments
408  integer(I4B), dimension(:), intent(inout) :: indx
409  real(DP), dimension(:), intent(inout) :: v
410  logical, intent(in), optional :: reverse
411  ! -- local variables
412  logical :: lrev
413  integer(I4B) :: nsizei
414  integer(I4B) :: nsizev
415  integer(I4B) :: i
416  integer(I4B) :: j
417  integer(I4B) :: k
418  integer(I4B) :: n
419  !integer(I4B) :: iidx
420  real(DP), dimension(:), allocatable :: vv
421  ! -- functions
422  ! -- code
423  !
424  ! -- process optional dummy variables
425  if (present(reverse)) then
426  lrev = reverse
427  else
428  lrev = .false.
429  end if
430  !
431  ! -- initialize heap
432  nsizev = size(v)
433  nsizei = min(nsizev, size(indx))
434  allocate (vv(nsizei))
435  !
436  ! -- initialize heap index (indx) and heap (vv)
437  do n = 1, nsizei
438  vv(n) = v(n)
439  indx(n) = n
440  end do
441  !
442  ! -- initial sort
443  call qsort(indx, vv)
444  !
445  ! -- evaluate the remaining elements in v
446  do i = nsizei + 1, nsizev
447  !
448  ! -- put the current value on the heap
449  if (v(i) > vv(1)) then
450  vv(1) = v(i)
451  indx(1) = i
452  j = 1
453  do
454  k = 2 * j
455  if (k > nsizei) then
456  exit
457  end if
458  if (k /= nsizei) then
459  if (vv(k) > vv(k + 1)) then
460  k = k + 1
461  end if
462  end if
463  if (vv(j) <= vv(k)) then
464  exit
465  end if
466  call rswap(vv(k), vv(j))
467  call iswap(indx(k), indx(j))
468  j = k
469  end do
470  end if
471  end do
472  !
473  ! -- final sort
474  call qsort(indx, vv)
475  !
476  ! -- reverse order of the heap index
477  if (lrev) then
478  j = nsizei
479  do i = 1, nsizei / 2
480  call iswap(indx(i), indx(j))
481  j = j - 1
482  end do
483  end if
484  !
485  ! -- return
486  return
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 345 of file sort.f90.

346  ! - dummy arguments
347  real(DP), dimension(:), allocatable, intent(in) :: a
348  real(DP), dimension(:), allocatable, intent(inout) :: b
349  ! -- local variables
350  integer(I4B) :: count
351  integer(I4B) :: n
352  integer(I4B), dimension(:), allocatable :: indxarr
353  real(DP), dimension(:), allocatable :: tarr
354  ! -- functions
355  ! -- code
356  !
357  ! -- allocate tarr and create idxarr
358  allocate (tarr(size(a)))
359  allocate (indxarr(size(a)))
360  !
361  ! -- fill tarr with a and create index
362  do n = 1, size(a)
363  tarr(n) = a(n)
364  indxarr(n) = n
365  end do
366  !
367  ! -- sort a in increasing order
368  call qsort(indxarr, tarr, reverse=.true.)
369  !
370  ! -- determine the number of unique values
371  count = 1
372  do n = 2, size(tarr)
373  if (tarr(n) > tarr(n - 1)) count = count + 1
374  end do
375  !
376  ! -- allocate b for unique values
377  if (allocated(b)) then
378  deallocate (b)
379  end if
380  allocate (b(count))
381  !
382  ! -- fill b with unique values
383  b(1) = tarr(1)
384  count = 1
385  do n = 2, size(a)
386  if (tarr(n) > b(count)) then
387  count = count + 1
388  b(count) = tarr(n)
389  end if
390  end do
391  !
392  ! -- allocate tarr and create idxarr
393  deallocate (tarr)
394  deallocate (indxarr)
395  !
396  ! -- return
397  return

◆ 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 290 of file sort.f90.

291  ! - dummy arguments
292  integer(I4B), dimension(:), allocatable, intent(in) :: a
293  integer(I4B), dimension(:), allocatable, intent(inout) :: b
294  ! -- local variables
295  integer(I4B) :: count
296  integer(I4B) :: n
297  integer(I4B), dimension(:), allocatable :: indxarr
298  integer(I4B), dimension(:), allocatable :: tarr
299  ! -- functions
300  ! -- code
301  !
302  ! -- allocate tarr and create idxarr
303  allocate (tarr(size(a)))
304  allocate (indxarr(size(a)))
305  !
306  ! -- fill tarr with a and create index
307  do n = 1, size(a)
308  tarr(n) = a(n)
309  indxarr(n) = n
310  end do
311  !
312  ! -- sort a in increasing order
313  call qsort(indxarr, tarr, reverse=.true.)
314  !
315  ! -- determine the number of unique values
316  count = 1
317  do n = 2, size(tarr)
318  if (tarr(n) > tarr(n - 1)) count = count + 1
319  end do
320  !
321  ! -- allocate b for unique values
322  if (allocated(b)) then
323  deallocate (b)
324  end if
325  allocate (b(count))
326  !
327  ! -- fill b with unique values
328  b(1) = tarr(1)
329  count = 1
330  do n = 2, size(a)
331  if (tarr(n) > b(count)) then
332  count = count + 1
333  b(count) = tarr(n)
334  end if
335  end do
336  !
337  ! -- allocate tarr and create idxarr
338  deallocate (tarr)
339  deallocate (indxarr)
340  !
341  ! -- return
342  return