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