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