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