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 = 128) :: 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 write(errmsg, 100) &
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), &
71 la_array_size_error)
72 return
73 end if
74 if (n <= 1) return
75
76 ! Process
77 call qsort_dbl_ind(dir, x, ind)
78
79 ! Formatting
80100 format(a, i0, a, i0, a)
81 end subroutine
82
83! ------------------------------------------------------------------------------
84 module subroutine sort_cmplx_array(x, ascend)
85 ! Arguments
86 complex(real64), intent(inout), dimension(:) :: x
87 logical, intent(in), optional :: ascend
88
89 ! Local Variables
90 logical :: dir
91
92 ! Initialization
93 if (present(ascend)) then
94 dir = ascend
95 else
96 dir = .true.
97 end if
98
99 ! Process
100 call qsort_cmplx(dir, x)
101 end subroutine
102
103! ------------------------------------------------------------------------------
104 module subroutine sort_cmplx_array_ind(x, ind, ascend, err)
105 ! Arguments
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
110
111 ! Local Variables
112 class(errors), pointer :: errmgr
113 type(errors), target :: deferr
114 character(len = 128) :: errmsg
115 integer(int32) :: n
116 logical :: dir
117
118 ! Initialization
119 n = size(x)
120 if (present(err)) then
121 errmgr => err
122 else
123 errmgr => deferr
124 end if
125 if (present(ascend)) then
126 dir = ascend
127 else
128 dir = .true. ! Ascend == true
129 end if
130
131 ! Input Check
132 if (size(ind) /= n) then
133 write(errmsg, 100) &
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), &
137 la_array_size_error)
138 return
139 end if
140 if (n <= 1) return
141
142 ! Process
143 call qsort_cmplx_ind(dir, x, ind)
144
145 ! Formatting
146100 format(a, i0, a, i0, a)
147 end subroutine
148
149! ------------------------------------------------------------------------------
150 module subroutine sort_eigen_cmplx(vals, vecs, ascend, err)
151 ! Arguments
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
156
157 ! Local Variables
158 class(errors), pointer :: errmgr
159 type(errors), target :: deferr
160 character(len = 128) :: errmsg
161 integer(int32) :: i, n, flag
162 logical :: dir
163 integer(int32), allocatable, dimension(:) :: ind
164
165 ! Initialization
166 if (present(err)) then
167 errmgr => err
168 else
169 errmgr => deferr
170 end if
171 if (present(ascend)) then
172 dir = ascend
173 else
174 dir = .true. ! Ascend == true
175 end if
176
177 ! Ensure the eigenvector matrix is sized appropriately
178 n = size(vals)
179 if (size(vecs, 1) /= n .or. size(vecs, 2) /= n) then
180 ! ARRAY SIZE ERROR
181 write(errmsg, 100) &
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), &
186 la_array_size_error)
187 end if
188
189 ! Allocate memory for the tracking array
190 allocate(ind(n), stat = flag)
191 if (flag /= 0) then
192 call errmgr%report_error("sort_eigen_cmplx", &
193 "Insufficient memory available.", la_out_of_memory_error)
194 return
195 end if
196 do i = 1, n
197 ind(i) = i
198 end do
199
200 ! Sort
201 call qsort_cmplx_ind(dir, vals, ind)
202
203 ! Shift the eigenvectors around to keep them associated with the
204 ! appropriate eigenvalue
205 vecs = vecs(:,ind)
206
207 ! Formatting
208100 format(a, i0, a, i0, a, i0, a, i0, a)
209 end subroutine
210
211! ------------------------------------------------------------------------------
212 module subroutine sort_eigen_dbl(vals, vecs, ascend, err)
213 ! Arguments
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
218
219 ! Local Variables
220 class(errors), pointer :: errmgr
221 type(errors), target :: deferr
222 character(len = 128) :: errmsg
223 integer(int32) :: i, n, flag
224 logical :: dir
225 integer(int32), allocatable, dimension(:) :: ind
226
227 ! Initialization
228 if (present(err)) then
229 errmgr => err
230 else
231 errmgr => deferr
232 end if
233 if (present(ascend)) then
234 dir = ascend
235 else
236 dir = .true. ! Ascend == true
237 end if
238
239 ! Ensure the eigenvector matrix is sized appropriately
240 n = size(vals)
241 if (size(vecs, 1) /= n .or. size(vecs, 2) /= n) then
242 ! ARRAY SIZE ERROR
243 write(errmsg, 100) &
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), &
248 la_array_size_error)
249 end if
250
251 ! Allocate memory for the tracking array
252 allocate(ind(n), stat = flag)
253 if (flag /= 0) then
254 call errmgr%report_error("sort_eigen_dbl", &
255 "Insufficient memory available.", la_out_of_memory_error)
256 return
257 end if
258 do i = 1, n
259 ind(i) = i
260 end do
261
262 ! Sort
263 call qsort_dbl_ind(dir, vals, ind)
264
265 ! Shift the eigenvectors around to keep them associated with the
266 ! appropriate eigenvalue
267 vecs = vecs(:,ind)
268
269 ! Formatting
270100 format(a, i0, a, i0, a, i0, a, i0, a)
271 end subroutine
272
273! ******************************************************************************
274! PRIVATE HELPER ROUTINES
275! ------------------------------------------------------------------------------
289 recursive subroutine qsort_dbl_ind(ascend, x, ind)
290 ! Arguments
291 logical, intent(in) :: ascend
292 real(real64), intent(inout), dimension(:) :: x
293 integer(int32), intent(inout), dimension(:) :: ind
294
295 ! Local Variables
296 integer(int32) :: iq
297
298 ! Process
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:))
303 end if
304 end subroutine
305
306! ------------------------------------------------------------------------------
322 subroutine dbl_partition_ind(ascend, x, ind, marker)
323 ! Arguments
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
328
329 ! Local Variables
330 integer(int32) :: i, j, itemp
331 real(real64) :: temp, pivot
332
333 ! Process
334 pivot = x(1)
335 i = 0
336 j = size(x) + 1
337 if (ascend) then
338 ! Ascending Sort
339 do
340 j = j - 1
341 do
342 if (x(j) <= pivot) exit
343 j = j - 1
344 end do
345 i = i + 1
346 do
347 if (x(i) >= pivot) exit
348 i = i + 1
349 end do
350 if (i < j) then
351 ! Exchage X(I) and X(J)
352 temp = x(i)
353 x(i) = x(j)
354 x(j) = temp
355
356 itemp = ind(i)
357 ind(i) = ind(j)
358 ind(j) = itemp
359 else if (i == j) then
360 marker = i + 1
361 return
362 else
363 marker = i
364 return
365 end if
366 end do
367 else
368 ! Descending Sort
369 do
370 j = j - 1
371 do
372 if (x(j) >= pivot) exit
373 j = j - 1
374 end do
375 i = i + 1
376 do
377 if (x(i) <= pivot) exit
378 i = i + 1
379 end do
380 if (i < j) then
381 ! Exchage X(I) and X(J)
382 temp = x(i)
383 x(i) = x(j)
384 x(j) = temp
385
386 itemp = ind(i)
387 ind(i) = ind(j)
388 ind(j) = itemp
389 else if (i == j) then
390 marker = i + 1
391 return
392 else
393 marker = i
394 return
395 end if
396 end do
397 end if
398 end subroutine
399
400! ------------------------------------------------------------------------------
415 recursive subroutine qsort_cmplx(ascend, x)
416 ! Arguments
417 logical, intent(in) :: ascend
418 complex(real64), intent(inout), dimension(:) :: x
419
420 ! Local Variables
421 integer(int32) :: iq
422
423 ! Process
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:))
428 end if
429 end subroutine
430
431! ------------------------------------------------------------------------------
448 subroutine cmplx_partition(ascend, x, marker)
449 ! Arguments
450 logical, intent(in) :: ascend
451 complex(real64), intent(inout), dimension(:) :: x
452 integer(int32), intent(out) :: marker
453
454 ! Local Variables
455 integer(int32) :: i, j
456 complex(real64) :: temp
457 real(real64) :: pivot
458
459 ! Process
460 pivot = real(x(1), real64)
461 i = 0
462 j = size(x) + 1
463 if (ascend) then
464 ! Ascending Sort
465 do
466 j = j - 1
467 do
468 if (real(x(j), real64) <= pivot) exit
469 j = j - 1
470 end do
471 i = i + 1
472 do
473 if (real(x(i), real64) >= pivot) exit
474 i = i + 1
475 end do
476 if (i < j) then
477 ! Exchage X(I) and X(J)
478 temp = x(i)
479 x(i) = x(j)
480 x(j) = temp
481 else if (i == j) then
482 marker = i + 1
483 return
484 else
485 marker = i
486 return
487 end if
488 end do
489 else
490 ! Descending Sort
491 do
492 j = j - 1
493 do
494 if (real(x(j), real64) >= pivot) exit
495 j = j - 1
496 end do
497 i = i + 1
498 do
499 if (real(x(i), real64) <= pivot) exit
500 i = i + 1
501 end do
502 if (i < j) then
503 ! Exchage X(I) and X(J)
504 temp = x(i)
505 x(i) = x(j)
506 x(j) = temp
507 else if (i == j) then
508 marker = i + 1
509 return
510 else
511 marker = i
512 return
513 end if
514 end do
515 end if
516 end subroutine
517
518! ------------------------------------------------------------------------------
536 recursive subroutine qsort_cmplx_ind(ascend, x, ind)
537 ! Arguments
538 logical, intent(in) :: ascend
539 complex(real64), intent(inout), dimension(:) :: x
540 integer(int32), intent(inout), dimension(:) :: ind
541
542 ! Local Variables
543 integer(int32) :: iq
544
545 ! Process
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:))
550 end if
551 end subroutine
552
553! ------------------------------------------------------------------------------
573 subroutine cmplx_partition_ind(ascend, x, ind, marker)
574 ! Arguments
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
579
580 ! Local Variables
581 integer(int32) :: i, j, itemp
582 complex(real64) :: temp
583 real(real64) :: pivot
584
585 ! Process
586 pivot = real(x(1), real64)
587 i = 0
588 j = size(x) + 1
589 if (ascend) then
590 ! Ascending Sort
591 do
592 j = j - 1
593 do
594 if (real(x(j), real64) <= pivot) exit
595 j = j - 1
596 end do
597 i = i + 1
598 do
599 if (real(x(i), real64) >= pivot) exit
600 i = i + 1
601 end do
602 if (i < j) then
603 ! Exchage X(I) and X(J)
604 temp = x(i)
605 x(i) = x(j)
606 x(j) = temp
607
608 itemp = ind(i)
609 ind(i) = ind(j)
610 ind(j) = itemp
611 else if (i == j) then
612 marker = i + 1
613 return
614 else
615 marker = i
616 return
617 end if
618 end do
619 else
620 ! Descending Sort
621 do
622 j = j - 1
623 do
624 if (real(x(j), real64) >= pivot) exit
625 j = j - 1
626 end do
627 i = i + 1
628 do
629 if (real(x(i), real64) <= pivot) exit
630 i = i + 1
631 end do
632 if (i < j) then
633 ! Exchage X(I) and X(J)
634 temp = x(i)
635 x(i) = x(j)
636 x(j) = temp
637
638 itemp = ind(i)
639 ind(i) = ind(j)
640 ind(j) = itemp
641 else if (i == j) then
642 marker = i + 1
643 return
644 else
645 marker = i
646 return
647 end if
648 end do
649 end if
650 end subroutine
651
652! ------------------------------------------------------------------------------
653end submodule
Provides a set of common linear algebra routines.
Definition: linalg.f90:145