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