linalg 1.7.4
A linear algebra library that provides a user-friendly interface to several BLAS and LAPACK routines.
Loading...
Searching...
No Matches
linalg_sorting.f90
1! linalg_sorting.f90
2
7submodule(linalg) linalg_sorting
8 use lapack
9 implicit none
10contains
11! ******************************************************************************
12! SORTING ROUTINES
13! ------------------------------------------------------------------------------
14 module subroutine sort_dbl_array(x, ascend)
15 ! Arguments
16 real(real64), intent(inout), dimension(:) :: x
17 logical, intent(in), optional :: ascend
18
19 ! Local Variables
20 character :: id
21 integer(int32) :: n, info
22
23 ! Initialization
24 if (present(ascend)) then
25 if (ascend) then
26 id = 'I'
27 else
28 id = 'D'
29 end if
30 else
31 id = 'I'
32 end if
33 n = size(x)
34
35 ! Process
36 call dlasrt(id, n, x, info)
37 end subroutine
38
39! ------------------------------------------------------------------------------
40 module subroutine sort_dbl_array_ind(x, ind, ascend, err)
41 ! Arguments
42 real(real64), intent(inout), dimension(:) :: x
43 integer(int32), intent(inout), dimension(:) :: ind
44 logical, intent(in), optional :: ascend
45 class(errors), intent(inout), optional, target :: err
46
47 ! Local Variables
48 class(errors), pointer :: errmgr
49 type(errors), target :: deferr
50 character(len = :), allocatable :: errmsg
51 integer(int32) :: n
52 logical :: dir
53
54 ! Initialization
55 n = size(x)
56 if (present(err)) then
57 errmgr => err
58 else
59 errmgr => deferr
60 end if
61 if (present(ascend)) then
62 dir = ascend
63 else
64 dir = .true. ! Ascend == true
65 end if
66
67 ! Input Check
68 if (size(ind) /= n) then
69 allocate(character(len = 256) :: errmsg)
70 write(errmsg, 100) &
71 "Expected the tracking array to be of size ", n, &
72 ", but found an array of size ", size(ind), "."
73 call errmgr%report_error("sort_dbl_array_ind", trim(errmsg), &
74 la_array_size_error)
75 return
76 end if
77 if (n <= 1) return
78
79 ! Process
80 call qsort_dbl_ind(dir, x, ind)
81
82 ! Formatting
83100 format(a, i0, a, i0, a)
84 end subroutine
85
86! ------------------------------------------------------------------------------
87 module subroutine sort_cmplx_array(x, ascend)
88 ! Arguments
89 complex(real64), intent(inout), dimension(:) :: x
90 logical, intent(in), optional :: ascend
91
92 ! Local Variables
93 logical :: dir
94
95 ! Initialization
96 if (present(ascend)) then
97 dir = ascend
98 else
99 dir = .true.
100 end if
101
102 ! Process
103 call qsort_cmplx(dir, x)
104 end subroutine
105
106! ------------------------------------------------------------------------------
107 module subroutine sort_cmplx_array_ind(x, ind, ascend, err)
108 ! Arguments
109 complex(real64), intent(inout), dimension(:) :: x
110 integer(int32), intent(inout), dimension(:) :: ind
111 logical, intent(in), optional :: ascend
112 class(errors), intent(inout), optional, target :: err
113
114 ! Local Variables
115 class(errors), pointer :: errmgr
116 type(errors), target :: deferr
117 character(len = :), allocatable :: errmsg
118 integer(int32) :: n
119 logical :: dir
120
121 ! Initialization
122 n = size(x)
123 if (present(err)) then
124 errmgr => err
125 else
126 errmgr => deferr
127 end if
128 if (present(ascend)) then
129 dir = ascend
130 else
131 dir = .true. ! Ascend == true
132 end if
133
134 ! Input Check
135 if (size(ind) /= n) then
136 allocate(character(len = 256) :: errmsg)
137 write(errmsg, 100) &
138 "Expected the tracking array to be of size ", n, &
139 ", but found an array of size ", size(ind), "."
140 call errmgr%report_error("sort_cmplx_array_ind", trim(errmsg), &
141 la_array_size_error)
142 return
143 end if
144 if (n <= 1) return
145
146 ! Process
147 call qsort_cmplx_ind(dir, x, ind)
148
149 ! Formatting
150100 format(a, i0, a, i0, a)
151 end subroutine
152
153! ------------------------------------------------------------------------------
154 module subroutine sort_eigen_cmplx(vals, vecs, ascend, err)
155 ! Arguments
156 complex(real64), intent(inout), dimension(:) :: vals
157 complex(real64), intent(inout), dimension(:,:) :: vecs
158 logical, intent(in), optional :: ascend
159 class(errors), intent(inout), optional, target :: err
160
161 ! Local Variables
162 class(errors), pointer :: errmgr
163 type(errors), target :: deferr
164 character(len = :), allocatable :: errmsg
165 integer(int32) :: i, n, flag
166 logical :: dir
167 integer(int32), allocatable, dimension(:) :: ind
168
169 ! Initialization
170 if (present(err)) then
171 errmgr => err
172 else
173 errmgr => deferr
174 end if
175 if (present(ascend)) then
176 dir = ascend
177 else
178 dir = .true. ! Ascend == true
179 end if
180
181 ! Ensure the eigenvector matrix is sized appropriately
182 n = size(vals)
183 if (size(vecs, 1) /= n .or. size(vecs, 2) /= n) then
184 ! ARRAY SIZE ERROR
185 allocate(character(len = 256) :: errmsg)
186 write(errmsg, 100) &
187 "Expected the eigenvector matrix to be of size ", n, &
188 "-by-", n, ", but found a matrix of size ", size(vecs, 1), &
189 "-by-", size(vecs, 2), "."
190 call errmgr%report_error("sort_eigen_cmplx", trim(errmsg), &
191 la_array_size_error)
192 end if
193
194 ! Allocate memory for the tracking array
195 allocate(ind(n), stat = flag)
196 if (flag /= 0) then
197 call errmgr%report_error("sort_eigen_cmplx", &
198 "Insufficient memory available.", la_out_of_memory_error)
199 return
200 end if
201 do i = 1, n
202 ind(i) = i
203 end do
204
205 ! Sort
206 call qsort_cmplx_ind(dir, vals, ind)
207
208 ! Shift the eigenvectors around to keep them associated with the
209 ! appropriate eigenvalue
210 vecs = vecs(:,ind)
211
212 ! Formatting
213100 format(a, i0, a, i0, a, i0, a, i0, a)
214 end subroutine
215
216! ------------------------------------------------------------------------------
217 module subroutine sort_eigen_dbl(vals, vecs, ascend, err)
218 ! Arguments
219 real(real64), intent(inout), dimension(:) :: vals
220 real(real64), intent(inout), dimension(:,:) :: vecs
221 logical, intent(in), optional :: ascend
222 class(errors), intent(inout), optional, target :: err
223
224 ! Local Variables
225 class(errors), pointer :: errmgr
226 type(errors), target :: deferr
227 character(len = :), allocatable :: errmsg
228 integer(int32) :: i, n, flag
229 logical :: dir
230 integer(int32), allocatable, dimension(:) :: ind
231
232 ! Initialization
233 if (present(err)) then
234 errmgr => err
235 else
236 errmgr => deferr
237 end if
238 if (present(ascend)) then
239 dir = ascend
240 else
241 dir = .true. ! Ascend == true
242 end if
243
244 ! Ensure the eigenvector matrix is sized appropriately
245 n = size(vals)
246 if (size(vecs, 1) /= n .or. size(vecs, 2) /= n) then
247 ! ARRAY SIZE ERROR
248 allocate(character(len = 256) :: errmsg)
249 write(errmsg, 100) &
250 "Expected the eigenvector matrix to be of size ", n, &
251 "-by-", n, ", but found a matrix of size ", size(vecs, 1), &
252 "-by-", size(vecs, 2), "."
253 call errmgr%report_error("sort_eigen_dbl", trim(errmsg), &
254 la_array_size_error)
255 end if
256
257 ! Allocate memory for the tracking array
258 allocate(ind(n), stat = flag)
259 if (flag /= 0) then
260 call errmgr%report_error("sort_eigen_dbl", &
261 "Insufficient memory available.", la_out_of_memory_error)
262 return
263 end if
264 do i = 1, n
265 ind(i) = i
266 end do
267
268 ! Sort
269 call qsort_dbl_ind(dir, vals, ind)
270
271 ! Shift the eigenvectors around to keep them associated with the
272 ! appropriate eigenvalue
273 vecs = vecs(:,ind)
274
275 ! Formatting
276100 format(a, i0, a, i0, a, i0, a, i0, a)
277 end subroutine
278
279! ******************************************************************************
280! PRIVATE HELPER ROUTINES
281! ------------------------------------------------------------------------------
295 recursive subroutine qsort_dbl_ind(ascend, x, ind)
296 ! Arguments
297 logical, intent(in) :: ascend
298 real(real64), intent(inout), dimension(:) :: x
299 integer(int32), intent(inout), dimension(:) :: ind
300
301 ! Local Variables
302 integer(int32) :: iq
303
304 ! Process
305 if (size(x) > 1) then
306 call dbl_partition_ind(ascend, x, ind, iq)
307 call qsort_dbl_ind(ascend, x(:iq-1), ind(:iq-1))
308 call qsort_dbl_ind(ascend, x(iq:), ind(iq:))
309 end if
310 end subroutine
311
312! ------------------------------------------------------------------------------
328 subroutine dbl_partition_ind(ascend, x, ind, marker)
329 ! Arguments
330 logical, intent(in) :: ascend
331 real(real64), intent(inout), dimension(:) :: x
332 integer(int32), intent(inout), dimension(:) :: ind
333 integer(int32), intent(out) :: marker
334
335 ! Local Variables
336 integer(int32) :: i, j, itemp
337 real(real64) :: temp, pivot
338
339 ! Process
340 pivot = x(1)
341 i = 0
342 j = size(x) + 1
343 if (ascend) then
344 ! Ascending Sort
345 do
346 j = j - 1
347 do
348 if (x(j) <= pivot) exit
349 j = j - 1
350 end do
351 i = i + 1
352 do
353 if (x(i) >= pivot) exit
354 i = i + 1
355 end do
356 if (i < j) then
357 ! Exchage X(I) and X(J)
358 temp = x(i)
359 x(i) = x(j)
360 x(j) = temp
361
362 itemp = ind(i)
363 ind(i) = ind(j)
364 ind(j) = itemp
365 else if (i == j) then
366 marker = i + 1
367 return
368 else
369 marker = i
370 return
371 end if
372 end do
373 else
374 ! Descending Sort
375 do
376 j = j - 1
377 do
378 if (x(j) >= pivot) exit
379 j = j - 1
380 end do
381 i = i + 1
382 do
383 if (x(i) <= pivot) exit
384 i = i + 1
385 end do
386 if (i < j) then
387 ! Exchage X(I) and X(J)
388 temp = x(i)
389 x(i) = x(j)
390 x(j) = temp
391
392 itemp = ind(i)
393 ind(i) = ind(j)
394 ind(j) = itemp
395 else if (i == j) then
396 marker = i + 1
397 return
398 else
399 marker = i
400 return
401 end if
402 end do
403 end if
404 end subroutine
405
406! ------------------------------------------------------------------------------
421 recursive subroutine qsort_cmplx(ascend, x)
422 ! Arguments
423 logical, intent(in) :: ascend
424 complex(real64), intent(inout), dimension(:) :: x
425
426 ! Local Variables
427 integer(int32) :: iq
428
429 ! Process
430 if (size(x) > 1) then
431 call cmplx_partition(ascend, x, iq)
432 call qsort_cmplx(ascend, x(:iq-1))
433 call qsort_cmplx(ascend, x(iq:))
434 end if
435 end subroutine
436
437! ------------------------------------------------------------------------------
454 subroutine cmplx_partition(ascend, x, marker)
455 ! Arguments
456 logical, intent(in) :: ascend
457 complex(real64), intent(inout), dimension(:) :: x
458 integer(int32), intent(out) :: marker
459
460 ! Local Variables
461 integer(int32) :: i, j
462 complex(real64) :: temp
463 real(real64) :: pivot
464
465 ! Process
466 pivot = real(x(1), real64)
467 i = 0
468 j = size(x) + 1
469 if (ascend) then
470 ! Ascending Sort
471 do
472 j = j - 1
473 do
474 if (real(x(j), real64) <= pivot) exit
475 j = j - 1
476 end do
477 i = i + 1
478 do
479 if (real(x(i), real64) >= pivot) exit
480 i = i + 1
481 end do
482 if (i < j) then
483 ! Exchage X(I) and X(J)
484 temp = x(i)
485 x(i) = x(j)
486 x(j) = temp
487 else if (i == j) then
488 marker = i + 1
489 return
490 else
491 marker = i
492 return
493 end if
494 end do
495 else
496 ! Descending Sort
497 do
498 j = j - 1
499 do
500 if (real(x(j), real64) >= pivot) exit
501 j = j - 1
502 end do
503 i = i + 1
504 do
505 if (real(x(i), real64) <= pivot) exit
506 i = i + 1
507 end do
508 if (i < j) then
509 ! Exchage X(I) and X(J)
510 temp = x(i)
511 x(i) = x(j)
512 x(j) = temp
513 else if (i == j) then
514 marker = i + 1
515 return
516 else
517 marker = i
518 return
519 end if
520 end do
521 end if
522 end subroutine
523
524! ------------------------------------------------------------------------------
542 recursive subroutine qsort_cmplx_ind(ascend, x, ind)
543 ! Arguments
544 logical, intent(in) :: ascend
545 complex(real64), intent(inout), dimension(:) :: x
546 integer(int32), intent(inout), dimension(:) :: ind
547
548 ! Local Variables
549 integer(int32) :: iq
550
551 ! Process
552 if (size(x) > 1) then
553 call cmplx_partition_ind(ascend, x, ind, iq)
554 call qsort_cmplx_ind(ascend, x(:iq-1), ind(:iq-1))
555 call qsort_cmplx_ind(ascend, x(iq:), ind(iq:))
556 end if
557 end subroutine
558
559! ------------------------------------------------------------------------------
579 subroutine cmplx_partition_ind(ascend, x, ind, marker)
580 ! Arguments
581 logical, intent(in) :: ascend
582 complex(real64), intent(inout), dimension(:) :: x
583 integer(int32), intent(inout), dimension(:) :: ind
584 integer(int32), intent(out) :: marker
585
586 ! Local Variables
587 integer(int32) :: i, j, itemp
588 complex(real64) :: temp
589 real(real64) :: pivot
590
591 ! Process
592 pivot = real(x(1), real64)
593 i = 0
594 j = size(x) + 1
595 if (ascend) then
596 ! Ascending Sort
597 do
598 j = j - 1
599 do
600 if (real(x(j), real64) <= pivot) exit
601 j = j - 1
602 end do
603 i = i + 1
604 do
605 if (real(x(i), real64) >= pivot) exit
606 i = i + 1
607 end do
608 if (i < j) then
609 ! Exchage X(I) and X(J)
610 temp = x(i)
611 x(i) = x(j)
612 x(j) = temp
613
614 itemp = ind(i)
615 ind(i) = ind(j)
616 ind(j) = itemp
617 else if (i == j) then
618 marker = i + 1
619 return
620 else
621 marker = i
622 return
623 end if
624 end do
625 else
626 ! Descending Sort
627 do
628 j = j - 1
629 do
630 if (real(x(j), real64) >= pivot) exit
631 j = j - 1
632 end do
633 i = i + 1
634 do
635 if (real(x(i), real64) <= pivot) exit
636 i = i + 1
637 end do
638 if (i < j) then
639 ! Exchage X(I) and X(J)
640 temp = x(i)
641 x(i) = x(j)
642 x(j) = temp
643
644 itemp = ind(i)
645 ind(i) = ind(j)
646 ind(j) = itemp
647 else if (i == j) then
648 marker = i + 1
649 return
650 else
651 marker = i
652 return
653 end if
654 end do
655 end if
656 end subroutine
657
658! ------------------------------------------------------------------------------
659end submodule
Provides a set of common linear algebra routines.
Definition linalg.f90:145