MODFLOW 6  version 6.6.0.dev0
USGS Modular Hydrologic Model
sortmodule::qsort Interface Reference
Collaboration diagram for sortmodule::qsort:
Collaboration graph

Private Member Functions

subroutine qsort_int1d (indx, v, reverse)
 
subroutine qsort_dbl1d (indx, v, reverse)
 

Detailed Description

Definition at line 11 of file sort.f90.

Member Function/Subroutine Documentation

◆ qsort_dbl1d()

subroutine sortmodule::qsort::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
Here is the call graph for this function:

◆ qsort_int1d()

subroutine sortmodule::qsort::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
Here is the call graph for this function:

The documentation for this interface was generated from the following file: