Octopus
sort.F90
Go to the documentation of this file.
1!! Copyright (C) 2002-2006 M. Marques, A. Castro, A. Rubio, G. Bertsch
2!! Copyright (C) 2021 S. Ohlmann
3!!
4!! This program is free software; you can redistribute it and/or modify
5!! it under the terms of the GNU General Public License as published by
6!! the Free Software Foundation; either version 2, or (at your option)
7!! any later version.
8!!
9!! This program is distributed in the hope that it will be useful,
10!! but WITHOUT ANY WARRANTY; without even the implied warranty of
11!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12!! GNU General Public License for more details.
13!!
14!! You should have received a copy of the GNU General Public License
15!! along with this program; if not, write to the Free Software
16!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
17!! 02110-1301, USA.
18!!
19
20#include "global.h"
21
24module sort_oct_m
25 use debug_oct_m
26 use global_oct_m
27
28 implicit none
29
30 private
31 public :: &
32 sort, &
34
35 ! ------------------------------------------------------------------------------
57 interface sort
58 module procedure dsort, isort, lsort
59 module procedure dshellsort1, zshellsort1, ishellsort1
60 module procedure dshellsort2, zshellsort2, ishellsort2
61 end interface sort
62
64 interface
65 subroutine dsort1(size, array)
66 use, intrinsic :: iso_fortran_env
67 implicit none
68 integer, intent(in) :: size
69 real(real64), intent(inout) :: array(*)
70 end subroutine dsort1
71
72 subroutine dsort2(size, array, indices)
73 use, intrinsic :: iso_fortran_env
74 implicit none
75 integer, intent(in) :: size
76 real(real64), intent(inout) :: array(*)
77 integer, intent(out) :: indices(*)
78 end subroutine dsort2
79
80 subroutine isort1(size, array)
81 implicit none
82 integer, intent(in) :: size
83 integer, intent(inout) :: array(*)
84 end subroutine isort1
85
86 subroutine isort2(size, array, indices)
87 implicit none
88 integer, intent(in) :: size
89 integer, intent(inout) :: array(*)
90 integer, intent(out) :: indices(*)
91 end subroutine isort2
92
93 subroutine lsort1(size, array)
94 use, intrinsic :: iso_fortran_env
95 implicit none
96 integer, intent(in) :: size
97 integer(int64), intent(inout) :: array(*)
98 end subroutine lsort1
99
100 subroutine lsort2(size, array, indices)
101 use, intrinsic :: iso_fortran_env
102 implicit none
103 integer, intent(in) :: size
104 integer(int64), intent(inout) :: array(*)
105 integer, intent(out) :: indices(*)
106 end subroutine lsort2
107 end interface
108
109contains
110
111 ! ---------------------------------------------------------
112 subroutine dsort(a, ind)
113 real(real64), intent(inout) :: a(:)
114 integer, optional, intent(out) :: ind(:)
115
116 push_sub(dsort)
117
118 if (size(a) > 0) then
120 if (.not. present(ind)) then
121 call dsort1(size(a), a)
122 else
123 call dsort2(size(a), a, ind)
124 end if
125
126 end if
127
128 pop_sub(dsort)
129 end subroutine dsort
130
131
132 ! ---------------------------------------------------------
134 subroutine isort(a, ind)
135 integer, intent(inout) :: a(:)
136 integer, optional, intent(out) :: ind(:)
137
138 push_sub(isort)
139
140 if (size(a) > 0) then
141
142 if (.not. present(ind)) then
143 call isort1(size(a), a)
144 else
145 call isort2(size(a), a, ind)
146 end if
147
148 end if
149
150 pop_sub(isort)
151 end subroutine isort
153 ! ---------------------------------------------------------
155 subroutine lsort(a, ind)
156 integer(int64), intent(inout) :: a(:)
157 integer, optional, intent(out) :: ind(:)
158
159 push_sub(lsort)
161 if (size(a) > 0) then
162
163 if (.not. present(ind)) then
164 call lsort1(size(a), a)
165 else
166 call lsort2(size(a), a, ind)
167 end if
168
169 end if
170
171 pop_sub(lsort)
172 end subroutine lsort
173
174
176 pure logical function less_idx(i, j, off, kabs, ksgn) result(less)
177 integer, intent(in) :: i, j
178 integer, intent(in) :: off(:, :) ! (ndim, n)
179 real(int64), intent(in) :: kabs(:)
180 integer, intent(in) :: ksgn(:)
182 integer :: d, ndim
183
184 if (.not. abs(kabs(i) - kabs(j)) <= abs(kabs(j) * 1.0e-14_real64)) then
185 less = kabs(i) < kabs(j)
186 return
187 end if
189 if (ksgn(i) /= ksgn(j)) then
190 less = ksgn(i) < ksgn(j) ! 0 before 1
191 return
192 end if
193
194 ndim = size(off, 1)
195 do d = 1, ndim
196 if (off(d,i) /= off(d,j)) then
197 less = off(d,i) < off(d,j)
198 return
199 end if
200 end do
201
202 less = i < j ! final stability (should not matter if offsets are unique)
203 end function less_idx
204
206 recursive subroutine mergesort_perm(perm, tmp, l, r, off, kabs, ksgn)
207 integer, intent(inout) :: perm(:), tmp(:)
208 integer, intent(in) :: l, r
209 integer, intent(in) :: off(:, :)
210 real(real64), intent(in) :: kabs(:)
211 integer, intent(in) :: ksgn(:)
212
213 integer :: m, i, j, k
214
215 if (r <= l) return
216 m = (l + r) / 2
217
218 call mergesort_perm(perm, tmp, l, m, off, kabs, ksgn)
219 call mergesort_perm(perm, tmp, m+1, r, off, kabs, ksgn)
220
221 i = l; j = m+1; k = l
222 do while (i <= m .and. j <= r)
223 if (less_idx(perm(i), perm(j), off, kabs, ksgn)) then
224 tmp(k) = perm(i); i = i+1
225 else
226 tmp(k) = perm(j); j = j+1
227 end if
228 k = k+1
229 end do
230 do while (i <= m); tmp(k) = perm(i); i=i+1; k=k+1; end do
231 do while (j <= r); tmp(k) = perm(j); j=j+1; k=k+1; end do
232
233 perm(l:r) = tmp(l:r)
234 end subroutine mergesort_perm
235
242 subroutine robust_sort_by_abs(v, off, perm, negative_first)
243 real(real64), intent(in) :: v(:)
244 integer, intent(in) :: off(:, :) ! (ndim, n)
245 integer, intent(out) :: perm(size(v))
246 logical, intent(in), optional :: negative_first
247
248 integer :: n, i
249 logical :: neg_first
250 integer, allocatable :: tmp(:)
251 real(real64), allocatable :: kabs(:)
252 integer, allocatable :: ksgn(:)
253
254 push_sub(robust_sort_by_abs)
255
256 n = size(v)
257 assert(size(off, dim=2) == n)
258
259 neg_first = optional_default(negative_first, .true.)
260
261 allocate(tmp(n), kabs(n), ksgn(n))
262
263 do i = 1, n
264 perm(i) = i
265 kabs(i) = abs(v(i))
266
267 if (neg_first) then
268 ksgn(i) = merge(0, 1, v(i) < 0.0_real64) ! neg (0) then pos (1)
269 else
270 ksgn(i) = merge(0, 1, v(i) >= 0.0_real64) ! pos (0) then neg (1)
271 end if
272 end do
273
274 call mergesort_perm(perm, tmp, 1, n, off, kabs, ksgn)
275
276 deallocate(tmp, kabs, ksgn)
277 pop_sub(robust_sort_by_abs)
278 end subroutine robust_sort_by_abs
279
280#include "undef.F90"
281#include "complex.F90"
282#include "sort_inc.F90"
283
284#include "undef.F90"
285#include "real.F90"
286#include "sort_inc.F90"
287
288#include "undef.F90"
289#include "integer.F90"
290#include "sort_inc.F90"
291
292end module sort_oct_m
293
294!! Local Variables:
295!! mode: f90
296!! coding: utf-8
297!! End:
from sort_low.cc
Definition: sort.F90:160
This is the common interface to a sorting routine. It performs the shell algorithm,...
Definition: sort.F90:152
This module is intended to contain "only mathematical" functions and procedures.
Definition: sort.F90:119
subroutine ishellsort2(a, x)
Definition: sort.F90:822
subroutine dsort(a, ind)
Definition: sort.F90:208
subroutine dshellsort1(a, x)
Definition: sort.F90:610
subroutine zshellsort2(a, x)
Definition: sort.F90:490
subroutine isort(a, ind)
Shell sort for integer arrays.
Definition: sort.F90:230
subroutine ishellsort1(a, x)
Definition: sort.F90:776
subroutine dshellsort2(a, x)
Definition: sort.F90:656
recursive subroutine mergesort_perm(perm, tmp, l, r, off, kabs, ksgn)
Perform the permutations for the sorting.
Definition: sort.F90:302
pure logical function less_idx(i, j, off, kabs, ksgn)
Sorting criterium for the robust sorting below.
Definition: sort.F90:272
subroutine lsort(a, ind)
Shell sort for integer(int64) arrays.
Definition: sort.F90:251
subroutine, public robust_sort_by_abs(v, off, perm, negative_first)
Robbust sorting of floating point numbers by absolute values.
Definition: sort.F90:338
subroutine zshellsort1(a, x)
Definition: sort.F90:444
int true(void)