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_solve.f90
1! linalg_solve.f90
2
7submodule(linalg) linalg_solve
8contains
9! ******************************************************************************
10! TRIANGULAR MATRIX SOLUTION ROUTINES
11! ------------------------------------------------------------------------------
12 module subroutine solve_tri_mtx(lside, upper, trans, nounit, alpha, a, b, err)
13 ! Arguments
14 logical, intent(in) :: lside, upper, trans, nounit
15 real(real64), intent(in) :: alpha
16 real(real64), intent(in), dimension(:,:) :: a
17 real(real64), intent(inout), dimension(:,:) :: b
18 class(errors), intent(inout), optional, target :: err
19
20 ! Parameters
21 character :: side, uplo, transa, diag
22
23 ! Local Variables
24 integer(int32) :: m, n, nrowa
25 class(errors), pointer :: errmgr
26 type(errors), target :: deferr
27
28 ! Initialization
29 m = size(b, 1)
30 n = size(b, 2)
31 if (lside) then
32 nrowa = m
33 side = 'L'
34 else
35 nrowa = n
36 side = 'R'
37 end if
38 if (upper) then
39 uplo = 'U'
40 else
41 uplo = 'L'
42 end if
43 if (trans) then
44 transa = 'T'
45 else
46 transa = 'N'
47 end if
48 if (nounit) then
49 diag = 'N'
50 else
51 diag = 'U'
52 end if
53 if (present(err)) then
54 errmgr => err
55 else
56 errmgr => deferr
57 end if
58
59 ! Input Check - matrix A must be square
60 if (size(a, 1) /= nrowa .or. size(a, 2) /= nrowa) then
61 ! ERROR: A must be square
62 call errmgr%report_error("solve_tri_mtx", &
63 "The input matrix must be square.", la_array_size_error)
64 return
65 end if
66
67 ! Call DTRSM
68 call dtrsm(side, uplo, transa, diag, m, n, alpha, a, nrowa, b, m)
69 end subroutine
70
71! ------------------------------------------------------------------------------
72 module subroutine solve_tri_mtx_cmplx(lside, upper, trans, nounit, alpha, a, b, err)
73 ! Arguments
74 logical, intent(in) :: lside, upper, trans, nounit
75 complex(real64), intent(in) :: alpha
76 complex(real64), intent(in), dimension(:,:) :: a
77 complex(real64), intent(inout), dimension(:,:) :: b
78 class(errors), intent(inout), optional, target :: err
79
80 ! Parameters
81 character :: side, uplo, transa, diag
82
83 ! Local Variables
84 integer(int32) :: m, n, nrowa
85 class(errors), pointer :: errmgr
86 type(errors), target :: deferr
87
88 ! Initialization
89 m = size(b, 1)
90 n = size(b, 2)
91 if (lside) then
92 nrowa = m
93 side = 'L'
94 else
95 nrowa = n
96 side = 'R'
97 end if
98 if (upper) then
99 uplo = 'U'
100 else
101 uplo = 'L'
102 end if
103 if (trans) then
104 transa = 'C'
105 else
106 transa = 'N'
107 end if
108 if (nounit) then
109 diag = 'N'
110 else
111 diag = 'U'
112 end if
113 if (present(err)) then
114 errmgr => err
115 else
116 errmgr => deferr
117 end if
118
119 ! Input Check - matrix A must be square
120 if (size(a, 1) /= nrowa .or. size(a, 2) /= nrowa) then
121 ! ERROR: A must be square
122 call errmgr%report_error("solve_tri_mtx_cmplx", &
123 "The input matrix must be square.", la_array_size_error)
124 return
125 end if
126
127 ! Call ZTRSM
128 call ztrsm(side, uplo, transa, diag, m, n, alpha, a, nrowa, b, m)
129 end subroutine
130
131! ------------------------------------------------------------------------------
132 module subroutine solve_tri_vec(upper, trans, nounit, a, x, err)
133 ! Arguments
134 logical, intent(in) :: upper, trans, nounit
135 real(real64), intent(in), dimension(:,:) :: a
136 real(real64), intent(inout), dimension(:) :: x
137 class(errors), intent(inout), optional, target :: err
138
139 ! Parameters
140 real(real64), parameter :: zero = 0.0d0
141
142 ! Local Variables
143 character :: uplo, t, diag
144 integer(int32) :: n
145 class(errors), pointer :: errmgr
146 type(errors), target :: deferr
147
148 ! Initialization
149 n = size(a, 1)
150 if (upper) then
151 uplo = 'U'
152 else
153 uplo = 'L'
154 end if
155 if (trans) then
156 t = 'T'
157 else
158 t = 'N'
159 end if
160 if (nounit) then
161 diag = 'N'
162 else
163 diag = 'U'
164 end if
165 if (present(err)) then
166 errmgr => err
167 else
168 errmgr => deferr
169 end if
170
171 ! Input Check
172 if (size(a, 2) /= n) then
173 ! ERROR: A must be square
174 call errmgr%report_error("solve_tri_vec", &
175 "The input matrix must be square.", la_array_size_error)
176 return
177 else if (size(x) /= n) then
178 ! ERROR: Inner matrix dimensions must agree
179 call errmgr%report_error("solve_tri_vec", &
180 "The inner matrix dimensions must be equal.", &
181 la_array_size_error)
182 return
183 end if
184
185 ! Call DTRSV
186 call dtrsv(uplo, t, diag, n, a, n, x, 1)
187 end subroutine
188
189! ------------------------------------------------------------------------------
190 module subroutine solve_tri_vec_cmplx(upper, trans, nounit, a, x, err)
191 ! Arguments
192 logical, intent(in) :: upper, trans, nounit
193 complex(real64), intent(in), dimension(:,:) :: a
194 complex(real64), intent(inout), dimension(:) :: x
195 class(errors), intent(inout), optional, target :: err
196
197 ! Parameters
198 real(real64), parameter :: zero = 0.0d0
199
200 ! Local Variables
201 character :: uplo, t, diag
202 integer(int32) :: n
203 class(errors), pointer :: errmgr
204 type(errors), target :: deferr
205
206 ! Initialization
207 n = size(a, 1)
208 if (upper) then
209 uplo = 'U'
210 else
211 uplo = 'L'
212 end if
213 if (trans) then
214 t = 'C'
215 else
216 t = 'N'
217 end if
218 if (nounit) then
219 diag = 'N'
220 else
221 diag = 'U'
222 end if
223 if (present(err)) then
224 errmgr => err
225 else
226 errmgr => deferr
227 end if
228
229 ! Input Check
230 if (size(a, 2) /= n) then
231 ! ERROR: A must be square
232 call errmgr%report_error("solve_tri_vec_cmplx", &
233 "The input matrix must be square.", la_array_size_error)
234 return
235 else if (size(x) /= n) then
236 ! ERROR: Inner matrix dimensions must agree
237 call errmgr%report_error("solve_tri_vec_cmplx", &
238 "The inner matrix dimensions must be equal.", &
239 la_array_size_error)
240 return
241 end if
242
243 ! Call ZTRSV
244 call ztrsv(uplo, t, diag, n, a, n, x, 1)
245 end subroutine
246
247! ******************************************************************************
248! LU SOLUTION
249! ------------------------------------------------------------------------------
250 module subroutine solve_lu_mtx(a, ipvt, b, err)
251 ! Arguments
252 real(real64), intent(in), dimension(:,:) :: a
253 integer(int32), intent(in), dimension(:) :: ipvt
254 real(real64), intent(inout), dimension(:,:) :: b
255 class(errors), intent(inout), optional, target :: err
256
257 ! Local Variables
258 integer(int32) :: n, nrhs, flag
259 class(errors), pointer :: errmgr
260 type(errors), target :: deferr
261 character(len = :), allocatable :: errmsg
262
263 ! Initialization
264 n = size(a, 1)
265 nrhs = size(b, 2)
266 if (present(err)) then
267 errmgr => err
268 else
269 errmgr => deferr
270 end if
271
272 ! Input Check
273 flag = 0
274 if (size(a, 2) /= n) then
275 flag = 1
276 else if (size(ipvt) /= n) then
277 flag = 2
278 else if (size(b, 1) /= n) then
279 flag = 3
280 end if
281 if (flag /= 0) then
282 ! One of the input arrays is not sized correctly
283 allocate(character(len = 256) :: errmsg)
284 write(errmsg, 100) "Input number ", flag, &
285 " is not sized correctly."
286 call errmgr%report_error("solve_lu_mtx", trim(errmsg), &
287 la_array_size_error)
288 return
289 end if
290
291 ! Call DGETRS
292 call dgetrs("N", n, nrhs, a, n, ipvt, b, n, flag)
293
294 ! Formatting
295100 format(a, i0, a)
296 end subroutine
297
298! ------------------------------------------------------------------------------
299 module subroutine solve_lu_mtx_cmplx(a, ipvt, b, err)
300 ! Arguments
301 complex(real64), intent(in), dimension(:,:) :: a
302 integer(int32), intent(in), dimension(:) :: ipvt
303 complex(real64), intent(inout), dimension(:,:) :: b
304 class(errors), intent(inout), optional, target :: err
305
306 ! Local Variables
307 integer(int32) :: n, nrhs, flag
308 class(errors), pointer :: errmgr
309 type(errors), target :: deferr
310 character(len = :), allocatable :: errmsg
311
312 ! Initialization
313 n = size(a, 1)
314 nrhs = size(b, 2)
315 if (present(err)) then
316 errmgr => err
317 else
318 errmgr => deferr
319 end if
320
321 ! Input Check
322 flag = 0
323 if (size(a, 2) /= n) then
324 flag = 1
325 else if (size(ipvt) /= n) then
326 flag = 2
327 else if (size(b, 1) /= n) then
328 flag = 3
329 end if
330 if (flag /= 0) then
331 ! One of the input arrays is not sized correctly
332 allocate(character(len = 256) :: errmsg)
333 write(errmsg, 100) "Input number ", flag, &
334 " is not sized correctly."
335 call errmgr%report_error("solve_lu_mtx_cmplx", trim(errmsg), &
336 la_array_size_error)
337 return
338 end if
339
340 ! Call ZGETRS
341 call zgetrs("N", n, nrhs, a, n, ipvt, b, n, flag)
342
343 ! Formatting
344100 format(a, i0, a)
345 end subroutine
346
347! ------------------------------------------------------------------------------
348 module subroutine solve_lu_vec(a, ipvt, b, err)
349 ! Arguments
350 real(real64), intent(in), dimension(:,:) :: a
351 integer(int32), intent(in), dimension(:) :: ipvt
352 real(real64), intent(inout), dimension(:) :: b
353 class(errors), intent(inout), optional, target :: err
354
355 ! Local Variables
356 integer(int32) :: n, flag
357 class(errors), pointer :: errmgr
358 type(errors), target :: deferr
359 character(len = :), allocatable :: errmsg
360
361 ! Initialization
362 n = size(a, 1)
363 if (present(err)) then
364 errmgr => err
365 else
366 errmgr => deferr
367 end if
368
369 ! Input Check
370 flag = 0
371 if (size(a, 2) /= n) then
372 flag = 1
373 else if (size(ipvt) /= n) then
374 flag = 2
375 else if (size(b) /= n) then
376 flag = 3
377 end if
378 if (flag /= 0) then
379 ! One of the input arrays is not sized correctly
380 allocate(character(len = 256) :: errmsg)
381 write(errmsg, 100) "Input number ", flag, &
382 " is not sized correctly."
383 call errmgr%report_error("solve_lu_vec", trim(errmsg), &
384 la_array_size_error)
385 return
386 end if
387
388 ! Call DGETRS
389 call dgetrs("N", n, 1, a, n, ipvt, b, n, flag)
390
391 ! Formatting
392100 format(a, i0, a)
393 end subroutine
394
395! ------------------------------------------------------------------------------
396 module subroutine solve_lu_vec_cmplx(a, ipvt, b, err)
397 ! Arguments
398 complex(real64), intent(in), dimension(:,:) :: a
399 integer(int32), intent(in), dimension(:) :: ipvt
400 complex(real64), intent(inout), dimension(:) :: b
401 class(errors), intent(inout), optional, target :: err
402
403 ! Local Variables
404 integer(int32) :: n, flag
405 class(errors), pointer :: errmgr
406 type(errors), target :: deferr
407 character(len = :), allocatable :: errmsg
408
409 ! Initialization
410 n = size(a, 1)
411 if (present(err)) then
412 errmgr => err
413 else
414 errmgr => deferr
415 end if
416
417 ! Input Check
418 flag = 0
419 if (size(a, 2) /= n) then
420 flag = 1
421 else if (size(ipvt) /= n) then
422 flag = 2
423 else if (size(b) /= n) then
424 flag = 3
425 end if
426 if (flag /= 0) then
427 ! One of the input arrays is not sized correctly
428 allocate(character(len = 256) :: errmsg)
429 write(errmsg, 100) "Input number ", flag, &
430 " is not sized correctly."
431 call errmgr%report_error("solve_lu_vec_cmplx", trim(errmsg), &
432 la_array_size_error)
433 return
434 end if
435
436 ! Call ZGETRS
437 call zgetrs("N", n, 1, a, n, ipvt, b, n, flag)
438
439 ! Formatting
440100 format(a, i0, a)
441 end subroutine
442
443! ******************************************************************************
444! QR SOLUTION
445! ------------------------------------------------------------------------------
446 module subroutine solve_qr_no_pivot_mtx(a, tau, b, work, olwork, err)
447 ! Arguments
448 real(real64), intent(inout), dimension(:,:) :: a, b
449 real(real64), intent(in), dimension(:) :: tau
450 real(real64), intent(out), target, optional, dimension(:) :: work
451 integer(int32), intent(out), optional :: olwork
452 class(errors), intent(inout), optional, target :: err
453
454 ! Parameters
455 real(real64), parameter :: one = 1.0d0
456
457 ! Local Variables
458 integer(int32) :: m, n, nrhs, k, lwork, flag, istat
459 real(real64), pointer, dimension(:) :: wptr
460 real(real64), allocatable, target, dimension(:) :: wrk
461 class(errors), pointer :: errmgr
462 type(errors), target :: deferr
463 character(len = :), allocatable :: errmsg
464
465 ! Initialization
466 m = size(a, 1)
467 n = size(a, 2)
468 nrhs = size(b, 2)
469 k = min(m, n)
470 if (present(err)) then
471 errmgr => err
472 else
473 errmgr => deferr
474 end if
475
476 ! Input Check
477 flag = 0
478 if (m < n) then
479 flag = 1
480 else if (size(tau) /= k) then
481 flag = 2
482 else if (size(b, 1) /= m) then
483 flag = 3
484 end if
485 if (flag /= 0) then
486 ! ERROR: One of the input arrays is not sized correctly
487 allocate(character(len = 256) :: errmsg)
488 write(errmsg, 100) "Input number ", flag, &
489 " is not sized correctly."
490 call errmgr%report_error("solve_qr_no_pivot_mtx", trim(errmsg), &
491 la_array_size_error)
492 return
493 end if
494
495 ! Workspace Query
496 call mult_qr(.true., .true., a, tau, b, olwork = lwork)
497 if (present(olwork)) then
498 olwork = lwork
499 return
500 end if
501
502 ! Local Memory Allocation
503 if (present(work)) then
504 if (size(work) < lwork) then
505 ! ERROR: WORK not sized correctly
506 call errmgr%report_error("solve_qr_no_pivot_mtx", &
507 "Incorrectly sized input array WORK, argument 4.", &
508 la_array_size_error)
509 return
510 end if
511 wptr => work(1:lwork)
512 else
513 allocate(wrk(lwork), stat = istat)
514 if (istat /= 0) then
515 ! ERROR: Out of memory
516 call errmgr%report_error("solve_qr_no_pivot_mtx", &
517 "Insufficient memory available.", &
518 la_out_of_memory_error)
519 return
520 end if
521 wptr => wrk
522 end if
523
524 ! Compute Q**T * B, and store in B
525 call mult_qr(.true., .true., a, tau, b, wptr)
526
527 ! Solve the triangular system: A(1:N,1:N)*X = B(1:N,:)
528 call solve_triangular_system(.true., .true., .false., .true., one, &
529 a(1:n,1:n), b(1:n,:))
530
531 ! Formatting
532100 format(a, i0, a)
533 end subroutine
534
535! ------------------------------------------------------------------------------
536 module subroutine solve_qr_no_pivot_mtx_cmplx(a, tau, b, work, olwork, err)
537 ! Arguments
538 complex(real64), intent(inout), dimension(:,:) :: a, b
539 complex(real64), intent(in), dimension(:) :: tau
540 complex(real64), intent(out), target, optional, dimension(:) :: work
541 integer(int32), intent(out), optional :: olwork
542 class(errors), intent(inout), optional, target :: err
543
544 ! Parameters
545 complex(real64), parameter :: one = (1.0d0, 0.0d0)
546
547 ! Local Variables
548 integer(int32) :: m, n, nrhs, k, lwork, flag, istat
549 complex(real64), pointer, dimension(:) :: wptr
550 complex(real64), allocatable, target, dimension(:) :: wrk
551 class(errors), pointer :: errmgr
552 type(errors), target :: deferr
553 character(len = :), allocatable :: errmsg
554
555 ! Initialization
556 m = size(a, 1)
557 n = size(a, 2)
558 nrhs = size(b, 2)
559 k = min(m, n)
560 if (present(err)) then
561 errmgr => err
562 else
563 errmgr => deferr
564 end if
565
566 ! Input Check
567 flag = 0
568 if (m < n) then
569 flag = 1
570 else if (size(tau) /= k) then
571 flag = 2
572 else if (size(b, 1) /= m) then
573 flag = 3
574 end if
575 if (flag /= 0) then
576 ! ERROR: One of the input arrays is not sized correctly
577 allocate(character(len = 256) :: errmsg)
578 write(errmsg, 100) "Input number ", flag, &
579 " is not sized correctly."
580 call errmgr%report_error("solve_qr_no_pivot_mtx_cmplx", &
581 trim(errmsg), la_array_size_error)
582 return
583 end if
584
585 ! Workspace Query
586 call mult_qr(.true., .true., a, tau, b, olwork = lwork)
587 if (present(olwork)) then
588 olwork = lwork
589 return
590 end if
591
592 ! Local Memory Allocation
593 if (present(work)) then
594 if (size(work) < lwork) then
595 ! ERROR: WORK not sized correctly
596 call errmgr%report_error("solve_qr_no_pivot_mtx_cmplx", &
597 "Incorrectly sized input array WORK, argument 4.", &
598 la_array_size_error)
599 return
600 end if
601 wptr => work(1:lwork)
602 else
603 allocate(wrk(lwork), stat = istat)
604 if (istat /= 0) then
605 ! ERROR: Out of memory
606 call errmgr%report_error("solve_qr_no_pivot_mtx_cmplx", &
607 "Insufficient memory available.", &
608 la_out_of_memory_error)
609 return
610 end if
611 wptr => wrk
612 end if
613
614 ! Compute Q**T * B, and store in B
615 call mult_qr(.true., .true., a, tau, b, wptr)
616
617 ! Solve the triangular system: A(1:N,1:N)*X = B(1:N,:)
618 call solve_triangular_system(.true., .true., .false., .true., one, &
619 a(1:n,1:n), b(1:n,:))
620
621 ! Formatting
622100 format(a, i0, a)
623 end subroutine
624
625! ------------------------------------------------------------------------------
626 module subroutine solve_qr_no_pivot_vec(a, tau, b, work, olwork, err)
627 ! Arguments
628 real(real64), intent(inout), dimension(:,:) :: a
629 real(real64), intent(in), dimension(:) :: tau
630 real(real64), intent(inout), dimension(:) :: b
631 real(real64), intent(out), target, optional, dimension(:) :: work
632 integer(int32), intent(out), optional :: olwork
633 class(errors), intent(inout), optional, target :: err
634
635 ! Local Variables
636 integer(int32) :: m, n, k, flag, lwork, istat
637 real(real64), pointer, dimension(:) :: wptr
638 real(real64), allocatable, target, dimension(:) :: wrk
639 class(errors), pointer :: errmgr
640 type(errors), target :: deferr
641 character(len = :), allocatable :: errmsg
642
643 ! Initialization
644 m = size(a, 1)
645 n = size(a, 2)
646 k = min(m, n)
647 if (present(err)) then
648 errmgr => err
649 else
650 errmgr => deferr
651 end if
652
653 ! Input Check
654 flag = 0
655 if (m < n) then
656 flag = 1
657 else if (size(tau) /= k) then
658 flag = 2
659 else if (size(b) /= m) then
660 flag = 3
661 end if
662 if (flag /= 0) then
663 ! ERROR: One of the input arrays is not sized correctly
664 allocate(character(len = 256) :: errmsg)
665 write(errmsg, 100) "Input number ", flag, &
666 " is not sized correctly."
667 call errmgr%report_error("solve_qr_no_pivot_vec", trim(errmsg), &
668 la_array_size_error)
669 return
670 end if
671
672 ! Workspace Query
673 call mult_qr(.true., a, tau, b, olwork = lwork)
674 if (present(olwork)) then
675 olwork = lwork
676 return
677 end if
678
679 ! Local Memory Allocation
680 if (present(work)) then
681 if (size(work) < lwork) then
682 ! ERROR: WORK not sized correctly
683 call errmgr%report_error("solve_qr_no_pivot_vec", &
684 "Incorrectly sized input array WORK, argument 4.", &
685 la_array_size_error)
686 return
687 end if
688 wptr => work(1:lwork)
689 else
690 allocate(wrk(lwork), stat = istat)
691 if (istat /= 0) then
692 ! ERROR: Out of memory
693 call errmgr%report_error("solve_qr_no_pivot_vec", &
694 "Insufficient memory available.", &
695 la_out_of_memory_error)
696 return
697 end if
698 wptr => wrk
699 end if
700
701 ! Compute Q**T * B, and store in B
702 call mult_qr(.true., a, tau, b, work = wptr)
703
704 ! Solve the triangular system: A(1:N,1:N)*X = B(1:N)
705 call solve_triangular_system(.true., .false., .true., a(1:n,1:n), b)
706
707 ! Formatting
708100 format(a, i0, a)
709 end subroutine
710
711! ------------------------------------------------------------------------------
712 module subroutine solve_qr_no_pivot_vec_cmplx(a, tau, b, work, olwork, err)
713 ! Arguments
714 complex(real64), intent(inout), dimension(:,:) :: a
715 complex(real64), intent(in), dimension(:) :: tau
716 complex(real64), intent(inout), dimension(:) :: b
717 complex(real64), intent(out), target, optional, dimension(:) :: work
718 integer(int32), intent(out), optional :: olwork
719 class(errors), intent(inout), optional, target :: err
720
721 ! Local Variables
722 integer(int32) :: m, n, k, flag, lwork, istat
723 complex(real64), pointer, dimension(:) :: wptr
724 complex(real64), allocatable, target, dimension(:) :: wrk
725 class(errors), pointer :: errmgr
726 type(errors), target :: deferr
727 character(len = :), allocatable :: errmsg
728
729 ! Initialization
730 m = size(a, 1)
731 n = size(a, 2)
732 k = min(m, n)
733 if (present(err)) then
734 errmgr => err
735 else
736 errmgr => deferr
737 end if
738
739 ! Input Check
740 flag = 0
741 if (m < n) then
742 flag = 1
743 else if (size(tau) /= k) then
744 flag = 2
745 else if (size(b) /= m) then
746 flag = 3
747 end if
748 if (flag /= 0) then
749 ! ERROR: One of the input arrays is not sized correctly
750 allocate(character(len = 256) :: errmsg)
751 write(errmsg, 100) "Input number ", flag, &
752 " is not sized correctly."
753 call errmgr%report_error("solve_qr_no_pivot_vec_cmplx", &
754 trim(errmsg), la_array_size_error)
755 return
756 end if
757
758 ! Workspace Query
759 call mult_qr(.true., a, tau, b, olwork = lwork)
760 if (present(olwork)) then
761 olwork = lwork
762 return
763 end if
764
765 ! Local Memory Allocation
766 if (present(work)) then
767 if (size(work) < lwork) then
768 ! ERROR: WORK not sized correctly
769 call errmgr%report_error("solve_qr_no_pivot_vec_cmplx", &
770 "Incorrectly sized input array WORK, argument 4.", &
771 la_array_size_error)
772 return
773 end if
774 wptr => work(1:lwork)
775 else
776 allocate(wrk(lwork), stat = istat)
777 if (istat /= 0) then
778 ! ERROR: Out of memory
779 call errmgr%report_error("solve_qr_no_pivot_vec_cmplx", &
780 "Insufficient memory available.", &
781 la_out_of_memory_error)
782 return
783 end if
784 wptr => wrk
785 end if
786
787 ! Compute Q**T * B, and store in B
788 call mult_qr(.true., a, tau, b, work = wptr)
789
790 ! Solve the triangular system: A(1:N,1:N)*X = B(1:N)
791 call solve_triangular_system(.true., .false., .true., a(1:n,1:n), b)
792
793 ! Formatting
794100 format(a, i0, a)
795 end subroutine
796
797! ------------------------------------------------------------------------------
798 module subroutine solve_qr_pivot_mtx(a, tau, jpvt, b, work, olwork, err)
799 ! Arguments
800 real(real64), intent(inout), dimension(:,:) :: a
801 real(real64), intent(in), dimension(:) :: tau
802 integer(int32), intent(in), dimension(:) :: jpvt
803 real(real64), intent(inout), dimension(:,:) :: b
804 real(real64), intent(out), target, optional, dimension(:) :: work
805 integer(int32), intent(out), optional :: olwork
806 class(errors), intent(inout), optional, target :: err
807
808 ! Parameters
809 integer(int32), parameter :: imin = 2
810 integer(int32), parameter :: imax = 1
811 real(real64), parameter :: zero = 0.0d0
812 real(real64), parameter :: one = 1.0d0
813
814 ! Local Variables
815 integer(int32) :: i, j, m, n, mn, nrhs, lwork, ismin, ismax, &
816 rnk, maxmn, flag, istat, lwork1, lwork2, lwork3
817 real(real64) :: rcond, smax, smin, smaxpr, sminpr, s1, c1, s2, c2
818 real(real64), pointer, dimension(:) :: wptr, w, tau2
819 real(real64), allocatable, target, dimension(:) :: wrk
820 class(errors), pointer :: errmgr
821 type(errors), target :: deferr
822 character(len = :), allocatable :: errmsg
823
824 ! Initialization
825 m = size(a, 1)
826 n = size(a, 2)
827 mn = min(m, n)
828 maxmn = max(m, n)
829 nrhs = size(b, 2)
830 ismin = mn + 1
831 ismax = 2 * mn + 1
832 rcond = epsilon(rcond)
833 if (present(err)) then
834 errmgr => err
835 else
836 errmgr => deferr
837 end if
838
839 ! Input Check
840 flag = 0
841 if (size(tau) /= mn) then
842 flag = 2
843 else if (size(jpvt) /= n) then
844 flag = 3
845 else if (size(b, 1) /= maxmn) then
846 flag = 4
847 end if
848 if (flag /= 0) then
849 ! ERROR: One of the input arrays is not sized correctly
850 allocate(character(len = 256) :: errmsg)
851 write(errmsg, 100) "Input number ", flag, &
852 " is not sized correctly."
853 call errmgr%report_error("solve_qr_pivot_mtx", trim(errmsg), &
854 la_array_size_error)
855 return
856 end if
857
858 ! Workspace Query
859 call rz_factor(a(1:mn,:), a(1:mn,1), olwork = lwork1)
860 call mult_qr(.true., .true., a, tau, b(1:m,:), olwork = lwork2)
861 call mult_rz(.true., .true., n, a(1:mn,:), a(1:mn,1), b(1:n,:), &
862 olwork = lwork3)
863 lwork = max(lwork1, lwork2, lwork3, 2 * mn + 1) + mn
864 if (present(olwork)) then
865 olwork = lwork
866 return
867 end if
868
869 ! Local Memory Allocation
870 if (present(work)) then
871 if (size(work) < lwork) then
872 ! ERROR: WORK not sized correctly
873 call errmgr%report_error("solve_qr_no_pivot_mtx", &
874 "Incorrectly sized input array WORK, argument 5.", &
875 la_array_size_error)
876 return
877 end if
878 wptr => work(1:lwork)
879 else
880 allocate(wrk(lwork), stat = istat)
881 if (istat /= 0) then
882 ! ERROR: Out of memory
883 call errmgr%report_error("solve_qr_pivot_mtx", &
884 "Insufficient memory available.", &
885 la_out_of_memory_error)
886 return
887 end if
888 wptr => wrk
889 end if
890
891 ! Determine the rank of R11 using an incremental condition estimation
892 wptr(ismin) = one
893 wptr(ismax) = one
894 smax = abs(a(1,1))
895 smin = smax
896 if (abs(a(1,1)) == zero) then
897 rnk = 0
898 b(1:maxmn,:) = zero
899 return
900 else
901 rnk = 1
902 end if
903 do
904 if (rnk < mn) then
905 i = rnk + 1
906 call dlaic1(imin, rnk, wptr(ismin:ismin+rnk-1), smin, &
907 a(1:rnk-1,i), a(i,i), sminpr, s1, c1)
908 call dlaic1(imax, rnk, wptr(ismax:ismax+rnk-1), smax, &
909 a(1:rnk-1,i), a(i,i), smaxpr, s2, c2)
910 if (smaxpr * rcond <= sminpr) then
911 do i = 1, rnk
912 wptr(ismin+i-1) = s1 * wptr(ismin+i-1)
913 wptr(ismax+i-1) = s2 * wptr(ismax+i-1)
914 end do
915 wptr(ismin+rnk) = c1
916 wptr(ismax+rnk) = c2
917 smin = sminpr
918 smax = smaxpr
919 rnk = rnk + 1
920 cycle
921 end if
922 end if
923 exit
924 end do
925
926 ! Partition R = [R11 R12]
927 ! [ 0 R22]
928 tau2 => wptr(1:rnk)
929 w => wptr(rnk+1:lwork)
930 if (rnk < n) call rz_factor(a(1:rnk,:), tau2, w)
931
932 ! Compute B(1:m,1:NRHS) = Q**T * B(1:M,1:NRHS)
933 call mult_qr(.true., .true., a, tau, b(1:m,:), w)
934
935 ! Solve the triangular system T11 * B(1:rnk,1:nrhs) = B(1:rnk,1:nrhs)
936 call solve_triangular_system(.true., .true., .false., .true., one, &
937 a(1:rnk,1:rnk), b(1:rnk,:))
938 if (n > rnk) b(rnk+1:n,:) = zero
939
940 ! Compute B(1:n,1:nrhs) = Y**T * B(1:n,1:nrhs)
941 if (rnk < n) then
942 call mult_rz(.true., .true., n - rnk, a(1:rnk,:), tau2, b(1:n,:), w)
943 end if
944
945 ! Apply the pivoting: B(1:N,1:NRHS) = P * B(1:N,1:NRHS)
946 do j = 1, nrhs
947 do i = 1, n
948 wptr(jpvt(i)) = b(i,j)
949 end do
950 b(:,j) = wptr(1:n)
951 end do
952
953 ! Formatting
954100 format(a, i0, a)
955 end subroutine
956
957! ------------------------------------------------------------------------------
958 module subroutine solve_qr_pivot_mtx_cmplx(a, tau, jpvt, b, work, olwork, err)
959 ! Arguments
960 complex(real64), intent(inout), dimension(:,:) :: a
961 complex(real64), intent(in), dimension(:) :: tau
962 integer(int32), intent(in), dimension(:) :: jpvt
963 complex(real64), intent(inout), dimension(:,:) :: b
964 complex(real64), intent(out), target, optional, dimension(:) :: work
965 integer(int32), intent(out), optional :: olwork
966 class(errors), intent(inout), optional, target :: err
967
968 ! Parameters
969 integer(int32), parameter :: imin = 2
970 integer(int32), parameter :: imax = 1
971 complex(real64), parameter :: zero = (0.0d0, 0.0d0)
972 complex(real64), parameter :: one = (1.0d0, 0.0d0)
973
974 ! Local Variables
975 integer(int32) :: i, j, m, n, mn, nrhs, lwork, ismin, ismax, &
976 rnk, maxmn, flag, istat, lwork1, lwork2, lwork3
977 real(real64) :: rcond, smax, smin, smaxpr, sminpr
978 complex(real64) :: s1, c1, s2, c2
979 complex(real64), pointer, dimension(:) :: wptr, w, tau2
980 complex(real64), allocatable, target, dimension(:) :: wrk
981 class(errors), pointer :: errmgr
982 type(errors), target :: deferr
983 character(len = :), allocatable :: errmsg
984
985 ! Initialization
986 m = size(a, 1)
987 n = size(a, 2)
988 mn = min(m, n)
989 maxmn = max(m, n)
990 nrhs = size(b, 2)
991 ismin = mn + 1
992 ismax = 2 * mn + 1
993 rcond = epsilon(rcond)
994 if (present(err)) then
995 errmgr => err
996 else
997 errmgr => deferr
998 end if
999
1000 ! Input Check
1001 flag = 0
1002 if (size(tau) /= mn) then
1003 flag = 2
1004 else if (size(jpvt) /= n) then
1005 flag = 3
1006 else if (size(b, 1) /= maxmn) then
1007 flag = 4
1008 end if
1009 if (flag /= 0) then
1010 ! ERROR: One of the input arrays is not sized correctly
1011 allocate(character(len = 256) :: errmsg)
1012 write(errmsg, 100) "Input number ", flag, &
1013 " is not sized correctly."
1014 call errmgr%report_error("solve_qr_pivot_mtx_cmplx", &
1015 trim(errmsg), la_array_size_error)
1016 return
1017 end if
1018
1019 ! Workspace Query
1020 call rz_factor(a(1:mn,:), a(1:mn,1), olwork = lwork1)
1021 call mult_qr(.true., .true., a, tau, b(1:m,:), olwork = lwork2)
1022 call mult_rz(.true., .true., n, a(1:mn,:), a(1:mn,1), b(1:n,:), &
1023 olwork = lwork3)
1024 lwork = max(lwork1, lwork2, lwork3, 2 * mn + 1) + mn
1025 if (present(olwork)) then
1026 olwork = lwork
1027 return
1028 end if
1029
1030 ! Local Memory Allocation
1031 if (present(work)) then
1032 if (size(work) < lwork) then
1033 ! ERROR: WORK not sized correctly
1034 call errmgr%report_error("solve_qr_no_pivot_mtx_cmplx", &
1035 "Incorrectly sized input array WORK, argument 5.", &
1036 la_array_size_error)
1037 return
1038 end if
1039 wptr => work(1:lwork)
1040 else
1041 allocate(wrk(lwork), stat = istat)
1042 if (istat /= 0) then
1043 ! ERROR: Out of memory
1044 call errmgr%report_error("solve_qr_pivot_mtx_cmplx", &
1045 "Insufficient memory available.", &
1046 la_out_of_memory_error)
1047 return
1048 end if
1049 wptr => wrk
1050 end if
1051
1052 ! Determine the rank of R11 using an incremental condition estimation
1053 wptr(ismin) = one
1054 wptr(ismax) = one
1055 smax = abs(a(1,1))
1056 smin = smax
1057 if (abs(a(1,1)) == zero) then
1058 rnk = 0
1059 b(1:maxmn,:) = zero
1060 return
1061 else
1062 rnk = 1
1063 end if
1064 do
1065 if (rnk < mn) then
1066 i = rnk + 1
1067 call zlaic1(imin, rnk, wptr(ismin:ismin+rnk-1), smin, &
1068 a(1:rnk-1,i), a(i,i), sminpr, s1, c1)
1069 call zlaic1(imax, rnk, wptr(ismax:ismax+rnk-1), smax, &
1070 a(1:rnk-1,i), a(i,i), smaxpr, s2, c2)
1071 if (smaxpr * rcond <= sminpr) then
1072 do i = 1, rnk
1073 wptr(ismin+i-1) = s1 * wptr(ismin+i-1)
1074 wptr(ismax+i-1) = s2 * wptr(ismax+i-1)
1075 end do
1076 wptr(ismin+rnk) = c1
1077 wptr(ismax+rnk) = c2
1078 smin = sminpr
1079 smax = smaxpr
1080 rnk = rnk + 1
1081 cycle
1082 end if
1083 end if
1084 exit
1085 end do
1086
1087 ! Partition R = [R11 R12]
1088 ! [ 0 R22]
1089 tau2 => wptr(1:rnk)
1090 w => wptr(rnk+1:lwork)
1091 if (rnk < n) call rz_factor(a(1:rnk,:), tau2, w)
1092
1093 ! Compute B(1:m,1:NRHS) = Q**T * B(1:M,1:NRHS)
1094 call mult_qr(.true., .true., a, tau, b(1:m,:), w)
1095
1096 ! Solve the triangular system T11 * B(1:rnk,1:nrhs) = B(1:rnk,1:nrhs)
1097 call solve_triangular_system(.true., .true., .false., .true., one, &
1098 a(1:rnk,1:rnk), b(1:rnk,:))
1099 if (n > rnk) b(rnk+1:n,:) = zero
1100
1101 ! Compute B(1:n,1:nrhs) = Y**T * B(1:n,1:nrhs)
1102 if (rnk < n) then
1103 call mult_rz(.true., .true., n - rnk, a(1:rnk,:), tau2, b(1:n,:), w)
1104 end if
1105
1106 ! Apply the pivoting: B(1:N,1:NRHS) = P * B(1:N,1:NRHS)
1107 do j = 1, nrhs
1108 do i = 1, n
1109 wptr(jpvt(i)) = b(i,j)
1110 end do
1111 b(:,j) = wptr(1:n)
1112 end do
1113
1114 ! Formatting
1115100 format(a, i0, a)
1116 end subroutine
1117
1118! ------------------------------------------------------------------------------
1119 module subroutine solve_qr_pivot_vec(a, tau, jpvt, b, work, olwork, err)
1120 ! Arguments
1121 real(real64), intent(inout), dimension(:,:) :: a
1122 real(real64), intent(in), dimension(:) :: tau
1123 integer(int32), intent(in), dimension(:) :: jpvt
1124 real(real64), intent(inout), dimension(:) :: b
1125 real(real64), intent(out), target, optional, dimension(:) :: work
1126 integer(int32), intent(out), optional :: olwork
1127 class(errors), intent(inout), optional, target :: err
1128
1129 ! Parameters
1130 integer(int32), parameter :: imin = 2
1131 integer(int32), parameter :: imax = 1
1132 real(real64), parameter :: zero = 0.0d0
1133 real(real64), parameter :: one = 1.0d0
1134
1135 ! Local Variables
1136 integer(int32) :: i, m, n, mn, lwork, ismin, ismax, rnk, maxmn, flag, &
1137 istat, lwork1, lwork2
1138 real(real64) :: rcond, smax, smin, smaxpr, sminpr, s1, c1, s2, c2
1139 real(real64), pointer, dimension(:) :: wptr, w, tau2
1140 real(real64), allocatable, target, dimension(:) :: wrk
1141 class(errors), pointer :: errmgr
1142 type(errors), target :: deferr
1143 character(len = :), allocatable :: errmsg
1144
1145 ! Initialization
1146 m = size(a, 1)
1147 n = size(a, 2)
1148 mn = min(m, n)
1149 maxmn = max(m, n)
1150 ismin = mn + 1
1151 ismax = 2 * mn + 1
1152 rcond = epsilon(rcond)
1153 if (present(err)) then
1154 errmgr => err
1155 else
1156 errmgr => deferr
1157 end if
1158
1159 ! Input Check
1160 flag = 0
1161 if (size(tau) /= mn) then
1162 flag = 2
1163 else if (size(jpvt) /= n) then
1164 flag = 3
1165 else if (size(b) /= maxmn) then
1166 flag = 4
1167 end if
1168 if (flag /= 0) then
1169 ! ERROR: One of the input arrays is not sized correctly
1170 allocate(character(len = 256) :: errmsg)
1171 write(errmsg, 100) "Input number ", flag, &
1172 " is not sized correctly."
1173 call errmgr%report_error("solve_qr_pivot_vec", trim(errmsg), &
1174 la_array_size_error)
1175 return
1176 end if
1177
1178 ! Workspace Query
1179 call rz_factor(a(1:mn,:), a(1:mn,1), olwork = lwork1)
1180 call mult_rz(.true., n, a(1:mn,:), a(1:mn,1), b(1:n), olwork = lwork2)
1181 lwork = max(lwork1, lwork2, 2 * mn + 1) + mn
1182 if (present(olwork)) then
1183 olwork = lwork
1184 return
1185 end if
1186
1187 ! Local Memory Allocation
1188 if (present(work)) then
1189 if (size(work) < lwork) then
1190 ! ERROR: WORK not sized correctly
1191 call errmgr%report_error("solve_qr_no_pivot_mtx", &
1192 "Incorrectly sized input array WORK, argument 5.", &
1193 la_array_size_error)
1194 return
1195 end if
1196 wptr => work(1:lwork)
1197 else
1198 allocate(wrk(lwork), stat = istat)
1199 if (istat /= 0) then
1200 ! ERROR: Out of memory
1201 call errmgr%report_error("solve_qr_pivot_vec", &
1202 "Insufficient memory available.", &
1203 la_out_of_memory_error)
1204 return
1205 end if
1206 wptr => wrk
1207 end if
1208
1209 ! Determine the rank of R11 using an incremental condition estimation
1210 wptr(ismin) = one
1211 wptr(ismax) = one
1212 smax = abs(a(1,1))
1213 smin = smax
1214 if (abs(a(1,1)) == zero) then
1215 rnk = 0
1216 b(maxmn) = zero
1217 return
1218 else
1219 rnk = 1
1220 end if
1221 do
1222 if (rnk < mn) then
1223 i = rnk + 1
1224 call dlaic1(imin, rnk, wptr(ismin:ismin+rnk-1), smin, &
1225 a(1:rnk-1,i), a(i,i), sminpr, s1, c1)
1226 call dlaic1(imax, rnk, wptr(ismax:ismax+rnk-1), smax, &
1227 a(1:rnk-1,i), a(i,i), smaxpr, s2, c2)
1228 if (smaxpr * rcond <= sminpr) then
1229 do i = 1, rnk
1230 wptr(ismin+i-1) = s1 * wptr(ismin+i-1)
1231 wptr(ismax+i-1) = s2 * wptr(ismax+i-1)
1232 end do
1233 wptr(ismin+rnk) = c1
1234 wptr(ismax+rnk) = c2
1235 smin = sminpr
1236 smax = smaxpr
1237 rnk = rnk + 1
1238 cycle
1239 end if
1240 end if
1241 exit
1242 end do
1243
1244 ! Partition R = [R11 R12]
1245 ! [ 0 R22]
1246 tau2 => wptr(1:rnk)
1247 w => wptr(rnk+1:lwork)
1248 if (rnk < n) call rz_factor(a(1:rnk,:), tau2, w)
1249
1250 ! Compute B(1:m,1:NRHS) = Q**T * B(1:M,1:NRHS)
1251 call mult_qr(.true., a, tau, b(1:m))
1252
1253 ! Solve the triangular system T11 * B(1:rnk) = B(1:rnk)
1254 call solve_triangular_system(.true., .false., .true., a(1:rnk,1:rnk), &
1255 b(1:rnk))
1256 if (n > rnk) b(rnk+1:n) = zero
1257
1258 ! Compute B(1:n) = Y**T * B(1:n)
1259 if (rnk < n) then
1260 call mult_rz(.true., n - rnk, a(1:rnk,:), tau2, b(1:n), w)
1261 end if
1262
1263 ! Apply the pivoting: B(1:N) = P * B(1:N)
1264 do i = 1, n
1265 wptr(jpvt(i)) = b(i)
1266 end do
1267 b = wptr(1:n)
1268
1269 ! Formatting
1270100 format(a, i0, a)
1271 end subroutine
1272
1273! ------------------------------------------------------------------------------
1274 module subroutine solve_qr_pivot_vec_cmplx(a, tau, jpvt, b, work, olwork, err)
1275 ! Arguments
1276 complex(real64), intent(inout), dimension(:,:) :: a
1277 complex(real64), intent(in), dimension(:) :: tau
1278 integer(int32), intent(in), dimension(:) :: jpvt
1279 complex(real64), intent(inout), dimension(:) :: b
1280 complex(real64), intent(out), target, optional, dimension(:) :: work
1281 integer(int32), intent(out), optional :: olwork
1282 class(errors), intent(inout), optional, target :: err
1283
1284 ! Parameters
1285 integer(int32), parameter :: imin = 2
1286 integer(int32), parameter :: imax = 1
1287 complex(real64), parameter :: zero = (0.0d0, 0.0d0)
1288 complex(real64), parameter :: one = (1.0d0, 0.0d0)
1289
1290 ! Local Variables
1291 integer(int32) :: i, m, n, mn, lwork, ismin, ismax, rnk, maxmn, flag, &
1292 istat, lwork1, lwork2
1293 real(real64) :: rcond, smax, smin, smaxpr, sminpr
1294 complex(real64) :: s1, c1, s2, c2
1295 complex(real64), pointer, dimension(:) :: wptr, w, tau2
1296 complex(real64), allocatable, target, dimension(:) :: wrk
1297 class(errors), pointer :: errmgr
1298 type(errors), target :: deferr
1299 character(len = :), allocatable :: errmsg
1300
1301 ! Initialization
1302 m = size(a, 1)
1303 n = size(a, 2)
1304 mn = min(m, n)
1305 maxmn = max(m, n)
1306 ismin = mn + 1
1307 ismax = 2 * mn + 1
1308 rcond = epsilon(rcond)
1309 if (present(err)) then
1310 errmgr => err
1311 else
1312 errmgr => deferr
1313 end if
1314
1315 ! Input Check
1316 flag = 0
1317 if (size(tau) /= mn) then
1318 flag = 2
1319 else if (size(jpvt) /= n) then
1320 flag = 3
1321 else if (size(b) /= maxmn) then
1322 flag = 4
1323 end if
1324 if (flag /= 0) then
1325 ! ERROR: One of the input arrays is not sized correctly
1326 allocate(character(len = 256) :: errmsg)
1327 write(errmsg, 100) "Input number ", flag, &
1328 " is not sized correctly."
1329 call errmgr%report_error("solve_qr_pivot_vec_cmplx", trim(errmsg), &
1330 la_array_size_error)
1331 return
1332 end if
1333
1334 ! Workspace Query
1335 call rz_factor(a(1:mn,:), a(1:mn,1), olwork = lwork1)
1336 call mult_rz(.true., n, a(1:mn,:), a(1:mn,1), b(1:n), olwork = lwork2)
1337 lwork = max(lwork1, lwork2, 2 * mn + 1) + mn
1338 if (present(olwork)) then
1339 olwork = lwork
1340 return
1341 end if
1342
1343 ! Local Memory Allocation
1344 if (present(work)) then
1345 if (size(work) < lwork) then
1346 ! ERROR: WORK not sized correctly
1347 call errmgr%report_error("solve_qr_no_pivot_mtx_cmplx", &
1348 "Incorrectly sized input array WORK, argument 5.", &
1349 la_array_size_error)
1350 return
1351 end if
1352 wptr => work(1:lwork)
1353 else
1354 allocate(wrk(lwork), stat = istat)
1355 if (istat /= 0) then
1356 ! ERROR: Out of memory
1357 call errmgr%report_error("solve_qr_pivot_vec_cmplx", &
1358 "Insufficient memory available.", &
1359 la_out_of_memory_error)
1360 return
1361 end if
1362 wptr => wrk
1363 end if
1364
1365 ! Determine the rank of R11 using an incremental condition estimation
1366 wptr(ismin) = one
1367 wptr(ismax) = one
1368 smax = abs(a(1,1))
1369 smin = smax
1370 if (abs(a(1,1)) == zero) then
1371 rnk = 0
1372 b(maxmn) = zero
1373 return
1374 else
1375 rnk = 1
1376 end if
1377 do
1378 if (rnk < mn) then
1379 i = rnk + 1
1380 call zlaic1(imin, rnk, wptr(ismin:ismin+rnk-1), smin, &
1381 a(1:rnk-1,i), a(i,i), sminpr, s1, c1)
1382 call zlaic1(imax, rnk, wptr(ismax:ismax+rnk-1), smax, &
1383 a(1:rnk-1,i), a(i,i), smaxpr, s2, c2)
1384 if (smaxpr * rcond <= sminpr) then
1385 do i = 1, rnk
1386 wptr(ismin+i-1) = s1 * wptr(ismin+i-1)
1387 wptr(ismax+i-1) = s2 * wptr(ismax+i-1)
1388 end do
1389 wptr(ismin+rnk) = c1
1390 wptr(ismax+rnk) = c2
1391 smin = sminpr
1392 smax = smaxpr
1393 rnk = rnk + 1
1394 cycle
1395 end if
1396 end if
1397 exit
1398 end do
1399
1400 ! Partition R = [R11 R12]
1401 ! [ 0 R22]
1402 tau2 => wptr(1:rnk)
1403 w => wptr(rnk+1:lwork)
1404 if (rnk < n) call rz_factor(a(1:rnk,:), tau2, w)
1405
1406 ! Compute B(1:m,1:NRHS) = Q**T * B(1:M,1:NRHS)
1407 call mult_qr(.true., a, tau, b(1:m))
1408
1409 ! Solve the triangular system T11 * B(1:rnk) = B(1:rnk)
1410 call solve_triangular_system(.true., .false., .true., a(1:rnk,1:rnk), &
1411 b(1:rnk))
1412 if (n > rnk) b(rnk+1:n) = zero
1413
1414 ! Compute B(1:n) = Y**T * B(1:n)
1415 if (rnk < n) then
1416 call mult_rz(.true., n - rnk, a(1:rnk,:), tau2, b(1:n), w)
1417 end if
1418
1419 ! Apply the pivoting: B(1:N) = P * B(1:N)
1420 do i = 1, n
1421 wptr(jpvt(i)) = b(i)
1422 end do
1423 b = wptr(1:n)
1424
1425 ! Formatting
1426100 format(a, i0, a)
1427 end subroutine
1428
1429! ******************************************************************************
1430! CHOLESKY SOLVE
1431! ------------------------------------------------------------------------------
1432 module subroutine solve_cholesky_mtx(upper, a, b, err)
1433 ! Arguments
1434 logical, intent(in) :: upper
1435 real(real64), intent(in), dimension(:,:) :: a
1436 real(real64), intent(inout), dimension(:,:) :: b
1437 class(errors), intent(inout), optional, target :: err
1438
1439 ! Local Variables
1440 character :: uplo
1441 integer(int32) :: n, nrhs, flag
1442 class(errors), pointer :: errmgr
1443 type(errors), target :: deferr
1444 character(len = :), allocatable :: errmsg
1445
1446 ! Initialization
1447 n = size(a, 1)
1448 nrhs = size(b, 2)
1449 if (upper) then
1450 uplo = 'U'
1451 else
1452 uplo = 'L'
1453 end if
1454 if (present(err)) then
1455 errmgr => err
1456 else
1457 errmgr => deferr
1458 end if
1459
1460 ! Input Check
1461 flag = 0
1462 if (size(a, 2) /= n) then
1463 flag = 2
1464 else if (size(b, 1) /= n) then
1465 flag = 3
1466 end if
1467 if (flag /= 0) then
1468 allocate(character(len = 256) :: errmsg)
1469 write(errmsg, 100) "Input number ", flag, &
1470 " is not sized correctly."
1471 call errmgr%report_error("solve_cholesky_mtx", trim(errmsg), &
1472 la_array_size_error)
1473 return
1474 end if
1475
1476 ! Process
1477 call dpotrs(uplo, n, nrhs, a, n, b, n, flag)
1478
1479 ! Formatting
1480100 format(a, i0, a)
1481 end subroutine
1482
1483! ------------------------------------------------------------------------------
1484 module subroutine solve_cholesky_mtx_cmplx(upper, a, b, err)
1485 ! Arguments
1486 logical, intent(in) :: upper
1487 complex(real64), intent(in), dimension(:,:) :: a
1488 complex(real64), intent(inout), dimension(:,:) :: b
1489 class(errors), intent(inout), optional, target :: err
1490
1491 ! Local Variables
1492 character :: uplo
1493 integer(int32) :: n, nrhs, flag
1494 class(errors), pointer :: errmgr
1495 type(errors), target :: deferr
1496 character(len = :), allocatable :: errmsg
1497
1498 ! Initialization
1499 n = size(a, 1)
1500 nrhs = size(b, 2)
1501 if (upper) then
1502 uplo = 'U'
1503 else
1504 uplo = 'L'
1505 end if
1506 if (present(err)) then
1507 errmgr => err
1508 else
1509 errmgr => deferr
1510 end if
1511
1512 ! Input Check
1513 flag = 0
1514 if (size(a, 2) /= n) then
1515 flag = 2
1516 else if (size(b, 1) /= n) then
1517 flag = 3
1518 end if
1519 if (flag /= 0) then
1520 allocate(character(len = 256) :: errmsg)
1521 write(errmsg, 100) "Input number ", flag, &
1522 " is not sized correctly."
1523 call errmgr%report_error("solve_cholesky_mtx_cmplx", trim(errmsg), &
1524 la_array_size_error)
1525 return
1526 end if
1527
1528 ! Process
1529 call zpotrs(uplo, n, nrhs, a, n, b, n, flag)
1530
1531 ! Formatting
1532100 format(a, i0, a)
1533 end subroutine
1534
1535! ------------------------------------------------------------------------------
1536 module subroutine solve_cholesky_vec(upper, a, b, err)
1537 ! Arguments
1538 logical, intent(in) :: upper
1539 real(real64), intent(in), dimension(:,:) :: a
1540 real(real64), intent(inout), dimension(:) :: b
1541 class(errors), intent(inout), optional, target :: err
1542
1543 ! Local Variables
1544 character :: uplo
1545 integer(int32) :: n, flag
1546 class(errors), pointer :: errmgr
1547 type(errors), target :: deferr
1548 character(len = :), allocatable :: errmsg
1549
1550 ! Initialization
1551 n = size(a, 1)
1552 if (upper) then
1553 uplo = 'U'
1554 else
1555 uplo = 'L'
1556 end if
1557 if (present(err)) then
1558 errmgr => err
1559 else
1560 errmgr => deferr
1561 end if
1562
1563 ! Input Check
1564 flag = 0
1565 if (size(a, 2) /= n) then
1566 flag = 2
1567 else if (size(b) /= n) then
1568 flag = 3
1569 end if
1570 if (flag /= 0) then
1571 allocate(character(len = 256) :: errmsg)
1572 write(errmsg, 100) "Input number ", flag, &
1573 " is not sized correctly."
1574 call errmgr%report_error("solve_cholesky_vec", trim(errmsg), &
1575 la_array_size_error)
1576 return
1577 end if
1578
1579 ! Process
1580 call dpotrs(uplo, n, 1, a, n, b, n, flag)
1581
1582 ! Formatting
1583100 format(a, i0, a)
1584 end subroutine
1585
1586! ------------------------------------------------------------------------------
1587 module subroutine solve_cholesky_vec_cmplx(upper, a, b, err)
1588 ! Arguments
1589 logical, intent(in) :: upper
1590 complex(real64), intent(in), dimension(:,:) :: a
1591 complex(real64), intent(inout), dimension(:) :: b
1592 class(errors), intent(inout), optional, target :: err
1593
1594 ! Local Variables
1595 character :: uplo
1596 integer(int32) :: n, flag
1597 class(errors), pointer :: errmgr
1598 type(errors), target :: deferr
1599 character(len = :), allocatable :: errmsg
1600
1601 ! Initialization
1602 n = size(a, 1)
1603 if (upper) then
1604 uplo = 'U'
1605 else
1606 uplo = 'L'
1607 end if
1608 if (present(err)) then
1609 errmgr => err
1610 else
1611 errmgr => deferr
1612 end if
1613
1614 ! Input Check
1615 flag = 0
1616 if (size(a, 2) /= n) then
1617 flag = 2
1618 else if (size(b) /= n) then
1619 flag = 3
1620 end if
1621 if (flag /= 0) then
1622 allocate(character(len = 256) :: errmsg)
1623 write(errmsg, 100) "Input number ", flag, &
1624 " is not sized correctly."
1625 call errmgr%report_error("solve_cholesky_vec_cmplx", trim(errmsg), &
1626 la_array_size_error)
1627 return
1628 end if
1629
1630 ! Process
1631 call zpotrs(uplo, n, 1, a, n, b, n, flag)
1632
1633 ! Formatting
1634100 format(a, i0, a)
1635 end subroutine
1636
1637! ******************************************************************************
1638! MATRIX INVERSION ROUTINES
1639! ------------------------------------------------------------------------------
1640 module subroutine mtx_inverse_dbl(a, iwork, work, olwork, err)
1641 ! Arguments
1642 real(real64), intent(inout), dimension(:,:) :: a
1643 integer(int32), intent(out), target, optional, dimension(:) :: iwork
1644 real(real64), intent(out), target, optional, dimension(:) :: work
1645 integer(int32), intent(out), optional :: olwork
1646 class(errors), intent(inout), optional, target :: err
1647
1648 ! Local Variables
1649 integer(int32) :: n, liwork, lwork, istat, flag, itemp(1)
1650 integer(int32), pointer, dimension(:) :: iptr
1651 integer(int32), allocatable, target, dimension(:) :: iwrk
1652 real(real64), pointer, dimension(:) :: wptr
1653 real(real64), allocatable, target, dimension(:) :: wrk
1654 real(real64), dimension(1) :: temp
1655 class(errors), pointer :: errmgr
1656 type(errors), target :: deferr
1657
1658 ! Initialization
1659 n = size(a, 1)
1660 liwork = n
1661 if (present(err)) then
1662 errmgr => err
1663 else
1664 errmgr => deferr
1665 end if
1666
1667 ! Input Check
1668 if (size(a, 2) /= n) then
1669 call errmgr%report_error("mtx_inverse", &
1670 "The matrix must be squre to invert.", la_array_size_error)
1671 return
1672 end if
1673
1674 ! Workspace Query
1675 call dgetri(n, a, n, itemp, temp, -1, flag)
1676 lwork = int(temp(1), int32)
1677 if (present(olwork)) then
1678 olwork = lwork
1679 return
1680 end if
1681
1682 ! Workspace Allocation
1683 if (present(work)) then
1684 if (size(work) < lwork) then
1685 ! ERROR: WORK not sized correctly
1686 call errmgr%report_error("mtx_inverse_dbl", &
1687 "Incorrectly sized input array WORK, argument 3.", &
1688 la_array_size_error)
1689 return
1690 end if
1691 wptr => work(1:lwork)
1692 else
1693 allocate(wrk(lwork), stat = istat)
1694 if (istat /= 0) then
1695 ! ERROR: Out of memory
1696 call errmgr%report_error("mtx_inverse_dbl", &
1697 "Insufficient memory available.", &
1698 la_out_of_memory_error)
1699 return
1700 end if
1701 wptr => wrk
1702 end if
1703
1704 ! Integer Workspace Allocation
1705 if (present(iwork)) then
1706 if (size(iwork) < liwork) then
1707 ! ERROR: IWORK not sized correctly
1708 call errmgr%report_error("mtx_inverse_dbl", &
1709 "Incorrectly sized input array IWORK, argument 2.", &
1710 la_array_size_error)
1711 return
1712 end if
1713 iptr => iwork(1:liwork)
1714 else
1715 allocate(iwrk(liwork), stat = istat)
1716 if (istat /= 0) then
1717 ! ERROR: Out of memory
1718 call errmgr%report_error("mtx_inverse_dbl", &
1719 "Insufficient memory available.", &
1720 la_out_of_memory_error)
1721 return
1722 end if
1723 iptr => iwrk
1724 end if
1725
1726 ! Compute the LU factorization of A
1727 call dgetrf(n, n, a, n, iptr, flag)
1728
1729 ! Compute the inverse of the LU factored matrix
1730 call dgetri(n, a, n, iptr, wptr, lwork, flag)
1731
1732 ! Check for a singular matrix
1733 if (flag > 0) then
1734 call errmgr%report_error("mtx_inverse_dbl", &
1735 "The matrix is singular; therefore, the inverse could " // &
1736 "not be computed.", la_singular_matrix_error)
1737 end if
1738 end subroutine
1739
1740! ------------------------------------------------------------------------------
1741 module subroutine mtx_inverse_cmplx(a, iwork, work, olwork, err)
1742 ! Arguments
1743 complex(real64), intent(inout), dimension(:,:) :: a
1744 integer(int32), intent(out), target, optional, dimension(:) :: iwork
1745 complex(real64), intent(out), target, optional, dimension(:) :: work
1746 integer(int32), intent(out), optional :: olwork
1747 class(errors), intent(inout), optional, target :: err
1748
1749 ! Local Variables
1750 integer(int32) :: n, liwork, lwork, istat, flag, itemp(1)
1751 integer(int32), pointer, dimension(:) :: iptr
1752 integer(int32), allocatable, target, dimension(:) :: iwrk
1753 complex(real64), pointer, dimension(:) :: wptr
1754 complex(real64), allocatable, target, dimension(:) :: wrk
1755 complex(real64), dimension(1) :: temp
1756 class(errors), pointer :: errmgr
1757 type(errors), target :: deferr
1758
1759 ! Initialization
1760 n = size(a, 1)
1761 liwork = n
1762 if (present(err)) then
1763 errmgr => err
1764 else
1765 errmgr => deferr
1766 end if
1767
1768 ! Input Check
1769 if (size(a, 2) /= n) then
1770 call errmgr%report_error("mtx_inverse_cmplx", &
1771 "The matrix must be squre to invert.", la_array_size_error)
1772 return
1773 end if
1774
1775 ! Workspace Query
1776 call zgetri(n, a, n, itemp, temp, -1, flag)
1777 lwork = int(temp(1), int32)
1778 if (present(olwork)) then
1779 olwork = lwork
1780 return
1781 end if
1782
1783 ! Workspace Allocation
1784 if (present(work)) then
1785 if (size(work) < lwork) then
1786 ! ERROR: WORK not sized correctly
1787 call errmgr%report_error("mtx_inverse_cmplx", &
1788 "Incorrectly sized input array WORK, argument 3.", &
1789 la_array_size_error)
1790 return
1791 end if
1792 wptr => work(1:lwork)
1793 else
1794 allocate(wrk(lwork), stat = istat)
1795 if (istat /= 0) then
1796 ! ERROR: Out of memory
1797 call errmgr%report_error("mtx_inverse_cmplx", &
1798 "Insufficient memory available.", &
1799 la_out_of_memory_error)
1800 return
1801 end if
1802 wptr => wrk
1803 end if
1804
1805 ! Integer Workspace Allocation
1806 if (present(iwork)) then
1807 if (size(iwork) < liwork) then
1808 ! ERROR: IWORK not sized correctly
1809 call errmgr%report_error("mtx_inverse_cmplx", &
1810 "Incorrectly sized input array IWORK, argument 2.", &
1811 la_array_size_error)
1812 return
1813 end if
1814 iptr => iwork(1:liwork)
1815 else
1816 allocate(iwrk(liwork), stat = istat)
1817 if (istat /= 0) then
1818 ! ERROR: Out of memory
1819 call errmgr%report_error("mtx_inverse_cmplx", &
1820 "Insufficient memory available.", &
1821 la_out_of_memory_error)
1822 return
1823 end if
1824 iptr => iwrk
1825 end if
1826
1827 ! Compute the LU factorization of A
1828 call zgetrf(n, n, a, n, iptr, flag)
1829
1830 ! Compute the inverse of the LU factored matrix
1831 call zgetri(n, a, n, iptr, wptr, lwork, flag)
1832
1833 ! Check for a singular matrix
1834 if (flag > 0) then
1835 call errmgr%report_error("mtx_inverse_cmplx", &
1836 "The matrix is singular; therefore, the inverse could " // &
1837 "not be computed.", la_singular_matrix_error)
1838 end if
1839 end subroutine
1840
1841! ------------------------------------------------------------------------------
1842 module subroutine mtx_pinverse_dbl(a, ainv, tol, work, olwork, err)
1843 ! Arguments
1844 real(real64), intent(inout), dimension(:,:) :: a
1845 real(real64), intent(out), dimension(:,:) :: ainv
1846 real(real64), intent(in), optional :: tol
1847 real(real64), intent(out), target, dimension(:), optional :: work
1848 integer(int32), intent(out), optional :: olwork
1849 class(errors), intent(inout), optional, target :: err
1850
1851 ! External Function Interfaces
1852 interface
1853 function dlamch(cmach) result(x)
1854 use, intrinsic :: iso_fortran_env, only : real64
1855 character, intent(in) :: cmach
1856 real(real64) :: x
1857 end function
1858 end interface
1859
1860 ! Parameters
1861 real(real64), parameter :: zero = 0.0d0
1862 real(real64), parameter :: one = 1.0d0
1863
1864 ! Local Variables
1865 integer(int32) :: i, m, n, mn, lwork, istat, flag, i1, i2a, i2b, i3a, &
1866 i3b, i4
1867 real(real64), pointer, dimension(:) :: s, wptr, w
1868 real(real64), pointer, dimension(:,:) :: u, vt
1869 real(real64), allocatable, target, dimension(:) :: wrk
1870 real(real64), dimension(1) :: temp
1871 real(real64) :: t, tref, tolcheck
1872 class(errors), pointer :: errmgr
1873 type(errors), target :: deferr
1874 character(len = :), allocatable :: errmsg
1875
1876 ! Initialization
1877 m = size(a, 1)
1878 n = size(a, 2)
1879 mn = min(m, n)
1880 i1 = m * mn
1881 i2a = i1 + 1
1882 i2b = i2a + n * n - 1
1883 i3a = i2b + 1
1884 i3b = i3a + mn - 1
1885 i4 = i3b + 1
1886 tolcheck = dlamch('s')
1887 if (present(err)) then
1888 errmgr => err
1889 else
1890 errmgr => deferr
1891 end if
1892
1893 ! Input Check
1894 if (size(ainv, 1) /= n .or. size(ainv, 2) /= m) then
1895 allocate(character(len = 256) :: errmsg)
1896 write(errmsg, 100) &
1897 "The output matrix AINV is not sized appropriately. " // &
1898 "It is expected to be ", n, "-by-", m, "."
1899 call errmgr%report_error("mtx_pinverse", errmsg, &
1900 la_array_size_error)
1901 return
1902 end if
1903
1904 ! Workspace Query
1905 call dgesvd('S', 'A', m, n, a, m, a(1:mn,:), a, m, a, n, temp, -1, flag)
1906 lwork = int(temp(1), int32)
1907 lwork = lwork + m * mn + n * n + mn
1908 if (present(olwork)) then
1909 olwork = lwork
1910 return
1911 end if
1912
1913 ! Local Memory Allocation
1914 if (present(work)) then
1915 if (size(work) < lwork) then
1916 ! ERROR: WORK not sized correctly
1917 call errmgr%report_error("mtx_pinverse", &
1918 "Incorrectly sized input array WORK, argument 4.", &
1919 la_array_size_error)
1920 return
1921 end if
1922 wptr => work(1:lwork)
1923 else
1924 allocate(wrk(lwork), stat = istat)
1925 if (istat /= 0) then
1926 ! ERROR: Out of memory
1927 call errmgr%report_error("mtx_pinverse", &
1928 "Insufficient memory available.", &
1929 la_out_of_memory_error)
1930 return
1931 end if
1932 wptr => wrk
1933 end if
1934 u(1:m,1:mn) => wptr(1:i1)
1935 vt(1:n,1:n) => wptr(i2a:i2b)
1936 s => wptr(i3a:i3b)
1937 w => wptr(i4:lwork)
1938
1939 ! Compute the SVD of A
1940 call dgesvd('S', 'A', m, n, a, m, s, u, m, vt, n, w, size(w), flag)
1941
1942 ! Check for convergence
1943 if (flag > 0) then
1944 allocate(character(len = 256) :: errmsg)
1945 write(errmsg, 101) flag, " superdiagonals could not " // &
1946 "converge to zero as part of the QR iteration process."
1947 call errmgr%report_warning("mtx_pinverse", errmsg, &
1948 la_convergence_error)
1949 return
1950 end if
1951
1952 ! Determine the threshold tolerance for the singular values such that
1953 ! singular values less than the threshold result in zero when inverted.
1954 tref = max(m, n) * epsilon(t) * s(1)
1955 if (present(tol)) then
1956 t = tol
1957 else
1958 t = tref
1959 end if
1960 !if (t < safe_denom(t)) then
1961 if (t < tolcheck) then
1962 ! The supplied tolerance is too small, simply fall back to the
1963 ! default, but issue a warning to the user
1964 t = tref
1965 ! call errmgr%report_warning("pinverse_1", "The supplied tolerance was " // &
1966 ! "smaller than a value that would result in an overflow " // &
1967 ! "condition, or is negative; therefore, the tolerance has " // &
1968 ! "been reset to its default value.")
1969 end if
1970
1971 ! Compute the pseudoinverse such that pinv(A) = V * inv(S) * U**T by
1972 ! first computing V * inv(S) (result is N-by-M), and store in the first
1973 ! MN rows of VT in a transposed manner.
1974 do i = 1, mn
1975 ! Apply 1 / S(I) to VT(I,:)
1976 if (s(i) < t) then
1977 vt(i,:) = zero
1978 else
1979 call recip_mult_array(s(i), vt(i,1:n))
1980 end if
1981 end do
1982
1983 ! Compute (VT**T * inv(S)) * U**T
1984 call mtx_mult(.true., .true., one, vt(1:mn,:), u, zero, ainv)
1985
1986 ! Formatting
1987100 format(a, i0, a, i0, a)
1988101 format(i0, a)
1989 end subroutine
1990
1991! ------------------------------------------------------------------------------
1992 module subroutine mtx_pinverse_cmplx(a, ainv, tol, work, olwork, rwork, err)
1993 ! Arguments
1994 complex(real64), intent(inout), dimension(:,:) :: a
1995 complex(real64), intent(out), dimension(:,:) :: ainv
1996 real(real64), intent(in), optional :: tol
1997 complex(real64), intent(out), target, dimension(:), optional :: work
1998 integer(int32), intent(out), optional :: olwork
1999 real(real64), intent(out), target, dimension(:), optional :: rwork
2000 class(errors), intent(inout), optional, target :: err
2001
2002 ! External Function Interfaces
2003 interface
2004 function dlamch(cmach) result(x)
2005 use, intrinsic :: iso_fortran_env, only : real64
2006 character, intent(in) :: cmach
2007 real(real64) :: x
2008 end function
2009 end interface
2010
2011 ! Parameters
2012 complex(real64), parameter :: zero = (0.0d0, 0.0d0)
2013 complex(real64), parameter :: one = (1.0d0, 0.0d0)
2014
2015 ! Local Variables
2016 integer(int32) :: i, m, n, mn, lwork, istat, flag, i1, i2a, i2b, i3, &
2017 lrwork, j, k
2018 real(real64), pointer, dimension(:) :: s, rwptr, rw
2019 real(real64), allocatable, target, dimension(:) :: rwrk
2020 complex(real64), pointer, dimension(:) :: wptr, w
2021 complex(real64), pointer, dimension(:,:) :: u, vt
2022 complex(real64), allocatable, target, dimension(:) :: wrk
2023 complex(real64) :: temp(1), val
2024 real(real64) :: t, tref, tolcheck, rtemp(1)
2025 class(errors), pointer :: errmgr
2026 type(errors), target :: deferr
2027 character(len = :), allocatable :: errmsg
2028
2029 ! Initialization
2030 m = size(a, 1)
2031 n = size(a, 2)
2032 mn = min(m, n)
2033 lrwork = 6 * mn
2034 i1 = m * mn
2035 i2a = i1 + 1
2036 i2b = i2a + n * n - 1
2037 i3 = i2b + 1
2038 tolcheck = dlamch('s')
2039 if (present(err)) then
2040 errmgr => err
2041 else
2042 errmgr => deferr
2043 end if
2044
2045 ! Input Check
2046 if (size(ainv, 1) /= n .or. size(ainv, 2) /= m) then
2047 allocate(character(len = 256) :: errmsg)
2048 write(errmsg, 100) &
2049 "The output matrix AINV is not sized appropriately. " // &
2050 "It is expected to be ", n, "-by-", m, "."
2051 call errmgr%report_error("mtx_pinverse_cmplx", errmsg, &
2052 la_array_size_error)
2053 return
2054 end if
2055
2056 ! Workspace Query
2057 call zgesvd('S', 'A', m, n, a, m, rtemp, a, m, a, n, temp, -1, &
2058 rtemp, flag)
2059 lwork = int(temp(1), int32)
2060 lwork = lwork + m * mn + n * n
2061 if (present(olwork)) then
2062 olwork = lwork
2063 return
2064 end if
2065
2066 ! Local Memory Allocation
2067 if (present(work)) then
2068 if (size(work) < lwork) then
2069 ! ERROR: WORK not sized correctly
2070 call errmgr%report_error("mtx_pinverse_cmplx", &
2071 "Incorrectly sized input array WORK, argument 4.", &
2072 la_array_size_error)
2073 return
2074 end if
2075 wptr => work(1:lwork)
2076 else
2077 allocate(wrk(lwork), stat = istat)
2078 if (istat /= 0) then
2079 ! ERROR: Out of memory
2080 call errmgr%report_error("mtx_pinverse_cmplx", &
2081 "Insufficient memory available.", &
2082 la_out_of_memory_error)
2083 return
2084 end if
2085 wptr => wrk
2086 end if
2087
2088 if (present(rwork)) then
2089 if (size(rwork) < lrwork) then
2090 ! ERROR: WORK not sized correctly
2091 call errmgr%report_error("mtx_pinverse_cmplx", &
2092 "Incorrectly sized input array RWORK, argument 6.", &
2093 la_array_size_error)
2094 return
2095 end if
2096 rwptr => rwork(1:lrwork)
2097 else
2098 allocate(rwrk(lrwork), stat = istat)
2099 if (istat /= 0) then
2100 ! ERROR: Out of memory
2101 call errmgr%report_error("mtx_pinverse_cmplx", &
2102 "Insufficient memory available.", &
2103 la_out_of_memory_error)
2104 return
2105 end if
2106 rwptr => rwrk
2107 end if
2108 u(1:m,1:mn) => wptr(1:i1)
2109 vt(1:n,1:n) => wptr(i2a:i2b)
2110 w => wptr(i3:lwork)
2111 s => rwptr(1:mn)
2112 rw => rwptr(mn+1:lrwork)
2113
2114 ! Compute the SVD of A
2115 call zgesvd('S', 'A', m, n, a, m, s, u, m, vt, n, w, size(w), rw, flag)
2116
2117 ! Check for convergence
2118 if (flag > 0) then
2119 allocate(character(len = 256) :: errmsg)
2120 write(errmsg, 101) flag, " superdiagonals could not " // &
2121 "converge to zero as part of the QR iteration process."
2122 call errmgr%report_warning("mtx_pinverse_cmplx", errmsg, &
2123 la_convergence_error)
2124 return
2125 end if
2126
2127 ! Determine the threshold tolerance for the singular values such that
2128 ! singular values less than the threshold result in zero when inverted.
2129 tref = max(m, n) * epsilon(t) * s(1)
2130 if (present(tol)) then
2131 t = tol
2132 else
2133 t = tref
2134 end if
2135 !if (t < safe_denom(t)) then
2136 if (t < tolcheck) then
2137 ! The supplied tolerance is too small, simply fall back to the
2138 ! default, but issue a warning to the user
2139 t = tref
2140 ! call errmgr%report_warning("pinverse_1", "The supplied tolerance was " // &
2141 ! "smaller than a value that would result in an overflow " // &
2142 ! "condition, or is negative; therefore, the tolerance has " // &
2143 ! "been reset to its default value.")
2144 end if
2145
2146 ! Compute the pseudoinverse such that pinv(A) = V * inv(S) * U**T by
2147 ! first computing V * inv(S) (result is N-by-M), and store in the first
2148 ! MN rows of VT in a transposed manner.
2149 do i = 1, mn
2150 ! Apply 1 / S(I) to VT(I,:)
2151 if (s(i) < t) then
2152 vt(i,:) = zero
2153 else
2154 ! call recip_mult_array(s(i), vt(i,1:n))
2155 vt(i,1:n) = conjg(vt(i,1:n)) / s(i)
2156 end if
2157 end do
2158
2159 ! Compute (VT**T * inv(S)) * U**H
2160 ! ainv = n-by-m
2161 ! vt is n-by-n
2162 ! u is m-by-mn such that u**H = mn-by-m
2163 ! Compute ainv = vt**T * u**H
2164 do j = 1, m
2165 do i = 1, n
2166 val = zero
2167 do k = 1, mn
2168 val = val + vt(k,i) * conjg(u(j,k))
2169 end do
2170 ainv(i,j) = val
2171 end do
2172 end do
2173
2174 ! Formatting
2175100 format(a, i0, a, i0, a)
2176101 format(i0, a)
2177 end subroutine
2178
2179! ******************************************************************************
2180! LEAST SQUARES SOLUTION ROUTINES
2181! ------------------------------------------------------------------------------
2182 module subroutine solve_least_squares_mtx(a, b, work, olwork, err)
2183 ! Arguments
2184 real(real64), intent(inout), dimension(:,:) :: a, b
2185 real(real64), intent(out), target, optional, dimension(:) :: work
2186 integer(int32), intent(out), optional :: olwork
2187 class(errors), intent(inout), optional, target :: err
2188
2189 ! Local Variables
2190 integer(int32) :: m, n, maxmn, nrhs, lwork, istat, flag
2191 real(real64), pointer, dimension(:) :: wptr
2192 real(real64), allocatable, target, dimension(:) :: wrk
2193 real(real64), dimension(1) :: temp
2194 class(errors), pointer :: errmgr
2195 type(errors), target :: deferr
2196
2197 ! Initialization
2198 m = size(a, 1)
2199 n = size(a, 2)
2200 maxmn = max(m, n)
2201 nrhs = size(b, 2)
2202 if (present(err)) then
2203 errmgr => err
2204 else
2205 errmgr => deferr
2206 end if
2207
2208 ! Input Check
2209 if (size(b, 1) /= maxmn) then
2210 call errmgr%report_error("solve_least_squares_mtx", &
2211 "Input 2 is not sized correctly.", la_array_size_error)
2212 return
2213 end if
2214
2215 ! Workspace Query
2216 call dgels('N', m, n, nrhs, a, m, b, maxmn, temp, -1, flag)
2217 lwork = int(temp(1), int32)
2218 if (present(olwork)) then
2219 olwork = lwork
2220 return
2221 end if
2222
2223 ! Local Memory Allocation
2224 if (present(work)) then
2225 if (size(work) < lwork) then
2226 ! ERROR: WORK not sized correctly
2227 call errmgr%report_error("solve_least_squares_mtx", &
2228 "Incorrectly sized input array WORK, argument 3.", &
2229 la_array_size_error)
2230 return
2231 end if
2232 wptr => work(1:lwork)
2233 else
2234 allocate(wrk(lwork), stat = istat)
2235 if (istat /= 0) then
2236 ! ERROR: Out of memory
2237 call errmgr%report_error("solve_least_squares_mtx", &
2238 "Insufficient memory available.", &
2239 la_out_of_memory_error)
2240 return
2241 end if
2242 wptr => wrk
2243 end if
2244
2245 ! Process
2246 call dgels('N', m, n, nrhs, a, m, b, maxmn, wptr, lwork, flag)
2247 if (flag > 0) then
2248 call errmgr%report_error("solve_least_squares_mtx", &
2249 "The supplied matrix is not of full rank; therefore, " // &
2250 "the solution could not be computed via this routine. " // &
2251 "Try a routine that utilizes column pivoting.", &
2252 la_invalid_operation_error)
2253 end if
2254 end subroutine
2255
2256! ------------------------------------------------------------------------------
2257 module subroutine solve_least_squares_mtx_cmplx(a, b, work, olwork, err)
2258 ! Arguments
2259 complex(real64), intent(inout), dimension(:,:) :: a, b
2260 complex(real64), intent(out), target, optional, dimension(:) :: work
2261 integer(int32), intent(out), optional :: olwork
2262 class(errors), intent(inout), optional, target :: err
2263
2264 ! Local Variables
2265 integer(int32) :: m, n, maxmn, nrhs, lwork, istat, flag
2266 complex(real64), pointer, dimension(:) :: wptr
2267 complex(real64), allocatable, target, dimension(:) :: wrk
2268 complex(real64), dimension(1) :: temp
2269 class(errors), pointer :: errmgr
2270 type(errors), target :: deferr
2271
2272 ! Initialization
2273 m = size(a, 1)
2274 n = size(a, 2)
2275 maxmn = max(m, n)
2276 nrhs = size(b, 2)
2277 if (present(err)) then
2278 errmgr => err
2279 else
2280 errmgr => deferr
2281 end if
2282
2283 ! Input Check
2284 if (size(b, 1) /= maxmn) then
2285 call errmgr%report_error("solve_least_squares_mtx_cmplx", &
2286 "Input 2 is not sized correctly.", la_array_size_error)
2287 return
2288 end if
2289
2290 ! Workspace Query
2291 call zgels('N', m, n, nrhs, a, m, b, maxmn, temp, -1, flag)
2292 lwork = int(temp(1), int32)
2293 if (present(olwork)) then
2294 olwork = lwork
2295 return
2296 end if
2297
2298 ! Local Memory Allocation
2299 if (present(work)) then
2300 if (size(work) < lwork) then
2301 ! ERROR: WORK not sized correctly
2302 call errmgr%report_error("solve_least_squares_mtx_cmplx", &
2303 "Incorrectly sized input array WORK, argument 3.", &
2304 la_array_size_error)
2305 return
2306 end if
2307 wptr => work(1:lwork)
2308 else
2309 allocate(wrk(lwork), stat = istat)
2310 if (istat /= 0) then
2311 ! ERROR: Out of memory
2312 call errmgr%report_error("solve_least_squares_mtx_cmplx", &
2313 "Insufficient memory available.", &
2314 la_out_of_memory_error)
2315 return
2316 end if
2317 wptr => wrk
2318 end if
2319
2320 ! Process
2321 call zgels('N', m, n, nrhs, a, m, b, maxmn, wptr, lwork, flag)
2322 if (flag > 0) then
2323 call errmgr%report_error("solve_least_squares_mtx_cmplx", &
2324 "The supplied matrix is not of full rank; therefore, " // &
2325 "the solution could not be computed via this routine. " // &
2326 "Try a routine that utilizes column pivoting.", &
2327 la_invalid_operation_error)
2328 end if
2329 end subroutine
2330
2331! ------------------------------------------------------------------------------
2332 module subroutine solve_least_squares_vec(a, b, work, olwork, err)
2333 ! Arguments
2334 real(real64), intent(inout), dimension(:,:) :: a
2335 real(real64), intent(inout), dimension(:) :: b
2336 real(real64), intent(out), target, optional, dimension(:) :: work
2337 integer(int32), intent(out), optional :: olwork
2338 class(errors), intent(inout), optional, target :: err
2339
2340 ! Local Variables
2341 integer(int32) :: m, n, maxmn, lwork, istat, flag
2342 real(real64), pointer, dimension(:) :: wptr
2343 real(real64), allocatable, target, dimension(:) :: wrk
2344 real(real64), dimension(1) :: temp
2345 class(errors), pointer :: errmgr
2346 type(errors), target :: deferr
2347
2348 ! Initialization
2349 m = size(a, 1)
2350 n = size(a, 2)
2351 maxmn = max(m, n)
2352 if (present(err)) then
2353 errmgr => err
2354 else
2355 errmgr => deferr
2356 end if
2357
2358 ! Input Check
2359 if (size(b) /= maxmn) then
2360 call errmgr%report_error("solve_least_squares_vec", &
2361 "Input 2 is not sized correctly.", la_array_size_error)
2362 return
2363 end if
2364
2365 ! Workspace Query
2366 call dgels('N', m, n, 1, a, m, b, maxmn, temp, -1, flag)
2367 lwork = int(temp(1), int32)
2368 if (present(olwork)) then
2369 olwork = lwork
2370 return
2371 end if
2372
2373 ! Local Memory Allocation
2374 if (present(work)) then
2375 if (size(work) < lwork) then
2376 ! ERROR: WORK not sized correctly
2377 call errmgr%report_error("solve_least_squares_vec", &
2378 "Incorrectly sized input array WORK, argument 3.", &
2379 la_array_size_error)
2380 return
2381 end if
2382 wptr => work(1:lwork)
2383 else
2384 allocate(wrk(lwork), stat = istat)
2385 if (istat /= 0) then
2386 ! ERROR: Out of memory
2387 call errmgr%report_error("solve_least_squares_vec", &
2388 "Insufficient memory available.", &
2389 la_out_of_memory_error)
2390 return
2391 end if
2392 wptr => wrk
2393 end if
2394
2395 ! Process
2396 call dgels('N', m, n, 1, a, m, b, maxmn, wptr, lwork, flag)
2397 if (flag > 0) then
2398 call errmgr%report_error("solve_least_squares_mtx", &
2399 "The supplied matrix is not of full rank; therefore, " // &
2400 "the solution could not be computed via this routine. " // &
2401 "Try a routine that utilizes column pivoting.", &
2402 la_invalid_operation_error)
2403 end if
2404 end subroutine
2405
2406! ------------------------------------------------------------------------------
2407 module subroutine solve_least_squares_vec_cmplx(a, b, work, olwork, err)
2408 ! Arguments
2409 complex(real64), intent(inout), dimension(:,:) :: a
2410 complex(real64), intent(inout), dimension(:) :: b
2411 complex(real64), intent(out), target, optional, dimension(:) :: work
2412 integer(int32), intent(out), optional :: olwork
2413 class(errors), intent(inout), optional, target :: err
2414
2415 ! Local Variables
2416 integer(int32) :: m, n, maxmn, lwork, istat, flag
2417 complex(real64), pointer, dimension(:) :: wptr
2418 complex(real64), allocatable, target, dimension(:) :: wrk
2419 complex(real64), dimension(1) :: temp
2420 class(errors), pointer :: errmgr
2421 type(errors), target :: deferr
2422
2423 ! Initialization
2424 m = size(a, 1)
2425 n = size(a, 2)
2426 maxmn = max(m, n)
2427 if (present(err)) then
2428 errmgr => err
2429 else
2430 errmgr => deferr
2431 end if
2432
2433 ! Input Check
2434 if (size(b) /= maxmn) then
2435 call errmgr%report_error("solve_least_squares_vec_cmplx", &
2436 "Input 2 is not sized correctly.", la_array_size_error)
2437 return
2438 end if
2439
2440 ! Workspace Query
2441 call zgels('N', m, n, 1, a, m, b, maxmn, temp, -1, flag)
2442 lwork = int(temp(1), int32)
2443 if (present(olwork)) then
2444 olwork = lwork
2445 return
2446 end if
2447
2448 ! Local Memory Allocation
2449 if (present(work)) then
2450 if (size(work) < lwork) then
2451 ! ERROR: WORK not sized correctly
2452 call errmgr%report_error("solve_least_squares_vec_cmplx", &
2453 "Incorrectly sized input array WORK, argument 3.", &
2454 la_array_size_error)
2455 return
2456 end if
2457 wptr => work(1:lwork)
2458 else
2459 allocate(wrk(lwork), stat = istat)
2460 if (istat /= 0) then
2461 ! ERROR: Out of memory
2462 call errmgr%report_error("solve_least_squares_vec_cmplx", &
2463 "Insufficient memory available.", &
2464 la_out_of_memory_error)
2465 return
2466 end if
2467 wptr => wrk
2468 end if
2469
2470 ! Process
2471 call zgels('N', m, n, 1, a, m, b, maxmn, wptr, lwork, flag)
2472 if (flag > 0) then
2473 call errmgr%report_error("solve_least_squares_mtx_cmplx", &
2474 "The supplied matrix is not of full rank; therefore, " // &
2475 "the solution could not be computed via this routine. " // &
2476 "Try a routine that utilizes column pivoting.", &
2477 la_invalid_operation_error)
2478 end if
2479 end subroutine
2480
2481! ------------------------------------------------------------------------------
2482 module subroutine solve_least_squares_mtx_pvt(a, b, ipvt, arnk, work, olwork, err)
2483 ! Arguments
2484 real(real64), intent(inout), dimension(:,:) :: a, b
2485 integer(int32), intent(inout), target, optional, dimension(:) :: ipvt
2486 integer(int32), intent(out), optional :: arnk
2487 real(real64), intent(out), target, optional, dimension(:) :: work
2488 integer(int32), intent(out), optional :: olwork
2489 class(errors), intent(inout), optional, target :: err
2490
2491 ! Local Variables
2492 integer(int32) :: m, n, maxmn, nrhs, lwork, istat, flag, rnk
2493 real(real64), pointer, dimension(:) :: wptr
2494 real(real64), allocatable, target, dimension(:) :: wrk
2495 integer(int32), allocatable, target, dimension(:) :: iwrk
2496 integer(int32), pointer, dimension(:) :: iptr
2497 real(real64), dimension(1) :: temp
2498 integer(int32), dimension(1) :: itemp
2499 real(real64) :: rc
2500 class(errors), pointer :: errmgr
2501 type(errors), target :: deferr
2502 character(len = :), allocatable :: errmsg
2503
2504 ! Initialization
2505 m = size(a, 1)
2506 n = size(a, 2)
2507 maxmn = max(m, n)
2508 nrhs = size(b, 2)
2509 rc = epsilon(rc)
2510 if (present(arnk)) arnk = 0
2511 if (present(err)) then
2512 errmgr => err
2513 else
2514 errmgr => deferr
2515 end if
2516
2517 ! Input Check
2518 flag = 0
2519 if (size(b, 1) /= maxmn) then
2520 flag = 2
2521 end if
2522 if (flag /= 0) then
2523 allocate(character(len = 256) :: errmsg)
2524 write(errmsg, 100) "Input number ", flag, &
2525 " is not sized correctly."
2526 call errmgr%report_error("solve_least_squares_mtx_pvt", &
2527 trim(errmsg), la_array_size_error)
2528 return
2529 end if
2530
2531 ! Workspace Query
2532 call dgelsy(m, n, nrhs, a, m, b, maxmn, itemp, rc, rnk, temp, -1, flag)
2533 lwork = int(temp(1), int32)
2534 if (present(olwork)) then
2535 olwork = lwork
2536 return
2537 end if
2538
2539 ! Local Memory Allocation
2540 if (present(ipvt)) then
2541 if (size(ipvt) < n) then
2542 ! ERROR: IPVT is not big enough
2543 call errmgr%report_error("solve_least_squares_mtx_pvt", &
2544 "Incorrectly sized pivot array, argument 3.", &
2545 la_array_size_error)
2546 return
2547 end if
2548 iptr => ipvt(1:n)
2549 else
2550 allocate(iwrk(n), stat = istat)
2551 if (istat /= 0) then
2552 ! ERROR: Out of memory
2553 call errmgr%report_error("solve_least_squares_mtx_pvt", &
2554 "Insufficient memory available.", &
2555 la_out_of_memory_error)
2556 return
2557 end if
2558 iptr => iwrk
2559 iptr = 0
2560 end if
2561
2562 if (present(work)) then
2563 if (size(work) < lwork) then
2564 ! ERROR: WORK not sized correctly
2565 call errmgr%report_error("solve_least_squares_mtx_pvt", &
2566 "Incorrectly sized input array WORK, argument 5.", &
2567 la_array_size_error)
2568 return
2569 end if
2570 wptr => work(1:lwork)
2571 else
2572 allocate(wrk(lwork), stat = istat)
2573 if (istat /= 0) then
2574 ! ERROR: Out of memory
2575 call errmgr%report_error("solve_least_squares_mtx_pvt", &
2576 "Insufficient memory available.", &
2577 la_out_of_memory_error)
2578 return
2579 end if
2580 wptr => wrk
2581 end if
2582
2583 ! Process
2584 call dgelsy(m, n, nrhs, a, m, b, maxmn, iptr, rc, rnk, wptr, lwork, &
2585 flag)
2586 if (present(arnk)) arnk = rnk
2587
2588 ! Formatting
2589100 format(a, i0, a)
2590 end subroutine
2591
2592! ------------------------------------------------------------------------------
2593 module subroutine solve_least_squares_mtx_pvt_cmplx(a, b, ipvt, arnk, &
2594 work, olwork, rwork, err)
2595 ! Arguments
2596 complex(real64), intent(inout), dimension(:,:) :: a, b
2597 integer(int32), intent(inout), target, optional, dimension(:) :: ipvt
2598 integer(int32), intent(out), optional :: arnk
2599 complex(real64), intent(out), target, optional, dimension(:) :: work
2600 integer(int32), intent(out), optional :: olwork
2601 real(real64), intent(out), target, optional, dimension(:) :: rwork
2602 class(errors), intent(inout), optional, target :: err
2603
2604 ! Local Variables
2605 integer(int32) :: m, n, maxmn, nrhs, lwork, istat, flag, rnk, lrwork
2606 complex(real64), pointer, dimension(:) :: wptr
2607 complex(real64), allocatable, target, dimension(:) :: wrk
2608 real(real64), pointer, dimension(:) :: rwptr
2609 real(real64), allocatable, target, dimension(:) :: rwrk
2610 integer(int32), allocatable, target, dimension(:) :: iwrk
2611 integer(int32), pointer, dimension(:) :: iptr
2612 complex(real64), dimension(1) :: temp
2613 real(real64), dimension(1) :: rtemp
2614 integer(int32), dimension(1) :: itemp
2615 real(real64) :: rc
2616 class(errors), pointer :: errmgr
2617 type(errors), target :: deferr
2618 character(len = :), allocatable :: errmsg
2619
2620 ! Initialization
2621 m = size(a, 1)
2622 n = size(a, 2)
2623 maxmn = max(m, n)
2624 nrhs = size(b, 2)
2625 lrwork = 2 * n
2626 rc = epsilon(rc)
2627 if (present(arnk)) arnk = 0
2628 if (present(err)) then
2629 errmgr => err
2630 else
2631 errmgr => deferr
2632 end if
2633
2634 ! Input Check
2635 flag = 0
2636 if (size(b, 1) /= maxmn) then
2637 flag = 2
2638 end if
2639 if (flag /= 0) then
2640 allocate(character(len = 256) :: errmsg)
2641 write(errmsg, 100) "Input number ", flag, &
2642 " is not sized correctly."
2643 call errmgr%report_error("solve_least_squares_mtx_pvt_cmplx", &
2644 trim(errmsg), la_array_size_error)
2645 return
2646 end if
2647
2648 ! Workspace Query
2649 call zgelsy(m, n, nrhs, a, m, b, maxmn, itemp, rc, rnk, temp, -1, &
2650 rtemp, flag)
2651 lwork = int(temp(1), int32)
2652 if (present(olwork)) then
2653 olwork = lwork
2654 return
2655 end if
2656
2657 ! Local Memory Allocation
2658 if (present(ipvt)) then
2659 if (size(ipvt) < n) then
2660 ! ERROR: IPVT is not big enough
2661 call errmgr%report_error("solve_least_squares_mtx_pvt", &
2662 "Incorrectly sized pivot array, argument 3.", &
2663 la_array_size_error)
2664 return
2665 end if
2666 iptr => ipvt(1:n)
2667 else
2668 allocate(iwrk(n), stat = istat)
2669 if (istat /= 0) then
2670 ! ERROR: Out of memory
2671 call errmgr%report_error("solve_least_squares_mtx_pvt", &
2672 "Insufficient memory available.", &
2673 la_out_of_memory_error)
2674 return
2675 end if
2676 iptr => iwrk
2677 iptr = 0
2678 end if
2679
2680 if (present(work)) then
2681 if (size(work) < lwork) then
2682 ! ERROR: WORK not sized correctly
2683 call errmgr%report_error("solve_least_squares_mtx_pvt", &
2684 "Incorrectly sized input array WORK, argument 5.", &
2685 la_array_size_error)
2686 return
2687 end if
2688 wptr => work(1:lwork)
2689 else
2690 allocate(wrk(lwork), stat = istat)
2691 if (istat /= 0) then
2692 ! ERROR: Out of memory
2693 call errmgr%report_error("solve_least_squares_mtx_pvt", &
2694 "Insufficient memory available.", &
2695 la_out_of_memory_error)
2696 return
2697 end if
2698 wptr => wrk
2699 end if
2700
2701 if (present(rwork)) then
2702 if (size(rwork) < lrwork) then
2703 ! ERROR: RWORK not sized correctly
2704 call errmgr%report_error("solve_least_squares_mtx_pvt_cmplx", &
2705 "Incorrectly sized input array RWORK, argument 7.", &
2706 la_array_size_error)
2707 return
2708 end if
2709 rwptr => rwork(1:lrwork)
2710 else
2711 allocate(rwrk(lrwork), stat = istat)
2712 if (istat /= 0) then
2713 ! ERROR: Out of memory
2714 call errmgr%report_error("solve_least_squares_mtx_pvt_cmplx", &
2715 "Insufficient memory available.", &
2716 la_out_of_memory_error)
2717 return
2718 end if
2719 rwptr => rwrk
2720 end if
2721
2722 ! Process
2723 call zgelsy(m, n, nrhs, a, m, b, maxmn, iptr, rc, rnk, wptr, lwork, &
2724 rwptr, flag)
2725 if (present(arnk)) arnk = rnk
2726
2727 ! Formatting
2728100 format(a, i0, a)
2729 end subroutine
2730
2731! ------------------------------------------------------------------------------
2732 module subroutine solve_least_squares_vec_pvt(a, b, ipvt, arnk, work, olwork, err)
2733 ! Arguments
2734 real(real64), intent(inout), dimension(:,:) :: a
2735 real(real64), intent(inout), dimension(:) :: b
2736 integer(int32), intent(inout), target, optional, dimension(:) :: ipvt
2737 integer(int32), intent(out), optional :: arnk
2738 real(real64), intent(out), target, optional, dimension(:) :: work
2739 integer(int32), intent(out), optional :: olwork
2740 class(errors), intent(inout), optional, target :: err
2741
2742 ! Local Variables
2743 integer(int32) :: m, n, maxmn, lwork, istat, flag, rnk
2744 real(real64), pointer, dimension(:) :: wptr
2745 real(real64), allocatable, target, dimension(:) :: wrk
2746 integer(int32), allocatable, target, dimension(:) :: iwrk
2747 integer(int32), pointer, dimension(:) :: iptr
2748 real(real64), dimension(1) :: temp
2749 integer(int32), dimension(1) :: itemp
2750 real(real64) :: rc
2751 class(errors), pointer :: errmgr
2752 type(errors), target :: deferr
2753 character(len = :), allocatable :: errmsg
2754
2755 ! Initialization
2756 m = size(a, 1)
2757 n = size(a, 2)
2758 maxmn = max(m, n)
2759 rc = epsilon(rc)
2760 if (present(arnk)) arnk = 0
2761 if (present(err)) then
2762 errmgr => err
2763 else
2764 errmgr => deferr
2765 end if
2766
2767 ! Input Check
2768 flag = 0
2769 if (size(b, 1) /= maxmn) then
2770 flag = 2
2771 end if
2772 if (flag /= 0) then
2773 allocate(character(len = 256) :: errmsg)
2774 write(errmsg, 100) "Input number ", flag, &
2775 " is not sized correctly."
2776 call errmgr%report_error("solve_least_squares_vec_pvt", &
2777 trim(errmsg), la_array_size_error)
2778 return
2779 end if
2780
2781 ! Workspace Query
2782 call dgelsy(m, n, 1, a, m, b, maxmn, itemp, rc, rnk, temp, -1, flag)
2783 lwork = int(temp(1), int32)
2784 if (present(olwork)) then
2785 olwork = lwork
2786 return
2787 end if
2788
2789 ! Local Memory Allocation
2790 if (present(ipvt)) then
2791 if (size(ipvt) < n) then
2792 ! ERROR: IPVT is not big enough
2793 call errmgr%report_error("solve_least_squares_mtx_pvt", &
2794 "Incorrectly sized pivot array, argument 3.", &
2795 la_array_size_error)
2796 return
2797 end if
2798 iptr => ipvt(1:n)
2799 else
2800 allocate(iwrk(n), stat = istat)
2801 if (istat /= 0) then
2802 ! ERROR: Out of memory
2803 call errmgr%report_error("solve_least_squares_mtx_pvt", &
2804 "Insufficient memory available.", &
2805 la_out_of_memory_error)
2806 return
2807 end if
2808 iptr => iwrk
2809 iptr = 0
2810 end if
2811
2812 if (present(work)) then
2813 if (size(work) < lwork) then
2814 ! ERROR: WORK not sized correctly
2815 call errmgr%report_error("solve_least_squares_vec_pvt", &
2816 "Incorrectly sized input array WORK, argument 5.", &
2817 la_array_size_error)
2818 return
2819 end if
2820 wptr => work(1:lwork)
2821 else
2822 allocate(wrk(lwork), stat = istat)
2823 if (istat /= 0) then
2824 ! ERROR: Out of memory
2825 call errmgr%report_error("solve_least_squares_vec_pvt", &
2826 "Insufficient memory available.", &
2827 la_out_of_memory_error)
2828 return
2829 end if
2830 wptr => wrk
2831 end if
2832
2833 ! Process
2834 call dgelsy(m, n, 1, a, m, b, maxmn, iptr, rc, rnk, wptr, lwork, flag)
2835 if (present(arnk)) arnk = rnk
2836
2837 ! Formatting
2838100 format(a, i0, a)
2839 end subroutine
2840
2841! ------------------------------------------------------------------------------
2842 module subroutine solve_least_squares_vec_pvt_cmplx(a, b, ipvt, arnk, &
2843 work, olwork, rwork, err)
2844 ! Arguments
2845 complex(real64), intent(inout), dimension(:,:) :: a
2846 complex(real64), intent(inout), dimension(:) :: b
2847 integer(int32), intent(inout), target, optional, dimension(:) :: ipvt
2848 integer(int32), intent(out), optional :: arnk
2849 complex(real64), intent(out), target, optional, dimension(:) :: work
2850 integer(int32), intent(out), optional :: olwork
2851 real(real64), intent(out), target, optional, dimension(:) :: rwork
2852 class(errors), intent(inout), optional, target :: err
2853
2854 ! Local Variables
2855 integer(int32) :: m, n, maxmn, lwork, istat, flag, rnk
2856 complex(real64), pointer, dimension(:) :: wptr
2857 complex(real64), allocatable, target, dimension(:) :: wrk
2858 real(real64), pointer, dimension(:) :: rwptr
2859 real(real64), allocatable, target, dimension(:) :: rwrk
2860 integer(int32), allocatable, target, dimension(:) :: iwrk
2861 integer(int32), pointer, dimension(:) :: iptr
2862 complex(real64), dimension(1) :: temp
2863 real(real64), dimension(1) :: rtemp
2864 integer(int32), dimension(1) :: itemp
2865 real(real64) :: rc
2866 class(errors), pointer :: errmgr
2867 type(errors), target :: deferr
2868 character(len = :), allocatable :: errmsg
2869
2870 ! Initialization
2871 m = size(a, 1)
2872 n = size(a, 2)
2873 maxmn = max(m, n)
2874 lrwork = 2 * n
2875 rc = epsilon(rc)
2876 if (present(arnk)) arnk = 0
2877 if (present(err)) then
2878 errmgr => err
2879 else
2880 errmgr => deferr
2881 end if
2882
2883 ! Input Check
2884 flag = 0
2885 if (size(b, 1) /= maxmn) then
2886 flag = 2
2887 end if
2888 if (flag /= 0) then
2889 allocate(character(len = 256) :: errmsg)
2890 write(errmsg, 100) "Input number ", flag, &
2891 " is not sized correctly."
2892 call errmgr%report_error("solve_least_squares_vec_pvt_cmplx", &
2893 trim(errmsg), la_array_size_error)
2894 return
2895 end if
2896
2897 ! Workspace Query
2898 call zgelsy(m, n, 1, a, m, b, maxmn, itemp, rc, rnk, temp, -1, rtemp, &
2899 flag)
2900 lwork = int(temp(1), int32)
2901 if (present(olwork)) then
2902 olwork = lwork
2903 return
2904 end if
2905
2906 ! Local Memory Allocation
2907 if (present(ipvt)) then
2908 if (size(ipvt) < n) then
2909 ! ERROR: IPVT is not big enough
2910 call errmgr%report_error("solve_least_squares_mtx_pvt", &
2911 "Incorrectly sized pivot array, argument 3.", &
2912 la_array_size_error)
2913 return
2914 end if
2915 iptr => ipvt(1:n)
2916 else
2917 allocate(iwrk(n), stat = istat)
2918 if (istat /= 0) then
2919 ! ERROR: Out of memory
2920 call errmgr%report_error("solve_least_squares_mtx_pvt", &
2921 "Insufficient memory available.", &
2922 la_out_of_memory_error)
2923 return
2924 end if
2925 iptr => iwrk
2926 iptr = 0
2927 end if
2928
2929 if (present(work)) then
2930 if (size(work) < lwork) then
2931 ! ERROR: WORK not sized correctly
2932 call errmgr%report_error("solve_least_squares_vec_pvt_cmplx", &
2933 "Incorrectly sized input array WORK, argument 5.", &
2934 la_array_size_error)
2935 return
2936 end if
2937 wptr => work(1:lwork)
2938 else
2939 allocate(wrk(lwork), stat = istat)
2940 if (istat /= 0) then
2941 ! ERROR: Out of memory
2942 call errmgr%report_error("solve_least_squares_vec_pvt_cmplx", &
2943 "Insufficient memory available.", &
2944 la_out_of_memory_error)
2945 return
2946 end if
2947 wptr => wrk
2948 end if
2949
2950 if (present(rwork)) then
2951 if (size(rwork) < lrwork) then
2952 ! ERROR: WORK not sized correctly
2953 call errmgr%report_error("solve_least_squares_vec_pvt_cmplx", &
2954 "Incorrectly sized input array RWORK, argument 7.", &
2955 la_array_size_error)
2956 return
2957 end if
2958 rwptr => rwork(1:lrwork)
2959 else
2960 allocate(rwrk(lrwork), stat = istat)
2961 if (istat /= 0) then
2962 ! ERROR: Out of memory
2963 call errmgr%report_error("solve_least_squares_vec_pvt_cmplx", &
2964 "Insufficient memory available.", &
2965 la_out_of_memory_error)
2966 return
2967 end if
2968 rwptr => rwrk
2969 end if
2970
2971 ! Process
2972 call zgelsy(m, n, 1, a, m, b, maxmn, iptr, rc, rnk, wptr, lwork, &
2973 rwptr, flag)
2974 if (present(arnk)) arnk = rnk
2975
2976 ! Formatting
2977100 format(a, i0, a)
2978 end subroutine
2979
2980! ------------------------------------------------------------------------------
2981 module subroutine solve_least_squares_mtx_svd(a, b, s, arnk, work, olwork, err)
2982 ! Arguments
2983 real(real64), intent(inout), dimension(:,:) :: a, b
2984 integer(int32), intent(out), optional :: arnk
2985 real(real64), intent(out), target, optional, dimension(:) :: work, s
2986 integer(int32), intent(out), optional :: olwork
2987 class(errors), intent(inout), optional, target :: err
2988
2989 ! Local Variables
2990 integer(int32) :: m, n, nrhs, mn, maxmn, istat, flag, lwork, rnk
2991 real(real64), pointer, dimension(:) :: wptr, sptr
2992 real(real64), allocatable, target, dimension(:) :: wrk, sing
2993 real(real64), dimension(1) :: temp
2994 real(real64) :: rcond
2995 class(errors), pointer :: errmgr
2996 type(errors), target :: deferr
2997 character(len = :), allocatable :: errmsg
2998
2999 ! Initialization
3000 m = size(a, 1)
3001 n = size(a, 2)
3002 nrhs = size(b, 2)
3003 mn = min(m, n)
3004 maxmn = max(m, n)
3005 rcond = epsilon(rcond)
3006 if (present(arnk)) arnk = 0
3007 if (present(err)) then
3008 errmgr => err
3009 else
3010 errmgr => deferr
3011 end if
3012
3013 ! Input Check
3014 flag = 0
3015 if (size(b, 1) /= maxmn) then
3016 flag = 2
3017 end if
3018 if (flag /= 0) then
3019 ! ERROR: One of the input arrays is not sized correctly
3020 allocate(character(len = 256) :: errmsg)
3021 write(errmsg, 100) "Input number ", flag, &
3022 " is not sized correctly."
3023 call errmgr%report_error("solve_least_squares_mtx_svd", &
3024 trim(errmsg), la_array_size_error)
3025 return
3026 end if
3027
3028 ! Workspace Query
3029 call dgelss(m, n, nrhs, a, m, b, maxmn, temp, rcond, rnk, temp, -1, &
3030 flag)
3031 lwork = int(temp(1), int32)
3032 if (present(olwork)) then
3033 olwork = lwork
3034 return
3035 end if
3036
3037 ! Local Memory Allocation
3038 if (present(s)) then
3039 if (size(s) < mn) then
3040 ! ERROR: S not sized correctly
3041 call errmgr%report_error("solve_least_squares_mtx_svd", &
3042 "Incorrectly sized input array S, argument 3.", &
3043 la_array_size_error)
3044 return
3045 end if
3046 sptr => s(1:mn)
3047 else
3048 allocate(sing(mn), stat = istat)
3049 if (istat /= 0) then
3050 ! ERROR: Out of memory
3051 call errmgr%report_error("solve_least_squares_mtx_svd", &
3052 "Insufficient memory available.", &
3053 la_out_of_memory_error)
3054 return
3055 end if
3056 sptr => sing
3057 end if
3058
3059 if (present(work)) then
3060 if (size(work) < lwork) then
3061 ! ERROR: WORK not sized correctly
3062 call errmgr%report_error("solve_least_squares_mtx_svd", &
3063 "Incorrectly sized input array WORK, argument 5.", &
3064 la_array_size_error)
3065 return
3066 end if
3067 wptr => work(1:lwork)
3068 else
3069 allocate(wrk(lwork), stat = istat)
3070 if (istat /= 0) then
3071 ! ERROR: Out of memory
3072 call errmgr%report_error("solve_least_squares_mtx_svd", &
3073 "Insufficient memory available.", &
3074 la_out_of_memory_error)
3075 return
3076 end if
3077 wptr => wrk
3078 end if
3079
3080 ! Process
3081 call dgelss(m, n, nrhs, a, m, b, maxmn, sptr, rcond, rnk, wptr, lwork, &
3082 flag)
3083 if (present(arnk)) arnk = rnk
3084 if (flag > 0) then
3085 allocate(character(len = 256) :: errmsg)
3086 write(errmsg, 101) flag, " superdiagonals could not " // &
3087 "converge to zero as part of the QR iteration process."
3088 call errmgr%report_warning("solve_least_squares_mtx_svd", errmsg, &
3089 la_convergence_error)
3090 end if
3091
3092 ! Formatting
3093100 format(a, i0, a)
3094101 format(i0, a)
3095 end subroutine
3096
3097! ------------------------------------------------------------------------------
3098 module subroutine solve_least_squares_mtx_svd_cmplx(a, b, s, arnk, work, &
3099 olwork, rwork, err)
3100 ! Arguments
3101 complex(real64), intent(inout), dimension(:,:) :: a, b
3102 integer(int32), intent(out), optional :: arnk
3103 complex(real64), intent(out), target, optional, dimension(:) :: work
3104 real(real64), intent(out), target, optional, dimension(:) :: s, rwork
3105 integer(int32), intent(out), optional :: olwork
3106 class(errors), intent(inout), optional, target :: err
3107
3108 ! Local Variables
3109 integer(int32) :: m, n, nrhs, mn, maxmn, istat, flag, lwork, rnk, lrwork
3110 complex(real64), pointer, dimension(:) :: wptr
3111 complex(real64), allocatable, target, dimension(:) :: wrk
3112 real(real64), pointer, dimension(:) :: rwptr, sptr
3113 real(real64), allocatable, target, dimension(:) :: rwrk, sing
3114 complex(real64), dimension(1) :: temp
3115 real(real64), dimension(1) :: rtemp
3116 real(real64) :: rcond
3117 class(errors), pointer :: errmgr
3118 type(errors), target :: deferr
3119 character(len = :), allocatable :: errmsg
3120
3121 ! Initialization
3122 m = size(a, 1)
3123 n = size(a, 2)
3124 nrhs = size(b, 2)
3125 mn = min(m, n)
3126 lrwork = 5 * mn
3127 maxmn = max(m, n)
3128 rcond = epsilon(rcond)
3129 if (present(arnk)) arnk = 0
3130 if (present(err)) then
3131 errmgr => err
3132 else
3133 errmgr => deferr
3134 end if
3135
3136 ! Input Check
3137 flag = 0
3138 if (size(b, 1) /= maxmn) then
3139 flag = 2
3140 end if
3141 if (flag /= 0) then
3142 ! ERROR: One of the input arrays is not sized correctly
3143 allocate(character(len = 256) :: errmsg)
3144 write(errmsg, 100) "Input number ", flag, &
3145 " is not sized correctly."
3146 call errmgr%report_error("solve_least_squares_mtx_svd_cmplx", &
3147 trim(errmsg), la_array_size_error)
3148 return
3149 end if
3150
3151 ! Workspace Query
3152 call zgelss(m, n, nrhs, a, m, b, maxmn, rtemp, rcond, rnk, temp, -1, &
3153 rtemp, flag)
3154 lwork = int(temp(1), int32)
3155 if (present(olwork)) then
3156 olwork = lwork
3157 return
3158 end if
3159
3160 ! Local Memory Allocation
3161 if (present(s)) then
3162 if (size(s) < mn) then
3163 ! ERROR: S not sized correctly
3164 call errmgr%report_error("solve_least_squares_mtx_svd_cmplx", &
3165 "Incorrectly sized input array S, argument 3.", &
3166 la_array_size_error)
3167 return
3168 end if
3169 sptr => s(1:mn)
3170 else
3171 allocate(sing(mn), stat = istat)
3172 if (istat /= 0) then
3173 ! ERROR: Out of memory
3174 call errmgr%report_error("solve_least_squares_mtx_svd_cmplx", &
3175 "Insufficient memory available.", &
3176 la_out_of_memory_error)
3177 return
3178 end if
3179 sptr => sing
3180 end if
3181
3182 if (present(work)) then
3183 if (size(work) < lwork) then
3184 ! ERROR: WORK not sized correctly
3185 call errmgr%report_error("solve_least_squares_mtx_svd_cmplx", &
3186 "Incorrectly sized input array WORK, argument 5.", &
3187 la_array_size_error)
3188 return
3189 end if
3190 wptr => work(1:lwork)
3191 else
3192 allocate(wrk(lwork), stat = istat)
3193 if (istat /= 0) then
3194 ! ERROR: Out of memory
3195 call errmgr%report_error("solve_least_squares_mtx_svd_cmplx", &
3196 "Insufficient memory available.", &
3197 la_out_of_memory_error)
3198 return
3199 end if
3200 wptr => wrk
3201 end if
3202
3203 if (present(rwork)) then
3204 if (size(rwork) < lrwork) then
3205 ! ERROR: WORK not sized correctly
3206 call errmgr%report_error("solve_least_squares_mtx_svd_cmplx", &
3207 "Incorrectly sized input array RWORK, argument 7.", &
3208 la_array_size_error)
3209 return
3210 end if
3211 rwptr => rwork(1:lrwork)
3212 else
3213 allocate(rwrk(lrwork), stat = istat)
3214 if (istat /= 0) then
3215 ! ERROR: Out of memory
3216 call errmgr%report_error("solve_least_squares_mtx_svd_cmplx", &
3217 "Insufficient memory available.", &
3218 la_out_of_memory_error)
3219 return
3220 end if
3221 rwptr => rwrk
3222 end if
3223
3224 ! Process
3225 call zgelss(m, n, nrhs, a, m, b, maxmn, sptr, rcond, rnk, wptr, lwork, &
3226 rwptr, flag)
3227 if (present(arnk)) arnk = rnk
3228 if (flag > 0) then
3229 allocate(character(len = 256) :: errmsg)
3230 write(errmsg, 101) flag, " superdiagonals could not " // &
3231 "converge to zero as part of the QR iteration process."
3232 call errmgr%report_warning("solve_least_squares_mtx_svd_cmplx", &
3233 errmsg, la_convergence_error)
3234 end if
3235
3236 ! Formatting
3237100 format(a, i0, a)
3238101 format(i0, a)
3239 end subroutine
3240
3241! ------------------------------------------------------------------------------
3242 module subroutine solve_least_squares_vec_svd(a, b, s, arnk, work, olwork, err)
3243 ! Arguments
3244 real(real64), intent(inout), dimension(:,:) :: a
3245 real(real64), intent(inout), dimension(:) :: b
3246 integer(int32), intent(out), optional :: arnk
3247 real(real64), intent(out), target, optional, dimension(:) :: work, s
3248 integer(int32), intent(out), optional :: olwork
3249 class(errors), intent(inout), optional, target :: err
3250
3251 ! Local Variables
3252 integer(int32) :: m, n, mn, maxmn, istat, flag, lwork, rnk
3253 real(real64), pointer, dimension(:) :: wptr, sptr
3254 real(real64), allocatable, target, dimension(:) :: wrk, sing
3255 real(real64), dimension(1) :: temp
3256 real(real64) :: rcond
3257 class(errors), pointer :: errmgr
3258 type(errors), target :: deferr
3259 character(len = :), allocatable :: errmsg
3260
3261 ! Initialization
3262 m = size(a, 1)
3263 n = size(a, 2)
3264 mn = min(m, n)
3265 maxmn = max(m, n)
3266 rcond = epsilon(rcond)
3267 if (present(arnk)) arnk = 0
3268 if (present(err)) then
3269 errmgr => err
3270 else
3271 errmgr => deferr
3272 end if
3273
3274 ! Input Check
3275 flag = 0
3276 if (size(b) /= maxmn) then
3277 flag = 2
3278 end if
3279 if (flag /= 0) then
3280 ! ERROR: One of the input arrays is not sized correctly
3281 allocate(character(len = 256) :: errmsg)
3282 write(errmsg, 100) "Input number ", flag, &
3283 " is not sized correctly."
3284 call errmgr%report_error("solve_least_squares_vec_svd", &
3285 trim(errmsg), la_array_size_error)
3286 return
3287 end if
3288
3289 ! Workspace Query
3290 call dgelss(m, n, 1, a, m, b, maxmn, temp, rcond, rnk, temp, -1, flag)
3291 lwork = int(temp(1), int32)
3292 if (present(olwork)) then
3293 olwork = lwork
3294 return
3295 end if
3296
3297 ! Local Memory Allocation
3298 if (present(s)) then
3299 if (size(s) < mn) then
3300 ! ERROR: S not sized correctly
3301 call errmgr%report_error("solve_least_squares_vec_svd", &
3302 "Incorrectly sized input array S, argument 3.", &
3303 la_array_size_error)
3304 return
3305 end if
3306 sptr => s(1:mn)
3307 else
3308 allocate(sing(mn), stat = istat)
3309 if (istat /= 0) then
3310 ! ERROR: Out of memory
3311 call errmgr%report_error("solve_least_squares_vec_svd", &
3312 "Insufficient memory available.", &
3313 la_out_of_memory_error)
3314 return
3315 end if
3316 sptr => sing
3317 end if
3318
3319 if (present(work)) then
3320 if (size(work) < lwork) then
3321 ! ERROR: WORK not sized correctly
3322 call errmgr%report_error("solve_least_squares_vec_svd", &
3323 "Incorrectly sized input array WORK, argument 5.", &
3324 la_array_size_error)
3325 return
3326 end if
3327 wptr => work(1:lwork)
3328 else
3329 allocate(wrk(lwork), stat = istat)
3330 if (istat /= 0) then
3331 ! ERROR: Out of memory
3332 call errmgr%report_error("solve_least_squares_vec_svd", &
3333 "Insufficient memory available.", &
3334 la_out_of_memory_error)
3335 return
3336 end if
3337 wptr => wrk
3338 end if
3339
3340 ! Process
3341 call dgelss(m, n, 1, a, m, b, maxmn, sptr, rcond, rnk, wptr, lwork, &
3342 flag)
3343 if (present(arnk)) arnk = rnk
3344 if (flag > 0) then
3345 allocate(character(len = 256) :: errmsg)
3346 write(errmsg, 101) flag, " superdiagonals could not " // &
3347 "converge to zero as part of the QR iteration process."
3348 call errmgr%report_warning("solve_least_squares_vec_svd", errmsg, &
3349 la_convergence_error)
3350 end if
3351
3352 ! Formatting
3353100 format(a, i0, a)
3354101 format(i0, a)
3355 end subroutine
3356
3357! ------------------------------------------------------------------------------
3358 module subroutine solve_least_squares_vec_svd_cmplx(a, b, s, arnk, work, &
3359 olwork, rwork, err)
3360 ! Arguments
3361 complex(real64), intent(inout), dimension(:,:) :: a
3362 complex(real64), intent(inout), dimension(:) :: b
3363 integer(int32), intent(out), optional :: arnk
3364 complex(real64), intent(out), target, optional, dimension(:) :: work
3365 real(real64), intent(out), target, optional, dimension(:) :: rwork, s
3366 integer(int32), intent(out), optional :: olwork
3367 class(errors), intent(inout), optional, target :: err
3368
3369 ! Local Variables
3370 integer(int32) :: m, n, mn, maxmn, istat, flag, lwork, rnk, lrwork
3371 real(real64), pointer, dimension(:) :: rwptr, sptr
3372 real(real64), allocatable, target, dimension(:) :: rwrk, sing
3373 complex(real64), pointer, dimension(:) :: wptr
3374 complex(real64), allocatable, target, dimension(:) :: wrk
3375 complex(real64), dimension(1) :: temp
3376 real(real64), dimension(1) :: rtemp
3377 real(real64) :: rcond
3378 class(errors), pointer :: errmgr
3379 type(errors), target :: deferr
3380 character(len = :), allocatable :: errmsg
3381
3382 ! Initialization
3383 m = size(a, 1)
3384 n = size(a, 2)
3385 mn = min(m, n)
3386 lrwork = 5 * mn
3387 maxmn = max(m, n)
3388 rcond = epsilon(rcond)
3389 if (present(arnk)) arnk = 0
3390 if (present(err)) then
3391 errmgr => err
3392 else
3393 errmgr => deferr
3394 end if
3395
3396 ! Input Check
3397 flag = 0
3398 if (size(b) /= maxmn) then
3399 flag = 2
3400 end if
3401 if (flag /= 0) then
3402 ! ERROR: One of the input arrays is not sized correctly
3403 allocate(character(len = 256) :: errmsg)
3404 write(errmsg, 100) "Input number ", flag, &
3405 " is not sized correctly."
3406 call errmgr%report_error("solve_least_squares_vec_svd_cmplx", &
3407 trim(errmsg), la_array_size_error)
3408 return
3409 end if
3410
3411 ! Workspace Query
3412 call zgelss(m, n, 1, a, m, b, maxmn, rtemp, rcond, rnk, temp, -1, &
3413 rtemp, flag)
3414 lwork = int(temp(1), int32)
3415 if (present(olwork)) then
3416 olwork = lwork
3417 return
3418 end if
3419
3420 ! Local Memory Allocation
3421 if (present(s)) then
3422 if (size(s) < mn) then
3423 ! ERROR: S not sized correctly
3424 call errmgr%report_error("solve_least_squares_vec_svd_cmplx", &
3425 "Incorrectly sized input array S, argument 3.", &
3426 la_array_size_error)
3427 return
3428 end if
3429 sptr => s(1:mn)
3430 else
3431 allocate(sing(mn), stat = istat)
3432 if (istat /= 0) then
3433 ! ERROR: Out of memory
3434 call errmgr%report_error("solve_least_squares_vec_svd_cmplx", &
3435 "Insufficient memory available.", &
3436 la_out_of_memory_error)
3437 return
3438 end if
3439 sptr => sing
3440 end if
3441
3442 if (present(work)) then
3443 if (size(work) < lwork) then
3444 ! ERROR: WORK not sized correctly
3445 call errmgr%report_error("solve_least_squares_vec_svd_cmplx", &
3446 "Incorrectly sized input array WORK, argument 5.", &
3447 la_array_size_error)
3448 return
3449 end if
3450 wptr => work(1:lwork)
3451 else
3452 allocate(wrk(lwork), stat = istat)
3453 if (istat /= 0) then
3454 ! ERROR: Out of memory
3455 call errmgr%report_error("solve_least_squares_vec_svd_cmplx", &
3456 "Insufficient memory available.", &
3457 la_out_of_memory_error)
3458 return
3459 end if
3460 wptr => wrk
3461 end if
3462
3463 if (present(rwork)) then
3464 if (size(rwork) < lrwork) then
3465 ! ERROR: WORK not sized correctly
3466 call errmgr%report_error("solve_least_squares_vec_svd_cmplx", &
3467 "Incorrectly sized input array RWORK, argument 7.", &
3468 la_array_size_error)
3469 return
3470 end if
3471 rwptr => rwork(1:lrwork)
3472 else
3473 allocate(rwrk(lrwork), stat = istat)
3474 if (istat /= 0) then
3475 ! ERROR: Out of memory
3476 call errmgr%report_error("solve_least_squares_vec_svd_cmplx", &
3477 "Insufficient memory available.", &
3478 la_out_of_memory_error)
3479 return
3480 end if
3481 rwptr => rwrk
3482 end if
3483
3484 ! Process
3485 call zgelss(m, n, 1, a, m, b, maxmn, sptr, rcond, rnk, wptr, lwork, &
3486 rwptr, flag)
3487 if (present(arnk)) arnk = rnk
3488 if (flag > 0) then
3489 allocate(character(len = 256) :: errmsg)
3490 write(errmsg, 101) flag, " superdiagonals could not " // &
3491 "converge to zero as part of the QR iteration process."
3492 call errmgr%report_warning("solve_least_squares_vec_svd_cmplx", &
3493 errmsg, la_convergence_error)
3494 end if
3495
3496 ! Formatting
3497100 format(a, i0, a)
3498101 format(i0, a)
3499 end subroutine
3500
3501! ******************************************************************************
3502! LQ SOLUTION
3503! ------------------------------------------------------------------------------
3504 module subroutine solve_lq_mtx(a, tau, b, work, olwork, err)
3505 ! Arguments
3506 real(real64), intent(in), dimension(:,:) :: a
3507 real(real64), intent(in), dimension(:) :: tau
3508 real(real64), intent(inout), dimension(:,:) :: b
3509 real(real64), intent(out), target, optional, dimension(:) :: work
3510 integer(int32), intent(out), optional :: olwork
3511 class(errors), intent(inout), optional, target :: err
3512
3513 ! Parameters
3514 real(real64), parameter :: one = 1.0d0
3515
3516 ! Local Variables
3517 integer(int32) :: m, n, nrhs, k, lwork, flag, istat
3518 real(real64), pointer, dimension(:) :: wptr
3519 real(real64), allocatable, target, dimension(:) :: wrk
3520 class(errors), pointer :: errmgr
3521 type(errors), target :: deferr
3522 character(len = :), allocatable :: errmsg
3523
3524 ! Initialization
3525 m = size(a, 1)
3526 n = size(a, 2)
3527 nrhs = size(b, 2)
3528 k = min(m, n)
3529 if (present(err)) then
3530 errmgr => err
3531 else
3532 errmgr => deferr
3533 end if
3534
3535 ! Input Check
3536 flag = 0
3537 if (m > n) then
3538 flag = 1
3539 else if (size(tau) /= k) then
3540 flag = 2
3541 else if (size(b, 1) /= n) then
3542 flag = 3
3543 end if
3544
3545 if (flag /= 0) then
3546 ! ERROR: One of the input arrays is not sized correctly
3547 allocate(character(len = 256) :: errmsg)
3548 write(errmsg, 100) "Input number ", flag, &
3549 " is not sized correctly."
3550 call errmgr%report_error("solve_lq_mtx", trim(errmsg), &
3551 la_array_size_error)
3552 return
3553 end if
3554
3555 ! Workspace Query
3556 call mult_lq(.true., .true., a, tau, b, olwork = lwork)
3557
3558 if (present(olwork)) then
3559 olwork = lwork
3560 return
3561 end if
3562
3563 ! Local Memory Allocation
3564 if (present(work)) then
3565 if (size(work) < lwork) then
3566 ! ERROR: WORK not sized correctly
3567 call errmgr%report_error("solve_lq_mtx", &
3568 "Incorrectly sized input array WORK, argument 4.", &
3569 la_array_size_error)
3570 return
3571 end if
3572 wptr => work(1:lwork)
3573 else
3574 allocate(wrk(lwork), stat = istat)
3575 if (istat /= 0) then
3576 ! ERROR: Out of memory
3577 call errmgr%report_error("solve_lq_mtx", &
3578 "Insufficient memory available.", &
3579 la_out_of_memory_error)
3580 return
3581 end if
3582 wptr => wrk
3583 end if
3584
3585 ! Solve the lower triangular system L * Y = B for Y, where Y = Q * X.
3586 ! The lower triangular system is M-by-M and Y is M-by-NHRS.
3587 call solve_triangular_system(.true., .false., .false., .true., one, &
3588 a(1:m,1:m), b(1:m,:), errmgr)
3589 if (errmgr%has_error_occurred()) return
3590
3591 ! Compute Q**T * Y = X
3592 call mult_lq(.true., .true., a, tau, b, work = wptr, err = errmgr)
3593 if (errmgr%has_error_occurred()) return
3594
3595 ! Formatting
3596100 format(a, i0, a)
3597 end subroutine
3598
3599! ------------------------------------------------------------------------------
3600 module subroutine solve_lq_mtx_cmplx(a, tau, b, work, olwork, err)
3601 ! Arguments
3602 complex(real64), intent(in), dimension(:,:) :: a
3603 complex(real64), intent(in), dimension(:) :: tau
3604 complex(real64), intent(inout), dimension(:,:) :: b
3605 complex(real64), intent(out), target, optional, dimension(:) :: work
3606 integer(int32), intent(out), optional :: olwork
3607 class(errors), intent(inout), optional, target :: err
3608
3609 ! Parameters
3610 complex(real64), parameter :: one = (1.0d0, 0.0d0)
3611
3612 ! Local Variables
3613 integer(int32) :: m, n, nrhs, k, lwork, flag, istat
3614 complex(real64), pointer, dimension(:) :: wptr
3615 complex(real64), allocatable, target, dimension(:) :: wrk
3616 class(errors), pointer :: errmgr
3617 type(errors), target :: deferr
3618 character(len = :), allocatable :: errmsg
3619
3620 ! Initialization
3621 m = size(a, 1)
3622 n = size(a, 2)
3623 nrhs = size(b, 2)
3624 k = min(m, n)
3625 if (present(err)) then
3626 errmgr => err
3627 else
3628 errmgr => deferr
3629 end if
3630
3631 ! Input Check
3632 flag = 0
3633 if (m > n) then
3634 flag = 1
3635 else if (size(tau) /= k) then
3636 flag = 2
3637 else if (size(b, 1) /= n) then
3638 flag = 3
3639 end if
3640
3641 if (flag /= 0) then
3642 ! ERROR: One of the input arrays is not sized correctly
3643 allocate(character(len = 256) :: errmsg)
3644 write(errmsg, 100) "Input number ", flag, &
3645 " is not sized correctly."
3646 call errmgr%report_error("solve_lq_mtx_cmplx", trim(errmsg), &
3647 la_array_size_error)
3648 return
3649 end if
3650
3651 ! Workspace Query
3652 call mult_lq(.true., .true., a, tau, b, olwork = lwork)
3653
3654 if (present(olwork)) then
3655 olwork = lwork
3656 return
3657 end if
3658
3659 ! Local Memory Allocation
3660 if (present(work)) then
3661 if (size(work) < lwork) then
3662 ! ERROR: WORK not sized correctly
3663 call errmgr%report_error("solve_lq_mtx_cmplx", &
3664 "Incorrectly sized input array WORK, argument 4.", &
3665 la_array_size_error)
3666 return
3667 end if
3668 wptr => work(1:lwork)
3669 else
3670 allocate(wrk(lwork), stat = istat)
3671 if (istat /= 0) then
3672 ! ERROR: Out of memory
3673 call errmgr%report_error("solve_lq_mtx_cmplx", &
3674 "Insufficient memory available.", &
3675 la_out_of_memory_error)
3676 return
3677 end if
3678 wptr => wrk
3679 end if
3680
3681 ! Solve the lower triangular system L * Y = B for Y, where Y = Q * X.
3682 ! The lower triangular system is M-by-M and Y is M-by-NHRS.
3683 call solve_triangular_system(.true., .false., .false., .true., one, &
3684 a(1:m,1:m), b(1:m,:), errmgr)
3685 if (errmgr%has_error_occurred()) return
3686
3687 ! Compute Q**T * Y = X
3688 call mult_lq(.true., .true., a, tau, b, work = wptr, err = errmgr)
3689 if (errmgr%has_error_occurred()) return
3690
3691 ! Formatting
3692100 format(a, i0, a)
3693 end subroutine
3694
3695! ------------------------------------------------------------------------------
3696 module subroutine solve_lq_vec(a, tau, b, work, olwork, err)
3697 ! Arguments
3698 real(real64), intent(in), dimension(:,:) :: a
3699 real(real64), intent(in), dimension(:) :: tau
3700 real(real64), intent(inout), dimension(:) :: b
3701 real(real64), intent(out), target, optional, dimension(:) :: work
3702 integer(int32), intent(out), optional :: olwork
3703 class(errors), intent(inout), optional, target :: err
3704
3705 ! Local Variables
3706 integer(int32) :: m, n, k, lwork, flag, istat
3707 real(real64), pointer, dimension(:) :: wptr
3708 real(real64), allocatable, target, dimension(:) :: wrk
3709 class(errors), pointer :: errmgr
3710 type(errors), target :: deferr
3711 character(len = :), allocatable :: errmsg
3712
3713 ! Initialization
3714 m = size(a, 1)
3715 n = size(a, 2)
3716 k = min(m, n)
3717 if (present(err)) then
3718 errmgr => err
3719 else
3720 errmgr => deferr
3721 end if
3722
3723 ! Input Check
3724 flag = 0
3725 if (m > n) then
3726 flag = 1
3727 else if (size(tau) /= k) then
3728 flag = 2
3729 else if (size(b) /= n) then
3730 flag = 3
3731 end if
3732
3733 if (flag /= 0) then
3734 ! ERROR: One of the input arrays is not sized correctly
3735 allocate(character(len = 256) :: errmsg)
3736 write(errmsg, 100) "Input number ", flag, &
3737 " is not sized correctly."
3738 call errmgr%report_error("solve_lq_vec", trim(errmsg), &
3739 la_array_size_error)
3740 return
3741 end if
3742
3743 ! Workspace Query
3744 call mult_lq(.true., a, tau, b, olwork = lwork)
3745
3746 if (present(olwork)) then
3747 olwork = lwork
3748 return
3749 end if
3750
3751 ! Local Memory Allocation
3752 if (present(work)) then
3753 if (size(work) < lwork) then
3754 ! ERROR: WORK not sized correctly
3755 call errmgr%report_error("solve_lq_vec", &
3756 "Incorrectly sized input array WORK, argument 4.", &
3757 la_array_size_error)
3758 return
3759 end if
3760 wptr => work(1:lwork)
3761 else
3762 allocate(wrk(lwork), stat = istat)
3763 if (istat /= 0) then
3764 ! ERROR: Out of memory
3765 call errmgr%report_error("solve_lq_vec", &
3766 "Insufficient memory available.", &
3767 la_out_of_memory_error)
3768 return
3769 end if
3770 wptr => wrk
3771 end if
3772
3773 ! Solve the lower triangular system L * Y = B for Y, where Y = Q * X.
3774 ! The lower triangular system is M-by-M and Y is M-by-NHRS.
3775 call solve_triangular_system(.false., .false., .true., a(1:m,1:m), &
3776 b(1:m), errmgr)
3777 if (errmgr%has_error_occurred()) return
3778
3779 ! Compute Q**T * Y = X
3780 call mult_lq(.true., a, tau, b, work = wptr, err = errmgr)
3781 if (errmgr%has_error_occurred()) return
3782
3783 ! Formatting
3784100 format(a, i0, a)
3785 end subroutine
3786
3787! ------------------------------------------------------------------------------
3788 module subroutine solve_lq_vec_cmplx(a, tau, b, work, olwork, err)
3789 ! Arguments
3790 complex(real64), intent(in), dimension(:,:) :: a
3791 complex(real64), intent(in), dimension(:) :: tau
3792 complex(real64), intent(inout), dimension(:) :: b
3793 complex(real64), intent(out), target, optional, dimension(:) :: work
3794 integer(int32), intent(out), optional :: olwork
3795 class(errors), intent(inout), optional, target :: err
3796
3797 ! Local Variables
3798 integer(int32) :: m, n, k, lwork, flag, istat
3799 complex(real64), pointer, dimension(:) :: wptr
3800 complex(real64), allocatable, target, dimension(:) :: wrk
3801 class(errors), pointer :: errmgr
3802 type(errors), target :: deferr
3803 character(len = :), allocatable :: errmsg
3804
3805 ! Initialization
3806 m = size(a, 1)
3807 n = size(a, 2)
3808 k = min(m, n)
3809 if (present(err)) then
3810 errmgr => err
3811 else
3812 errmgr => deferr
3813 end if
3814
3815 ! Input Check
3816 flag = 0
3817 if (m > n) then
3818 flag = 1
3819 else if (size(tau) /= k) then
3820 flag = 2
3821 else if (size(b) /= n) then
3822 flag = 3
3823 end if
3824
3825 if (flag /= 0) then
3826 ! ERROR: One of the input arrays is not sized correctly
3827 allocate(character(len = 256) :: errmsg)
3828 write(errmsg, 100) "Input number ", flag, &
3829 " is not sized correctly."
3830 call errmgr%report_error("solve_lq_vec_cmplx", trim(errmsg), &
3831 la_array_size_error)
3832 return
3833 end if
3834
3835 ! Workspace Query
3836 call mult_lq(.true., a, tau, b, olwork = lwork)
3837
3838 if (present(olwork)) then
3839 olwork = lwork
3840 return
3841 end if
3842
3843 ! Local Memory Allocation
3844 if (present(work)) then
3845 if (size(work) < lwork) then
3846 ! ERROR: WORK not sized correctly
3847 call errmgr%report_error("solve_lq_vec_cmplx", &
3848 "Incorrectly sized input array WORK, argument 4.", &
3849 la_array_size_error)
3850 return
3851 end if
3852 wptr => work(1:lwork)
3853 else
3854 allocate(wrk(lwork), stat = istat)
3855 if (istat /= 0) then
3856 ! ERROR: Out of memory
3857 call errmgr%report_error("solve_lq_vec_cmplx", &
3858 "Insufficient memory available.", &
3859 la_out_of_memory_error)
3860 return
3861 end if
3862 wptr => wrk
3863 end if
3864
3865 ! Solve the lower triangular system L * Y = B for Y, where Y = Q * X.
3866 ! The lower triangular system is M-by-M and Y is M-by-NHRS.
3867 call solve_triangular_system(.false., .false., .true., a(1:m,1:m), &
3868 b(1:m), errmgr)
3869 if (errmgr%has_error_occurred()) return
3870
3871 ! Compute Q**T * Y = X
3872 call mult_lq(.true., a, tau, b, work = wptr, err = errmgr)
3873 if (errmgr%has_error_occurred()) return
3874
3875 ! Formatting
3876100 format(a, i0, a)
3877 end subroutine
3878
3879! ------------------------------------------------------------------------------
3880end submodule
Provides a set of common linear algebra routines.
Definition: linalg.f90:145