SCALE-RM
Data Types | Functions/Subroutines
scale_sort Module Reference

module SORT More...

Functions/Subroutines

subroutine, public sort_uniq_int_sorted (n, ary, c)
 
recursive subroutine, public sort_quicksort (n, array)
 
subroutine, public sort_heapsort (n, array)
 
recursive subroutine, public sort_quickselect (A, left, right, K)
 
recursive subroutine, public sort_quickselect_arg (A, X, left, right, K)
 
recursive subroutine, public sort_quickselect_desc (A, left, right, K)
 
recursive subroutine, public sort_quickselect_desc_arg (A, X, left, right, K)
 

Detailed Description

module SORT

Description
Sort data
Author
Team SCALE

Function/Subroutine Documentation

◆ sort_uniq_int_sorted()

subroutine, public scale_sort::sort_uniq_int_sorted ( integer(dp), intent(in)  n,
integer(dp), dimension(n), intent(inout)  ary,
integer(dp), intent(out)  c 
)

Definition at line 200 of file scale_sort.F90.

200  integer(DP), intent(in) :: n
201  integer(DP), intent(inout) :: ary(n)
202  integer(DP), intent(out) :: c
203  integer(DP) :: i, j
204 
205  c = 1
206  j = ary(1)
207  do i = 2, n
208  if(j .ne. ary(i)) then
209  j = ary(i)
210  c = c + 1
211  ary(c) = ary(i)
212  end if
213  end do
214 
215  return

Referenced by scale_letkf::radar_superobing().

Here is the caller graph for this function:

◆ sort_quicksort()

recursive subroutine, public scale_sort::sort_quicksort ( integer(dp), intent(in)  n,
integer(dp), dimension(n), intent(inout)  array 
)

Definition at line 219 of file scale_sort.F90.

219  implicit none
220  integer(DP), intent(in) :: n
221  integer(DP), intent(inout) :: array(n)
222  integer(DP) :: i, j, k, k1, kn2, kn, tmp
223  integer, parameter :: threshold = 4
224 
225  if(n <= 1) return
226 
227  k1 = array(1)
228  kn2 = array(n / 2)
229  kn = array(n)
230  if(k1 < kn2) then
231  if(kn2 < kn) then
232  k = kn2 ! 1, n / 2, n
233  else ! array(n) <= array(n / 2)
234  if(k1 < kn) then
235  k = kn ! 1, n, n / 2
236  else
237  k = k1 ! n, 1, n / 2
238  end if
239  end if
240  else ! array(n / 2) <= array(1)
241  if(k1 < kn) then
242  k = k1 ! n / 2, 1, n
243  else ! array(n) <= array(1)
244  if(kn2 < kn) then
245  k = kn ! n / 2, n, 1
246  else ! array(n) <= array(n / 2)
247  k = kn2 ! n, n / 2, 1
248  end if
249  end if
250  end if
251 
252  i = 1
253  j = n
254  do
255  do i = i, n, 1
256  if(array(i) .ge. k) exit
257  end do
258  do j = j, 1, -1
259  if(array(j) .le. k) exit
260  end do
261  if(i >= j) exit
262  tmp = array(i)
263  array(i) = array(j)
264  array(j) = tmp
265  i = i + 1
266  j = j - 1
267  end do
268 
269  if(i - 1 > threshold) then
270  call sort_quicksort(i - 1, array(1:(i - 1)))
271  else if(i - 1 > 1) then
272  call sort_heapsort(i - 1, array(1:(i - 1)))
273  end if
274  if(n - j > threshold) then
275  call sort_quicksort(n - j, array((j + 1):n))
276  else if(n - j > 1) then
277  call sort_heapsort(n - j, array((j + 1):n))
278  end if
279 
280  return

References sort_heapsort().

Referenced by com_gamma().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ sort_heapsort()

subroutine, public scale_sort::sort_heapsort ( integer(dp), intent(in)  n,
integer(dp), dimension(n), intent(inout)  array 
)

Definition at line 284 of file scale_sort.F90.

284  implicit none
285  integer(DP), intent(in) :: n
286  integer(DP), intent(inout) :: array(n)
287 
288  integer(DP) :: i, j, k, l
289  integer(DP) :: t
290 
291  l = n / 2 + 1
292  k = n
293  do while(k .ne. 1)
294  if(l .gt. 1) then
295  l = l - 1
296  t = array(l)
297  else
298  t = array(k)
299  array(k) = array(1)
300  k = k - 1
301  if(k .eq. 1) then
302  array(1) = t
303  exit
304  end if
305  end if
306  i = l
307  j = l + l
308  do while(j .le. k)
309  if(j .lt. k)then
310  if(array(j) .lt. array(j + 1)) j = j + 1
311  end if
312  if(t .lt. array(j))then
313  array(i) = array(j)
314  i = j
315  j = j + j
316  else
317  j = k + 1
318  endif
319  enddo
320  array(i) = t
321  end do
322 
323  return

Referenced by sort_quicksort().

Here is the caller graph for this function:

◆ sort_quickselect()

recursive subroutine, public scale_sort::sort_quickselect ( real(rp), dimension(:), intent(inout)  A,
integer, intent(in)  left,
integer, intent(in)  right,
integer, intent(in)  K 
)

Definition at line 327 of file scale_sort.F90.

327  implicit none
328  real(RP), intent(inout) :: A(:)
329  integer, intent(in) :: left, right
330  integer, intent(in) :: K
331  integer :: middle, pivot
332 
333  if (left < right) then
334  if ((right-left)/k >= 2) then
335  call sample_second_min(a, left, right, k, pivot)
336  else
337  middle = (left + right) / 2
338  call median_of_three(a, left, middle, right, pivot)
339  end if
340  call partition(a, left, right, pivot)
341  if (k < pivot) then
342  call sort_quickselect(a, left, pivot-1, k)
343  else if (k > pivot) then
344  call sort_quickselect(a, pivot+1, right, k)
345  end if
346  end if
347 
348  return

◆ sort_quickselect_arg()

recursive subroutine, public scale_sort::sort_quickselect_arg ( real(rp), dimension(:), intent(in)  A,
integer, dimension(:), intent(inout)  X,
integer, intent(in)  left,
integer, intent(in)  right,
integer, intent(in)  K 
)

Definition at line 352 of file scale_sort.F90.

352  implicit none
353  real(RP), intent(in) :: A(:)
354  integer, intent(inout) :: X(:)
355  integer, intent(in) :: left, right
356  integer, intent(in) :: K
357  integer :: middle, pivot
358 
359  if (left < right) then
360  if ((right-left)/k >= 2) then
361  call sample_second_min_arg(a, x, left, right, k, pivot)
362  else
363  middle = (left + right) / 2
364  call median_of_three_arg(a, x, left, middle, right, pivot)
365  end if
366  call partition_arg(a, x, left, right, pivot)
367  if (k < pivot) then
368  call sort_quickselect_arg(a, x, left, pivot-1, k)
369  else if (k > pivot) then
370  call sort_quickselect_arg(a, x, pivot+1, right, k)
371  end if
372  end if
373 
374  return

Referenced by com_gamma(), and scale_letkf::letkf_add_inflation_setup().

Here is the caller graph for this function:

◆ sort_quickselect_desc()

recursive subroutine, public scale_sort::sort_quickselect_desc ( real(rp), dimension(:), intent(inout)  A,
integer, intent(in)  left,
integer, intent(in)  right,
integer, intent(in)  K 
)

Definition at line 378 of file scale_sort.F90.

378  implicit none
379  real(RP), intent(inout) :: A(:)
380  integer, intent(in) :: left, right
381  integer, intent(in) :: K
382  integer :: middle, pivot
383 
384  if (left < right) then
385  if ((right-left)/k >= 2) then
386  call sample_second_max(a, left, right, k, pivot)
387  else
388  middle = (left + right) / 2
389  call median_of_three(a, left, middle, right, pivot)
390  end if
391  call partition_desc(a, left, right, pivot)
392  if (k < pivot) then
393  call sort_quickselect_desc(a, left, pivot-1, k)
394  else if (k > pivot) then
395  call sort_quickselect_desc(a, pivot+1, right, k)
396  end if
397  end if
398 
399  return

◆ sort_quickselect_desc_arg()

recursive subroutine, public scale_sort::sort_quickselect_desc_arg ( real(rp), dimension(:), intent(in)  A,
integer, dimension(:), intent(inout)  X,
integer, intent(in)  left,
integer, intent(in)  right,
integer, intent(in)  K 
)

Definition at line 403 of file scale_sort.F90.

403  implicit none
404  real(RP), intent(in) :: A(:)
405  integer, intent(inout) :: X(:)
406  integer, intent(in) :: left, right
407  integer, intent(in) :: K
408  integer :: middle, pivot
409 
410  if (left < right) then
411  if ((right-left)/k >= 2) then
412  call sample_second_max_arg(a, x, left, right, k, pivot)
413  else
414  middle = (left + right) / 2
415  call median_of_three_arg(a, x, left, middle, right, pivot)
416  end if
417  call partition_desc_arg(a, x, left, right, pivot)
418  if (k < pivot) then
419  call sort_quickselect_desc_arg(a, x, left, pivot-1, k)
420  else if (k > pivot) then
421  call sort_quickselect_desc_arg(a, x, pivot+1, right, k)
422  end if
423  end if
424 
425  return

Referenced by com_gamma(), and scale_letkf::letkf_add_inflation_setup().

Here is the caller graph for this function: