7submodule(
linalg) linalg_sorting
12 module subroutine sort_dbl_array(x, ascend)
14 real(real64),
intent(inout),
dimension(:) :: x
15 logical,
intent(in),
optional :: ascend
19 integer(int32) :: n, info
22 if (
present(ascend))
then
34 call dlasrt(id, n, x, info)
38 module subroutine sort_dbl_array_ind(x, ind, ascend, err)
40 real(real64),
intent(inout),
dimension(:) :: x
41 integer(int32),
intent(inout),
dimension(:) :: ind
42 logical,
intent(in),
optional :: ascend
43 class(errors),
intent(inout),
optional,
target :: err
46 class(errors),
pointer :: errmgr
47 type(errors),
target :: deferr
48 character(len = 128) :: errmsg
54 if (
present(err))
then
59 if (
present(ascend))
then
66 if (
size(ind) /= n)
then
68 "Expected the tracking array to be of size ", n, &
69 ", but found an array of size ",
size(ind),
"."
70 call errmgr%report_error(
"sort_dbl_array_ind", trim(errmsg), &
77 call qsort_dbl_ind(dir, x, ind)
80100
format(a, i0, a, i0, a)
84 module subroutine sort_cmplx_array(x, ascend)
86 complex(real64),
intent(inout),
dimension(:) :: x
87 logical,
intent(in),
optional :: ascend
93 if (
present(ascend))
then
100 call qsort_cmplx(dir, x)
104 module subroutine sort_cmplx_array_ind(x, ind, ascend, err)
106 complex(real64),
intent(inout),
dimension(:) :: x
107 integer(int32),
intent(inout),
dimension(:) :: ind
108 logical,
intent(in),
optional :: ascend
109 class(errors),
intent(inout),
optional,
target :: err
112 class(errors),
pointer :: errmgr
113 type(errors),
target :: deferr
114 character(len = 128) :: errmsg
120 if (
present(err))
then
125 if (
present(ascend))
then
132 if (
size(ind) /= n)
then
134 "Expected the tracking array to be of size ", n, &
135 ", but found an array of size ",
size(ind),
"."
136 call errmgr%report_error(
"sort_cmplx_array_ind", trim(errmsg), &
143 call qsort_cmplx_ind(dir, x, ind)
146100
format(a, i0, a, i0, a)
150 module subroutine sort_eigen_cmplx(vals, vecs, ascend, err)
152 complex(real64),
intent(inout),
dimension(:) :: vals
153 complex(real64),
intent(inout),
dimension(:,:) :: vecs
154 logical,
intent(in),
optional :: ascend
155 class(errors),
intent(inout),
optional,
target :: err
158 class(errors),
pointer :: errmgr
159 type(errors),
target :: deferr
160 character(len = 128) :: errmsg
161 integer(int32) :: i, n, flag
163 integer(int32),
allocatable,
dimension(:) :: ind
166 if (
present(err))
then
171 if (
present(ascend))
then
179 if (
size(vecs, 1) /= n .or.
size(vecs, 2) /= n)
then
182 "Expected the eigenvector matrix to be of size ", n, &
183 "-by-", n,
", but found a matrix of size ",
size(vecs, 1), &
184 "-by-",
size(vecs, 2),
"."
185 call errmgr%report_error(
"sort_eigen_cmplx", trim(errmsg), &
190 allocate(ind(n), stat = flag)
192 call errmgr%report_error(
"sort_eigen_cmplx", &
193 "Insufficient memory available.", la_out_of_memory_error)
201 call qsort_cmplx_ind(dir, vals, ind)
208100
format(a, i0, a, i0, a, i0, a, i0, a)
212 module subroutine sort_eigen_dbl(vals, vecs, ascend, err)
214 real(real64),
intent(inout),
dimension(:) :: vals
215 real(real64),
intent(inout),
dimension(:,:) :: vecs
216 logical,
intent(in),
optional :: ascend
217 class(errors),
intent(inout),
optional,
target :: err
220 class(errors),
pointer :: errmgr
221 type(errors),
target :: deferr
222 character(len = 128) :: errmsg
223 integer(int32) :: i, n, flag
225 integer(int32),
allocatable,
dimension(:) :: ind
228 if (
present(err))
then
233 if (
present(ascend))
then
241 if (
size(vecs, 1) /= n .or.
size(vecs, 2) /= n)
then
244 "Expected the eigenvector matrix to be of size ", n, &
245 "-by-", n,
", but found a matrix of size ",
size(vecs, 1), &
246 "-by-",
size(vecs, 2),
"."
247 call errmgr%report_error(
"sort_eigen_dbl", trim(errmsg), &
252 allocate(ind(n), stat = flag)
254 call errmgr%report_error(
"sort_eigen_dbl", &
255 "Insufficient memory available.", la_out_of_memory_error)
263 call qsort_dbl_ind(dir, vals, ind)
270100
format(a, i0, a, i0, a, i0, a, i0, a)
289 recursive subroutine qsort_dbl_ind(ascend, x, ind)
291 logical,
intent(in) :: ascend
292 real(real64),
intent(inout),
dimension(:) :: x
293 integer(int32),
intent(inout),
dimension(:) :: ind
299 if (
size(x) > 1)
then
300 call dbl_partition_ind(ascend, x, ind, iq)
301 call qsort_dbl_ind(ascend, x(:iq-1), ind(:iq-1))
302 call qsort_dbl_ind(ascend, x(iq:), ind(iq:))
322 subroutine dbl_partition_ind(ascend, x, ind, marker)
324 logical,
intent(in) :: ascend
325 real(real64),
intent(inout),
dimension(:) :: x
326 integer(int32),
intent(inout),
dimension(:) :: ind
327 integer(int32),
intent(out) :: marker
330 integer(int32) :: i, j, itemp
331 real(real64) :: temp, pivot
342 if (x(j) <= pivot)
exit
347 if (x(i) >= pivot)
exit
359 else if (i == j)
then
372 if (x(j) >= pivot)
exit
377 if (x(i) <= pivot)
exit
389 else if (i == j)
then
415 recursive subroutine qsort_cmplx(ascend, x)
417 logical,
intent(in) :: ascend
418 complex(real64),
intent(inout),
dimension(:) :: x
424 if (
size(x) > 1)
then
425 call cmplx_partition(ascend, x, iq)
426 call qsort_cmplx(ascend, x(:iq-1))
427 call qsort_cmplx(ascend, x(iq:))
448 subroutine cmplx_partition(ascend, x, marker)
450 logical,
intent(in) :: ascend
451 complex(real64),
intent(inout),
dimension(:) :: x
452 integer(int32),
intent(out) :: marker
455 integer(int32) :: i, j
456 complex(real64) :: temp
457 real(real64) :: pivot
460 pivot = real(x(1), real64)
468 if (real(x(j), real64) <= pivot)
exit
473 if (real(x(i), real64) >= pivot)
exit
481 else if (i == j)
then
494 if (real(x(j), real64) >= pivot)
exit
499 if (real(x(i), real64) <= pivot)
exit
507 else if (i == j)
then
536 recursive subroutine qsort_cmplx_ind(ascend, x, ind)
538 logical,
intent(in) :: ascend
539 complex(real64),
intent(inout),
dimension(:) :: x
540 integer(int32),
intent(inout),
dimension(:) :: ind
546 if (
size(x) > 1)
then
547 call cmplx_partition_ind(ascend, x, ind, iq)
548 call qsort_cmplx_ind(ascend, x(:iq-1), ind(:iq-1))
549 call qsort_cmplx_ind(ascend, x(iq:), ind(iq:))
573 subroutine cmplx_partition_ind(ascend, x, ind, marker)
575 logical,
intent(in) :: ascend
576 complex(real64),
intent(inout),
dimension(:) :: x
577 integer(int32),
intent(inout),
dimension(:) :: ind
578 integer(int32),
intent(out) :: marker
581 integer(int32) :: i, j, itemp
582 complex(real64) :: temp
583 real(real64) :: pivot
586 pivot = real(x(1), real64)
594 if (real(x(j), real64) <= pivot)
exit
599 if (real(x(i), real64) >= pivot)
exit
611 else if (i == j)
then
624 if (real(x(j), real64) >= pivot)
exit
629 if (real(x(i), real64) <= pivot)
exit
641 else if (i == j)
then
Provides a set of common linear algebra routines.