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_basic.f90
1! linalg_basic.f90
2
3submodule(linalg) linalg_basic
4contains
5! ******************************************************************************
6! MATRIX MULTIPLICATION ROUTINES
7! ------------------------------------------------------------------------------
8 module subroutine mtx_mult_mtx(transa, transb, alpha, a, b, beta, c, err)
9 ! Arguments
10 logical, intent(in) :: transa, transb
11 real(real64), intent(in) :: alpha, beta
12 real(real64), intent(in), dimension(:,:) :: a, b
13 real(real64), intent(inout), dimension(:,:) :: c
14 class(errors), intent(inout), optional, target :: err
15
16 ! Parameters
17 real(real64), parameter :: zero = 0.0d0
18 real(real64), parameter :: one = 1.0d0
19
20 ! Local Variables
21 character :: ta, tb
22 integer(int32) :: m, n, k, lda, ldb, flag
23 class(errors), pointer :: errmgr
24 type(errors), target :: deferr
25 character(len = 128) :: errmsg
26
27 ! Initialization
28 m = size(c, 1)
29 n = size(c, 2)
30 if (transa) then ! K = # of columns in op(A) (# of rows in op(B))
31 k = size(a, 1)
32 ta = 'T'
33 lda = k
34 else
35 k = size(a, 2)
36 ta = 'N'
37 lda = m
38 end if
39 if (transb) then
40 tb = 'T'
41 ldb = n
42 else
43 tb = 'N'
44 ldb = k
45 end if
46 if (present(err)) then
47 errmgr => err
48 else
49 errmgr => deferr
50 end if
51
52 ! Input Check
53 flag = 0
54 if (transa) then
55 if (size(a, 2) /= m) flag = 4
56 else
57 if (size(a, 1) /= m) flag = 4
58 end if
59 if (transb) then
60 if (size(b, 2) /= k .or. size(b, 1) /= n) flag = 5
61 else
62 if (size(b, 1) /= k .or. size(b, 2) /= n) flag = 5
63 end if
64 if (flag /= 0) then
65 ! ERROR: Matrix dimensions mismatch
66 write(errmsg, 100) &
67 "Matrix dimension mismatch. Input number ", flag, &
68 " was not sized correctly."
69 call errmgr%report_error("mtx_mult_mtx", errmsg, &
70 la_array_size_error)
71 return
72 end if
73
74 ! Call DGEMM
75 call dgemm(ta, tb, m, n, k, alpha, a, lda, b, ldb, beta, c, m)
76
77 ! Formatting
78100 format(a, i0, a)
79 end subroutine
80
81! ------------------------------------------------------------------------------
82 module subroutine mtx_mult_vec(trans, alpha, a, b, beta, c, err)
83 ! Arguments
84 logical, intent(in) :: trans
85 real(real64), intent(in) :: alpha, beta
86 real(real64), intent(in), dimension(:,:) :: a
87 real(real64), intent(in), dimension(:) :: b
88 real(real64), intent(inout), dimension(:) :: c
89 class(errors), intent(inout), optional, target :: err
90
91 ! Local Variables
92 character :: t
93 integer(int32) :: m, n, flag
94 class(errors), pointer :: errmgr
95 type(errors), target :: deferr
96 character(len = 128) :: errmsg
97
98 ! Initialization
99 m = size(a, 1)
100 n = size(a, 2)
101 t = 'N'
102 if (trans) t = 'T'
103 if (present(err)) then
104 errmgr => err
105 else
106 errmgr => deferr
107 end if
108
109 ! Input Check
110 flag = 0
111 if (trans) then
112 if (size(b) /= m) then
113 flag = 4
114 else if (size(c) /= n) then
115 flag = 6
116 end if
117 else
118 if (size(b) /= n) then
119 flag = 4
120 else if (size(c) /= m) then
121 flag = 6
122 end if
123 end if
124 if (flag /= 0) then
125 ! ERROR: Matrix dimensions mismatch
126 write(errmsg, 100) &
127 "Matrix dimension mismatch. Input number ", flag, &
128 " was not sized correctly."
129 call errmgr%report_error("mtx_mult_vec", errmsg, &
130 la_array_size_error)
131 return
132 end if
133
134 ! Call DGEMV
135 call dgemv(t, m, n, alpha, a, m, b, 1, beta, c, 1)
136
137 ! Formatting
138100 format(a, i0, a)
139 end subroutine
140
141! xxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxx !
142! COMPLEX VALUED VERSIONS !
143! xxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxx !
144 module subroutine cmtx_mult_mtx(opa, opb, alpha, a, b, beta, c, err)
145 ! Arguments
146 integer(int32), intent(in) :: opa, opb
147 complex(real64), intent(in) :: alpha, beta
148 complex(real64), intent(in), dimension(:,:) :: a, b
149 complex(real64), intent(inout), dimension(:,:) :: c
150 class(errors), intent(inout), optional, target :: err
151
152 ! Parameters
153 real(real64), parameter :: zero = 0.0d0
154 real(real64), parameter :: one = 1.0d0
155
156 ! Local Variables
157 character :: ta, tb
158 integer(int32) :: m, n, k, lda, ldb, flag
159 class(errors), pointer :: errmgr
160 type(errors), target :: deferr
161 character(len = 128) :: errmsg
162
163 ! Initialization
164 m = size(c, 1)
165 n = size(c, 2)
166 if (opa == la_transpose) then ! K = # of columns in op(A) (# of rows in op(B))
167 k = size(a, 1)
168 ta = 'T'
169 lda = k
170 else if (opa == la_hermitian_transpose) then
171 k = size(a, 1)
172 ta = 'H'
173 lda = k
174 else
175 k = size(a, 2)
176 ta = 'N'
177 lda = m
178 end if
179 if (opb == la_transpose) then
180 tb = 'T'
181 ldb = n
182 else if (opb == la_hermitian_transpose) then
183 tb = 'H'
184 ldb = n
185 else
186 tb = 'N'
187 ldb = k
188 end if
189 if (present(err)) then
190 errmgr => err
191 else
192 errmgr => deferr
193 end if
194
195 ! Input Check
196 flag = 0
197 if (opa == la_transpose .or. opa == la_hermitian_transpose) then
198 if (size(a, 2) /= m) flag = 4
199 else
200 if (size(a, 1) /= m) flag = 4
201 end if
202 if (opb == la_transpose .or. opb == la_hermitian_transpose) then
203 if (size(b, 2) /= k .or. size(b, 1) /= n) flag = 5
204 else
205 if (size(b, 1) /= k .or. size(b, 2) /= n) flag = 5
206 end if
207 if (flag /= 0) then
208 ! ERROR: Matrix dimensions mismatch
209 write(errmsg, 100) &
210 "Matrix dimension mismatch. Input number ", flag, &
211 " was not sized correctly."
212 call errmgr%report_error("cmtx_mult_mtx", errmsg, &
213 la_array_size_error)
214 return
215 end if
216
217 ! Call ZGEMM
218 call zgemm(ta, tb, m, n, k, alpha, a, lda, b, ldb, beta, c, m)
219
220 ! Formatting
221100 format(a, i0, a)
222 end subroutine
223
224! ------------------------------------------------------------------------------
225 module subroutine cmtx_mult_vec(opa, alpha, a, b, beta, c, err)
226 ! Arguments
227 integer(int32), intent(in) :: opa
228 complex(real64), intent(in) :: alpha, beta
229 complex(real64), intent(in), dimension(:,:) :: a
230 complex(real64), intent(in), dimension(:) :: b
231 complex(real64), intent(inout), dimension(:) :: c
232 class(errors), intent(inout), optional, target :: err
233
234 ! Local Variables
235 character :: t
236 integer(int32) :: m, n, flag
237 class(errors), pointer :: errmgr
238 type(errors), target :: deferr
239 character(len = 128) :: errmsg
240
241 ! Initialization
242 m = size(a, 1)
243 n = size(a, 2)
244 if (opa == la_transpose) then
245 t = 'T'
246 else if (opa == la_hermitian_transpose) then
247 t = 'H'
248 else
249 t = 'N'
250 end if
251 if (present(err)) then
252 errmgr => err
253 else
254 errmgr => deferr
255 end if
256
257 ! Input Check
258 flag = 0
259 if (opa == la_transpose .or. opa == la_hermitian_transpose) then
260 if (size(b) /= m) then
261 flag = 4
262 else if (size(c) /= n) then
263 flag = 6
264 end if
265 else
266 if (size(b) /= n) then
267 flag = 4
268 else if (size(c) /= m) then
269 flag = 6
270 end if
271 end if
272 if (flag /= 0) then
273 ! ERROR: Matrix dimensions mismatch
274 write(errmsg, 100) &
275 "Matrix dimension mismatch. Input number ", flag, &
276 " was not sized correctly."
277 call errmgr%report_error("cmtx_mult_vec", errmsg, &
278 la_array_size_error)
279 return
280 end if
281
282 ! Call ZGEMV
283 call zgemv(t, m, n, alpha, a, m, b, 1, beta, c, 1)
284
285 ! Formatting
286100 format(a, i0, a)
287 end subroutine
288
289! ******************************************************************************
290! RANK 1 UPDATE
291! ------------------------------------------------------------------------------
292 module subroutine rank1_update_dbl(alpha, x, y, a, err)
293 ! Arguments
294 real(real64), intent(in) :: alpha
295 real(real64), intent(in), dimension(:) :: x, y
296 real(real64), intent(inout), dimension(:,:) :: a
297 class(errors), intent(inout), optional, target :: err
298
299 ! Parameters
300 real(real64), parameter :: zero = 0.0d0
301
302 ! Local Variables
303 integer(int32) :: j, m, n
304 real(real64) :: temp
305 class(errors), pointer :: errmgr
306 type(errors), target :: deferr
307
308 ! Initialization
309 m = size(x)
310 n = size(y)
311 if (present(err)) then
312 errmgr => err
313 else
314 errmgr => deferr
315 end if
316
317 ! Input Check
318 if (size(a, 1) /= m .or. size(a, 2) /= n) then
319 ! ERROR: Matrix dimension array
320 call errmgr%report_error("rank1_update_dbl", &
321 "Matrix dimension mismatch.", la_array_size_error)
322 return
323 end if
324
325 ! Process
326 do j = 1, n
327 if (y(j) /= zero) then
328 temp = alpha * y(j)
329 a(:,j) = a(:,j) + temp * x
330 end if
331 end do
332 end subroutine
333
334! xxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxx !
335! COMPLEX VALUED VERSIONS !
336! xxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxx !
337 module subroutine rank1_update_cmplx(alpha, x, y, a, err)
338 ! Arguments
339 complex(real64), intent(in) :: alpha
340 complex(real64), intent(in), dimension(:) :: x, y
341 complex(real64), intent(inout), dimension(:,:) :: a
342 class(errors), intent(inout), optional, target :: err
343
344 ! Parameters
345 complex(real64), parameter :: zero = (0.0d0, 0.0d0)
346
347 ! Local Variables
348 integer(int32) :: j, m, n
349 complex(real64) :: temp
350 class(errors), pointer :: errmgr
351 type(errors), target :: deferr
352
353 ! Initialization
354 m = size(x)
355 n = size(y)
356 if (present(err)) then
357 errmgr => err
358 else
359 errmgr => deferr
360 end if
361
362 ! Input Check
363 if (size(a, 1) /= m .or. size(a, 2) /= n) then
364 ! ERROR: Matrix dimension array
365 call errmgr%report_error("rank1_update_cmplx", &
366 "Matrix dimension mismatch.", la_array_size_error)
367 return
368 end if
369
370 ! Process
371 do j = 1, n
372 if (y(j) /= zero) then
373 temp = alpha * conjg(y(j))
374 a(:,j) = a(:,j) + temp * x
375 end if
376 end do
377 end subroutine
378
379! ******************************************************************************
380! DIAGONAL MATRIX MULTIPLICATION ROUTINES
381! ------------------------------------------------------------------------------
382 module subroutine diag_mtx_mult_mtx(lside, trans, alpha, a, b, beta, c, err)
383 ! Arguments
384 logical, intent(in) :: lside, trans
385 real(real64) :: alpha, beta
386 real(real64), intent(in), dimension(:) :: a
387 real(real64), intent(in), dimension(:,:) :: b
388 real(real64), intent(inout), dimension(:,:) :: c
389 class(errors), intent(inout), optional, target :: err
390
391 ! Parameters
392 real(real64), parameter :: zero = 0.0d0
393 real(real64), parameter :: one = 1.0d0
394
395 ! Local Variables
396 integer(int32) :: i, m, n, k, nrowb, ncolb, flag
397 real(real64) :: temp
398 class(errors), pointer :: errmgr
399 type(errors), target :: deferr
400 character(len = 128) :: errmsg
401
402 ! Initialization
403 m = size(c, 1)
404 n = size(c, 2)
405 k = size(a)
406 nrowb = size(b, 1)
407 ncolb = size(b, 2)
408 if (present(err)) then
409 errmgr => err
410 else
411 errmgr => deferr
412 end if
413
414 ! Input Check
415 flag = 0
416 if (lside) then
417 if (k > m) then
418 flag = 4
419 else
420 if (trans) then
421 ! Compute C = alpha * A * B**T + beta * C
422 if (nrowb /= n .or. ncolb < k) flag = 5
423 else
424 ! Compute C = alpha * A * B + beta * C
425 if (nrowb < k .or. ncolb /= n) flag = 5
426 end if
427 end if
428 else
429 if (k > n) then
430 flag = 4
431 else
432 if (trans) then
433 ! Compute C = alpha * B**T * A + beta * C
434 if (ncolb /= m .or. nrowb < k) flag = 5
435 else
436 ! Compute C = alpha * B * A + beta * C
437 if (nrowb /= m .or. ncolb < k) flag = 5
438 end if
439 end if
440 end if
441 if (flag /= 0) then
442 ! ERROR: One of the input arrays is not sized correctly
443 write(errmsg, 100) "Input number ", flag, &
444 " is not sized correctly."
445 call errmgr%report_error("diag_mtx_mult_mtx", trim(errmsg), &
446 la_array_size_error)
447 return
448 end if
449
450 ! Deal with ALPHA == 0
451 if (alpha == 0) then
452 if (beta == zero) then
453 c = zero
454 else if (beta /= one) then
455 c = beta * c
456 end if
457 return
458 end if
459
460 ! Process
461 if (lside) then
462 if (trans) then
463 ! Compute C = alpha * A * B**T + beta * C
464 do i = 1, k
465 if (beta == zero) then
466 c(i,:) = zero
467 else if (beta /= one) then
468 c(i,:) = beta * c(i,:)
469 end if
470 temp = alpha * a(i)
471 if (temp /= one) c(i,:) = c(i,:) + temp * b(:,i)
472 end do
473 else
474 ! Compute C = alpha * A * B + beta * C
475 do i = 1, k
476 if (beta == zero) then
477 c(i,:) = zero
478 else if (beta /= one) then
479 c(i,:) = beta * c(i,:)
480 end if
481 temp = alpha * a(i)
482 if (temp /= one) c(i,:) = c(i,:) + temp * b(i,:)
483 end do
484 end if
485
486 ! Handle extra rows
487 if (m > k) then
488 if (beta == zero) then
489 c(k+1:m,:) = zero
490 else
491 c(k+1:m,:) = beta * c(k+1:m,:)
492 end if
493 end if
494 else
495 if (trans) then
496 ! Compute C = alpha * B**T * A + beta * C
497 do i = 1, k
498 if (beta == zero) then
499 c(:,i) = zero
500 else if (beta /= one) then
501 c(:,i) = beta * c(:,i)
502 end if
503 temp = alpha * a(i)
504 if (temp /= one) c(:,i) = c(:,i) + temp * b(i,:)
505 end do
506 else
507 ! Compute C = alpha * B * A + beta * C
508 do i = 1, k
509 if (beta == zero) then
510 c(:,i) = zero
511 else if (beta /= one) then
512 c(:,i) = beta * c(:,i)
513 end if
514 temp = alpha * a(i)
515 if (temp /= one) c(:,i) = c(:,i) + temp * b(:,i)
516 end do
517 end if
518
519 ! Handle extra columns
520 if (n > k) then
521 if (beta == zero) then
522 c(:,k+1:m) = zero
523 else if (beta /= one) then
524 c(:,k+1:m) = beta * c(:,k+1:m)
525 end if
526 end if
527 end if
528
529 ! Formatting
530100 format(a, i0, a)
531 end subroutine
532
533! ------------------------------------------------------------------------------
534 module subroutine diag_mtx_mult_mtx2(lside, alpha, a, b, err)
535 ! Arguments
536 logical, intent(in) :: lside
537 real(real64), intent(in) :: alpha
538 real(real64), intent(in), dimension(:) :: a
539 real(real64), intent(inout), dimension(:,:) :: b
540 class(errors), intent(inout), optional, target :: err
541
542 ! Parameters
543 real(real64), parameter :: zero = 0.0d0
544 real(real64), parameter :: one = 1.0d0
545
546 ! Local Variables
547 integer(int32) :: i, m, n, k
548 real(real64) :: temp
549 class(errors), pointer :: errmgr
550 type(errors), target :: deferr
551
552 ! Initialization
553 m = size(b, 1)
554 n = size(b, 2)
555 k = size(a)
556 if (present(err)) then
557 errmgr => err
558 else
559 errmgr => deferr
560 end if
561
562 ! Input Check
563 if ((lside .and. k > m) .or. (.not.lside .and. k > n)) then
564 ! ERROR: One of the input arrays is not sized correctly
565 call errmgr%report_error("diag_mtx_mult_mtx2", &
566 "Input number 3 is not sized correctly.", &
567 la_array_size_error)
568 return
569 end if
570
571 ! Process
572 if (lside) then
573 ! Compute B = alpha * A * B
574 do i = 1, k
575 temp = alpha * a(i)
576 if (temp /= one) b(i,:) = temp * b(i,:)
577 end do
578 if (m > k) b(k+1:m,:) = zero
579 else
580 ! Compute B = alpha * B * A
581 do i = 1, k
582 temp = alpha * a(i)
583 if (temp /= one) b(:,i) = temp * b(:,i)
584 end do
585 if (n > k) b(:,k+1:n) = zero
586 end if
587 end subroutine
588
589! ------------------------------------------------------------------------------
590 module subroutine diag_mtx_mult_mtx3(lside, trans, alpha, a, b, beta, c, err)
591 ! Arguments
592 logical, intent(in) :: lside, trans
593 real(real64) :: alpha, beta
594 complex(real64), intent(in), dimension(:) :: a
595 real(real64), intent(in), dimension(:,:) :: b
596 complex(real64), intent(inout), dimension(:,:) :: c
597 class(errors), intent(inout), optional, target :: err
598
599 ! Parameters
600 complex(real64), parameter :: zero = (0.0d0, 0.0d0)
601 complex(real64), parameter :: one = (1.0d0, 0.0d0)
602
603 ! Local Variables
604 integer(int32) :: i, m, n, k, nrowb, ncolb, flag
605 complex(real64) :: temp
606 class(errors), pointer :: errmgr
607 type(errors), target :: deferr
608 character(len = 128) :: errmsg
609
610 ! Initialization
611 m = size(c, 1)
612 n = size(c, 2)
613 k = size(a)
614 nrowb = size(b, 1)
615 ncolb = size(b, 2)
616 if (present(err)) then
617 errmgr => err
618 else
619 errmgr => deferr
620 end if
621
622 ! Input Check
623 flag = 0
624 if (lside) then
625 if (k > m) then
626 flag = 4
627 else
628 if (trans) then
629 ! Compute C = alpha * A * B**T + beta * C
630 if (nrowb /= n .or. ncolb < k) flag = 5
631 else
632 ! Compute C = alpha * A * B + beta * C
633 if (nrowb < k .or. ncolb /= n) flag = 5
634 end if
635 end if
636 else
637 if (k > n) then
638 flag = 4
639 else
640 if (trans) then
641 ! Compute C = alpha * B**T * A + beta * C
642 if (ncolb /= m .or. nrowb < k) flag = 5
643 else
644 ! Compute C = alpha * B * A + beta * C
645 if (nrowb /= m .or. ncolb < k) flag = 5
646 end if
647 end if
648 end if
649 if (flag /= 0) then
650 ! ERROR: One of the input arrays is not sized correctly
651 write(errmsg, 100) "Input number ", flag, &
652 " is not sized correctly."
653 call errmgr%report_error("diag_mtx_mult_mtx3", trim(errmsg), &
654 la_array_size_error)
655 return
656 end if
657
658 ! Deal with ALPHA == 0
659 if (alpha == 0) then
660 if (beta == zero) then
661 c = zero
662 else if (beta /= one) then
663 c = beta * c
664 end if
665 return
666 end if
667
668 ! Process
669 if (lside) then
670 if (trans) then
671 ! Compute C = alpha * A * B**T + beta * C
672 do i = 1, k
673 if (beta == zero) then
674 c(i,:) = zero
675 else if (beta /= one) then
676 c(i,:) = beta * c(i,:)
677 end if
678 temp = alpha * a(i)
679 if (temp /= one) c(i,:) = c(i,:) + temp * b(:,i)
680 end do
681 else
682 ! Compute C = alpha * A * B + beta * C
683 do i = 1, k
684 if (beta == zero) then
685 c(i,:) = zero
686 else if (beta /= one) then
687 c(i,:) = beta * c(i,:)
688 end if
689 temp = alpha * a(i)
690 if (temp /= one) c(i,:) = c(i,:) + temp * b(i,:)
691 end do
692 end if
693
694 ! Handle extra rows
695 if (m > k) then
696 if (beta == zero) then
697 c(k+1:m,:) = zero
698 else
699 c(k+1:m,:) = beta * c(k+1:m,:)
700 end if
701 end if
702 else
703 if (trans) then
704 ! Compute C = alpha * B**T * A + beta * C
705 do i = 1, k
706 if (beta == zero) then
707 c(:,i) = zero
708 else if (beta /= one) then
709 c(:,i) = beta * c(:,i)
710 end if
711 temp = alpha * a(i)
712 if (temp /= one) c(:,i) = c(:,i) + temp * b(i,:)
713 end do
714 else
715 ! Compute C = alpha * B * A + beta * C
716 do i = 1, k
717 if (beta == zero) then
718 c(:,i) = zero
719 else if (beta /= one) then
720 c(:,i) = beta * c(:,i)
721 end if
722 temp = alpha * a(i)
723 if (temp /= one) c(:,i) = c(:,i) + temp * b(:,i)
724 end do
725 end if
726
727 ! Handle extra columns
728 if (n > k) then
729 if (beta == zero) then
730 c(:,k+1:m) = zero
731 else if (beta /= one) then
732 c(:,k+1:m) = beta * c(:,k+1:m)
733 end if
734 end if
735 end if
736
737 ! Formatting
738100 format(a, i0, a)
739 end subroutine
740
741! ------------------------------------------------------------------------------
742 module subroutine diag_mtx_mult_mtx4(lside, opb, alpha, a, b, beta, c, err)
743 ! Arguments
744 logical, intent(in) :: lside
745 integer(int32), intent(in) :: opb
746 real(real64) :: alpha, beta
747 complex(real64), intent(in), dimension(:) :: a
748 complex(real64), intent(in), dimension(:,:) :: b
749 complex(real64), intent(inout), dimension(:,:) :: c
750 class(errors), intent(inout), optional, target :: err
751
752 ! Parameters
753 complex(real64), parameter :: zero = (0.0d0, 0.0d0)
754 complex(real64), parameter :: one = (1.0d0, 0.0d0)
755
756 ! Local Variables
757 integer(int32) :: i, m, n, k, nrowb, ncolb, flag
758 complex(real64) :: temp
759 class(errors), pointer :: errmgr
760 type(errors), target :: deferr
761 character(len = 128) :: errmsg
762
763 ! Initialization
764 m = size(c, 1)
765 n = size(c, 2)
766 k = size(a)
767 nrowb = size(b, 1)
768 ncolb = size(b, 2)
769 if (present(err)) then
770 errmgr => err
771 else
772 errmgr => deferr
773 end if
774
775 ! Input Check
776 flag = 0
777 if (lside) then
778 if (k > m) then
779 flag = 4
780 else
781 if (opb == la_transpose .or. opb == la_hermitian_transpose) then
782 ! Compute C = alpha * A * B**T + beta * C
783 if (nrowb /= n .or. ncolb < k) flag = 5
784 else
785 ! Compute C = alpha * A * B + beta * C
786 if (nrowb < k .or. ncolb /= n) flag = 5
787 end if
788 end if
789 else
790 if (k > n) then
791 flag = 4
792 else
793 if (opb == la_transpose .or. opb == la_hermitian_transpose) then
794 ! Compute C = alpha * B**T * A + beta * C
795 if (ncolb /= m .or. nrowb < k) flag = 5
796 else
797 ! Compute C = alpha * B * A + beta * C
798 if (nrowb /= m .or. ncolb < k) flag = 5
799 end if
800 end if
801 end if
802 if (flag /= 0) then
803 ! ERROR: One of the input arrays is not sized correctly
804 write(errmsg, 100) "Input number ", flag, &
805 " is not sized correctly."
806 call errmgr%report_error("diag_mtx_mult_mtx4", trim(errmsg), &
807 la_array_size_error)
808 return
809 end if
810
811 ! Deal with ALPHA == 0
812 if (alpha == 0) then
813 if (beta == zero) then
814 c = zero
815 else if (beta /= one) then
816 c = beta * c
817 end if
818 return
819 end if
820
821 ! Process
822 if (lside) then
823 if (opb == la_transpose) then
824 ! Compute C = alpha * A * B**T + beta * C
825 do i = 1, k
826 if (beta == zero) then
827 c(i,:) = zero
828 else if (beta /= one) then
829 c(i,:) = beta * c(i,:)
830 end if
831 temp = alpha * a(i)
832 if (temp /= one) c(i,:) = c(i,:) + temp * b(:,i)
833 end do
834 else if (opb == la_hermitian_transpose) then
835 ! Compute C = alpha * A * B**H + beta * C
836 do i = 1, k
837 if (beta == zero) then
838 c(i,:) = zero
839 else if (beta /= one) then
840 c(i,:) = beta * c(i,:)
841 end if
842 temp = alpha * a(i)
843 if (temp /= one) c(i,:) = c(i,:) + temp * conjg(b(:,i))
844 end do
845 else
846 ! Compute C = alpha * A * B + beta * C
847 do i = 1, k
848 if (beta == zero) then
849 c(i,:) = zero
850 else if (beta /= one) then
851 c(i,:) = beta * c(i,:)
852 end if
853 temp = alpha * a(i)
854 if (temp /= one) c(i,:) = c(i,:) + temp * b(i,:)
855 end do
856 end if
857
858 ! Handle extra rows
859 if (m > k) then
860 if (beta == zero) then
861 c(k+1:m,:) = zero
862 else
863 c(k+1:m,:) = beta * c(k+1:m,:)
864 end if
865 end if
866 else
867 if (opb == la_transpose) then
868 ! Compute C = alpha * B**T * A + beta * C
869 do i = 1, k
870 if (beta == zero) then
871 c(:,i) = zero
872 else if (beta /= one) then
873 c(:,i) = beta * c(:,i)
874 end if
875 temp = alpha * a(i)
876 if (temp /= one) c(:,i) = c(:,i) + temp * b(i,:)
877 end do
878 else if (opb == la_hermitian_transpose) then
879 ! Compute C = alpha * B**H * A + beta * C
880 do i = 1, k
881 if (beta == zero) then
882 c(:,i) = zero
883 else if (beta /= one) then
884 c(:,i) = beta * c(:,i)
885 end if
886 temp = alpha * a(i)
887 if (temp /= one) c(:,i) = c(:,i) + temp * conjg(b(i,:))
888 end do
889 else
890 ! Compute C = alpha * B * A + beta * C
891 do i = 1, k
892 if (beta == zero) then
893 c(:,i) = zero
894 else if (beta /= one) then
895 c(:,i) = beta * c(:,i)
896 end if
897 temp = alpha * a(i)
898 if (temp /= one) c(:,i) = c(:,i) + temp * b(:,i)
899 end do
900 end if
901
902 ! Handle extra columns
903 if (n > k) then
904 if (beta == zero) then
905 c(:,k+1:m) = zero
906 else if (beta /= one) then
907 c(:,k+1:m) = beta * c(:,k+1:m)
908 end if
909 end if
910 end if
911
912 ! Formatting
913100 format(a, i0, a)
914 end subroutine
915
916! ------------------------------------------------------------------------------
917 module subroutine diag_mtx_mult_mtx_cmplx(lside, opb, alpha, a, b, beta, c, err)
918 ! Arguments
919 logical, intent(in) :: lside
920 integer(int32), intent(in) :: opb
921 complex(real64) :: alpha, beta
922 complex(real64), intent(in), dimension(:) :: a
923 complex(real64), intent(in), dimension(:,:) :: b
924 complex(real64), intent(inout), dimension(:,:) :: c
925 class(errors), intent(inout), optional, target :: err
926
927 ! Parameters
928 complex(real64), parameter :: zero = (0.0d0, 0.0d0)
929 complex(real64), parameter :: one = (1.0d0, 0.0d0)
930
931 ! Local Variables
932 integer(int32) :: i, m, n, k, nrowb, ncolb, flag
933 complex(real64) :: temp
934 class(errors), pointer :: errmgr
935 type(errors), target :: deferr
936 character(len = 128) :: errmsg
937
938 ! Initialization
939 m = size(c, 1)
940 n = size(c, 2)
941 k = size(a)
942 nrowb = size(b, 1)
943 ncolb = size(b, 2)
944 if (present(err)) then
945 errmgr => err
946 else
947 errmgr => deferr
948 end if
949
950 ! Input Check
951 flag = 0
952 if (lside) then
953 if (k > m) then
954 flag = 4
955 else
956 if (opb == la_transpose .or. opb == la_hermitian_transpose) then
957 ! Compute C = alpha * A * B**T + beta * C
958 if (nrowb /= n .or. ncolb < k) flag = 5
959 else
960 ! Compute C = alpha * A * B + beta * C
961 if (nrowb < k .or. ncolb /= n) flag = 5
962 end if
963 end if
964 else
965 if (k > n) then
966 flag = 4
967 else
968 if (opb == la_transpose .or. opb == la_hermitian_transpose) then
969 ! Compute C = alpha * B**T * A + beta * C
970 if (ncolb /= m .or. nrowb < k) flag = 5
971 else
972 ! Compute C = alpha * B * A + beta * C
973 if (nrowb /= m .or. ncolb < k) flag = 5
974 end if
975 end if
976 end if
977 if (flag /= 0) then
978 ! ERROR: One of the input arrays is not sized correctly
979 write(errmsg, 100) "Input number ", flag, &
980 " is not sized correctly."
981 call errmgr%report_error("diag_mtx_mult_mtx_cmplx", trim(errmsg), &
982 la_array_size_error)
983 return
984 end if
985
986 ! Deal with ALPHA == 0
987 if (alpha == 0) then
988 if (beta == zero) then
989 c = zero
990 else if (beta /= one) then
991 c = beta * c
992 end if
993 return
994 end if
995
996 ! Process
997 if (lside) then
998 if (opb == la_transpose) then
999 ! Compute C = alpha * A * B**T + beta * C
1000 do i = 1, k
1001 if (beta == zero) then
1002 c(i,:) = zero
1003 else if (beta /= one) then
1004 c(i,:) = beta * c(i,:)
1005 end if
1006 temp = alpha * a(i)
1007 if (temp /= one) c(i,:) = c(i,:) + temp * b(:,i)
1008 end do
1009 else if (opb == la_hermitian_transpose) then
1010 ! Compute C = alpha * A * B**H + beta * C
1011 do i = 1, k
1012 if (beta == zero) then
1013 c(i,:) = zero
1014 else if (beta /= one) then
1015 c(i,:) = beta * c(i,:)
1016 end if
1017 temp = alpha * a(i)
1018 if (temp /= one) c(i,:) = c(i,:) + temp * conjg(b(:,i))
1019 end do
1020 else
1021 ! Compute C = alpha * A * B + beta * C
1022 do i = 1, k
1023 if (beta == zero) then
1024 c(i,:) = zero
1025 else if (beta /= one) then
1026 c(i,:) = beta * c(i,:)
1027 end if
1028 temp = alpha * a(i)
1029 if (temp /= one) c(i,:) = c(i,:) + temp * b(i,:)
1030 end do
1031 end if
1032
1033 ! Handle extra rows
1034 if (m > k) then
1035 if (beta == zero) then
1036 c(k+1:m,:) = zero
1037 else
1038 c(k+1:m,:) = beta * c(k+1:m,:)
1039 end if
1040 end if
1041 else
1042 if (opb == la_transpose) then
1043 ! Compute C = alpha * B**T * A + beta * C
1044 do i = 1, k
1045 if (beta == zero) then
1046 c(:,i) = zero
1047 else if (beta /= one) then
1048 c(:,i) = beta * c(:,i)
1049 end if
1050 temp = alpha * a(i)
1051 if (temp /= one) c(:,i) = c(:,i) + temp * b(i,:)
1052 end do
1053 else if (opb == la_hermitian_transpose) then
1054 ! Compute C = alpha * B**H * A + beta * C
1055 do i = 1, k
1056 if (beta == zero) then
1057 c(:,i) = zero
1058 else if (beta /= one) then
1059 c(:,i) = beta * c(:,i)
1060 end if
1061 temp = alpha * a(i)
1062 if (temp /= one) c(:,i) = c(:,i) + temp * conjg(b(i,:))
1063 end do
1064 else
1065 ! Compute C = alpha * B * A + beta * C
1066 do i = 1, k
1067 if (beta == zero) then
1068 c(:,i) = zero
1069 else if (beta /= one) then
1070 c(:,i) = beta * c(:,i)
1071 end if
1072 temp = alpha * a(i)
1073 if (temp /= one) c(:,i) = c(:,i) + temp * b(:,i)
1074 end do
1075 end if
1076
1077 ! Handle extra columns
1078 if (n > k) then
1079 if (beta == zero) then
1080 c(:,k+1:m) = zero
1081 else if (beta /= one) then
1082 c(:,k+1:m) = beta * c(:,k+1:m)
1083 end if
1084 end if
1085 end if
1086
1087 ! Formatting
1088100 format(a, i0, a)
1089 end subroutine
1090
1091! ------------------------------------------------------------------------------
1092 module subroutine diag_mtx_mult_mtx2_cmplx(lside, alpha, a, b, err)
1093 ! Arguments
1094 logical, intent(in) :: lside
1095 complex(real64), intent(in) :: alpha
1096 complex(real64), intent(in), dimension(:) :: a
1097 complex(real64), intent(inout), dimension(:,:) :: b
1098 class(errors), intent(inout), optional, target :: err
1099
1100 ! Parameters
1101 complex(real64), parameter :: zero = (0.0d0, 0.0d0)
1102 complex(real64), parameter :: one = (1.0d0, 0.0d0)
1103
1104 ! Local Variables
1105 integer(int32) :: i, m, n, k
1106 complex(real64) :: temp
1107 class(errors), pointer :: errmgr
1108 type(errors), target :: deferr
1109
1110 ! Initialization
1111 m = size(b, 1)
1112 n = size(b, 2)
1113 k = size(a)
1114 if (present(err)) then
1115 errmgr => err
1116 else
1117 errmgr => deferr
1118 end if
1119
1120 ! Input Check
1121 if ((lside .and. k > m) .or. (.not.lside .and. k > n)) then
1122 ! ERROR: One of the input arrays is not sized correctly
1123 call errmgr%report_error("diag_mtx_mult_mtx2_cmplx", &
1124 "Input number 3 is not sized correctly.", &
1125 la_array_size_error)
1126 return
1127 end if
1128
1129 ! Process
1130 if (lside) then
1131 ! Compute B = alpha * A * B
1132 do i = 1, k
1133 temp = alpha * a(i)
1134 if (temp /= one) b(i,:) = temp * b(i,:)
1135 end do
1136 if (m > k) b(k+1:m,:) = zero
1137 else
1138 ! Compute B = alpha * B * A
1139 do i = 1, k
1140 temp = alpha * a(i)
1141 if (temp /= one) b(:,i) = temp * b(:,i)
1142 end do
1143 if (n > k) b(:,k+1:n) = zero
1144 end if
1145 end subroutine
1146
1147! ------------------------------------------------------------------------------
1148 module subroutine diag_mtx_mult_mtx_mix(lside, opb, alpha, a, b, beta, c, err)
1149 ! Arguments
1150 logical, intent(in) :: lside
1151 integer(int32), intent(in) :: opb
1152 complex(real64) :: alpha, beta
1153 real(real64), intent(in), dimension(:) :: a
1154 complex(real64), intent(in), dimension(:,:) :: b
1155 complex(real64), intent(inout), dimension(:,:) :: c
1156 class(errors), intent(inout), optional, target :: err
1157
1158 ! Parameters
1159 complex(real64), parameter :: zero = (0.0d0, 0.0d0)
1160 complex(real64), parameter :: one = (1.0d0, 0.0d0)
1161
1162 ! Local Variables
1163 integer(int32) :: i, m, n, k, nrowb, ncolb, flag
1164 complex(real64) :: temp
1165 class(errors), pointer :: errmgr
1166 type(errors), target :: deferr
1167 character(len = 128) :: errmsg
1168
1169 ! Initialization
1170 m = size(c, 1)
1171 n = size(c, 2)
1172 k = size(a)
1173 nrowb = size(b, 1)
1174 ncolb = size(b, 2)
1175 if (present(err)) then
1176 errmgr => err
1177 else
1178 errmgr => deferr
1179 end if
1180
1181 ! Input Check
1182 flag = 0
1183 if (lside) then
1184 if (k > m) then
1185 flag = 4
1186 else
1187 if (opb == la_transpose .or. opb == la_hermitian_transpose) then
1188 ! Compute C = alpha * A * B**T + beta * C
1189 if (nrowb /= n .or. ncolb < k) flag = 5
1190 else
1191 ! Compute C = alpha * A * B + beta * C
1192 if (nrowb < k .or. ncolb /= n) flag = 5
1193 end if
1194 end if
1195 else
1196 if (k > n) then
1197 flag = 4
1198 else
1199 if (opb == la_transpose .or. opb == la_hermitian_transpose) then
1200 ! Compute C = alpha * B**T * A + beta * C
1201 if (ncolb /= m .or. nrowb < k) flag = 5
1202 else
1203 ! Compute C = alpha * B * A + beta * C
1204 if (nrowb /= m .or. ncolb < k) flag = 5
1205 end if
1206 end if
1207 end if
1208 if (flag /= 0) then
1209 ! ERROR: One of the input arrays is not sized correctly
1210 write(errmsg, 100) "Input number ", flag, &
1211 " is not sized correctly."
1212 call errmgr%report_error("diag_mtx_mult_mtx_mix", trim(errmsg), &
1213 la_array_size_error)
1214 return
1215 end if
1216
1217 ! Deal with ALPHA == 0
1218 if (alpha == 0) then
1219 if (beta == zero) then
1220 c = zero
1221 else if (beta /= one) then
1222 c = beta * c
1223 end if
1224 return
1225 end if
1226
1227 ! Process
1228 if (lside) then
1229 if (opb == la_transpose) then
1230 ! Compute C = alpha * A * B**T + beta * C
1231 do i = 1, k
1232 if (beta == zero) then
1233 c(i,:) = zero
1234 else if (beta /= one) then
1235 c(i,:) = beta * c(i,:)
1236 end if
1237 temp = alpha * a(i)
1238 if (temp /= one) c(i,:) = c(i,:) + temp * b(:,i)
1239 end do
1240 else if (opb == la_hermitian_transpose) then
1241 ! Compute C = alpha * A * B**H + beta * C
1242 do i = 1, k
1243 if (beta == zero) then
1244 c(i,:) = zero
1245 else if (beta /= one) then
1246 c(i,:) = beta * c(i,:)
1247 end if
1248 temp = alpha * a(i)
1249 if (temp /= one) c(i,:) = c(i,:) + temp * conjg(b(:,i))
1250 end do
1251 else
1252 ! Compute C = alpha * A * B + beta * C
1253 do i = 1, k
1254 if (beta == zero) then
1255 c(i,:) = zero
1256 else if (beta /= one) then
1257 c(i,:) = beta * c(i,:)
1258 end if
1259 temp = alpha * a(i)
1260 if (temp /= one) c(i,:) = c(i,:) + temp * b(i,:)
1261 end do
1262 end if
1263
1264 ! Handle extra rows
1265 if (m > k) then
1266 if (beta == zero) then
1267 c(k+1:m,:) = zero
1268 else
1269 c(k+1:m,:) = beta * c(k+1:m,:)
1270 end if
1271 end if
1272 else
1273 if (opb == la_transpose) then
1274 ! Compute C = alpha * B**T * A + beta * C
1275 do i = 1, k
1276 if (beta == zero) then
1277 c(:,i) = zero
1278 else if (beta /= one) then
1279 c(:,i) = beta * c(:,i)
1280 end if
1281 temp = alpha * a(i)
1282 if (temp /= one) c(:,i) = c(:,i) + temp * b(i,:)
1283 end do
1284 else if (opb == la_hermitian_transpose) then
1285 ! Compute C = alpha * B**H * A + beta * C
1286 do i = 1, k
1287 if (beta == zero) then
1288 c(:,i) = zero
1289 else if (beta /= one) then
1290 c(:,i) = beta * c(:,i)
1291 end if
1292 temp = alpha * a(i)
1293 if (temp /= one) c(:,i) = c(:,i) + temp * conjg(b(i,:))
1294 end do
1295 else
1296 ! Compute C = alpha * B * A + beta * C
1297 do i = 1, k
1298 if (beta == zero) then
1299 c(:,i) = zero
1300 else if (beta /= one) then
1301 c(:,i) = beta * c(:,i)
1302 end if
1303 temp = alpha * a(i)
1304 if (temp /= one) c(:,i) = c(:,i) + temp * b(:,i)
1305 end do
1306 end if
1307
1308 ! Handle extra columns
1309 if (n > k) then
1310 if (beta == zero) then
1311 c(:,k+1:m) = zero
1312 else if (beta /= one) then
1313 c(:,k+1:m) = beta * c(:,k+1:m)
1314 end if
1315 end if
1316 end if
1317
1318 ! Formatting
1319100 format(a, i0, a)
1320 end subroutine
1321
1322! ------------------------------------------------------------------------------
1323 module subroutine diag_mtx_mult_mtx2_mix(lside, alpha, a, b, err)
1324 ! Arguments
1325 logical, intent(in) :: lside
1326 complex(real64), intent(in) :: alpha
1327 real(real64), intent(in), dimension(:) :: a
1328 complex(real64), intent(inout), dimension(:,:) :: b
1329 class(errors), intent(inout), optional, target :: err
1330
1331 ! Parameters
1332 complex(real64), parameter :: zero = (0.0d0, 0.0d0)
1333 complex(real64), parameter :: one = (1.0d0, 0.0d0)
1334
1335 ! Local Variables
1336 integer(int32) :: i, m, n, k
1337 complex(real64) :: temp
1338 class(errors), pointer :: errmgr
1339 type(errors), target :: deferr
1340
1341 ! Initialization
1342 m = size(b, 1)
1343 n = size(b, 2)
1344 k = size(a)
1345 if (present(err)) then
1346 errmgr => err
1347 else
1348 errmgr => deferr
1349 end if
1350
1351 ! Input Check
1352 if ((lside .and. k > m) .or. (.not.lside .and. k > n)) then
1353 ! ERROR: One of the input arrays is not sized correctly
1354 call errmgr%report_error("diag_mtx_mult_mtx2_cmplx", &
1355 "Input number 3 is not sized correctly.", &
1356 la_array_size_error)
1357 return
1358 end if
1359
1360 ! Process
1361 if (lside) then
1362 ! Compute B = alpha * A * B
1363 do i = 1, k
1364 temp = alpha * a(i)
1365 if (temp /= one) b(i,:) = temp * b(i,:)
1366 end do
1367 if (m > k) b(k+1:m,:) = zero
1368 else
1369 ! Compute B = alpha * B * A
1370 do i = 1, k
1371 temp = alpha * a(i)
1372 if (temp /= one) b(:,i) = temp * b(:,i)
1373 end do
1374 if (n > k) b(:,k+1:n) = zero
1375 end if
1376 end subroutine
1377
1378! ******************************************************************************
1379! BASIC OPERATION ROUTINES
1380! ------------------------------------------------------------------------------
1381 pure module function trace_dbl(x) result(y)
1382 ! Arguments
1383 real(real64), intent(in), dimension(:,:) :: x
1384 real(real64) :: y
1385
1386 ! Parameters
1387 real(real64), parameter :: zero = 0.0d0
1388
1389 ! Local Variables
1390 integer(int32) :: i, m, n, mn
1391
1392 ! Initialization
1393 y = zero
1394 m = size(x, 1)
1395 n = size(x, 2)
1396 mn = min(m, n)
1397
1398 ! Process
1399 do i = 1, mn
1400 y = y + x(i,i)
1401 end do
1402 end function
1403
1404! ------------------------------------------------------------------------------
1405 pure module function trace_cmplx(x) result(y)
1406 ! Arguments
1407 complex(real64), intent(in), dimension(:,:) :: x
1408 complex(real64) :: y
1409
1410 ! Parameters
1411 complex(real64), parameter :: zero = (0.0d0, 0.0d0)
1412
1413 ! Local Variables
1414 integer(int32) :: i, m, n, mn
1415
1416 ! Initialization
1417 y = zero
1418 m = size(x, 1)
1419 n = size(x, 2)
1420 mn = min(m, n)
1421
1422 ! Process
1423 do i = 1, mn
1424 y = y + x(i,i)
1425 end do
1426 end function
1427
1428! ------------------------------------------------------------------------------
1429 module function mtx_rank_dbl(a, tol, work, olwork, err) result(rnk)
1430 ! Arguments
1431 real(real64), intent(inout), dimension(:,:) :: a
1432 real(real64), intent(in), optional :: tol
1433 real(real64), intent(out), target, optional, dimension(:) :: work
1434 integer(int32), intent(out), optional :: olwork
1435 class(errors), intent(inout), optional, target :: err
1436 integer(int32) :: rnk
1437
1438 ! External Function Interfaces
1439 interface
1440 function dlamch(cmach) result(x)
1441 use, intrinsic :: iso_fortran_env, only : real64
1442 character, intent(in) :: cmach
1443 real(real64) :: x
1444 end function
1445 end interface
1446
1447 ! Local Variables
1448 integer(int32) :: i, m, n, mn, istat, lwork, flag
1449 real(real64), pointer, dimension(:) :: wptr, s, w
1450 real(real64), allocatable, target, dimension(:) :: wrk
1451 real(real64) :: t, tref, smlnum
1452 real(real64), dimension(1) :: dummy, temp
1453 class(errors), pointer :: errmgr
1454 type(errors), target :: deferr
1455 character(len = 128) :: errmsg
1456
1457 ! Initialization
1458 m = size(a, 1)
1459 n = size(a, 2)
1460 mn = min(m, n)
1461 smlnum = dlamch('s')
1462 rnk = 0
1463 if (present(err)) then
1464 errmgr => err
1465 else
1466 errmgr => deferr
1467 end if
1468
1469 ! Workspace Query
1470 !call svd(a, a(1:mn,1), olwork = lwork)
1471 call dgesvd('N', 'N', m, n, a, m, dummy, dummy, m, dummy, n, temp, &
1472 -1, flag)
1473 lwork = int(temp(1), int32) + mn
1474 if (present(olwork)) then
1475 olwork = lwork
1476 return
1477 end if
1478
1479 ! Local Memory Allocation
1480 if (present(work)) then
1481 if (size(work) < lwork) then
1482 ! ERROR: WORK not sized correctly
1483 call errmgr%report_error("mtx_rank", &
1484 "Incorrectly sized input array WORK, argument 5.", &
1485 la_array_size_error)
1486 return
1487 end if
1488 wptr => work(1:lwork)
1489 else
1490 allocate(wrk(lwork), stat = istat)
1491 if (istat /= 0) then
1492 ! ERROR: Out of memory
1493 call errmgr%report_error("mtx_rank", &
1494 "Insufficient memory available.", &
1495 la_out_of_memory_error)
1496 return
1497 end if
1498 wptr => wrk
1499 end if
1500 s => wptr(1:mn)
1501 w => wptr(mn+1:lwork)
1502
1503 ! Compute the singular values of A
1504 call dgesvd('N', 'N', m, n, a, m, s, dummy, m, dummy, n, w, &
1505 lwork - mn, flag)
1506 if (flag > 0) then
1507 write(errmsg, 100) flag, " superdiagonals could not " // &
1508 "converge to zero as part of the QR iteration process."
1509 call errmgr%report_warning("mtx_rank", errmsg, la_convergence_error)
1510 end if
1511
1512 ! Determine the threshold tolerance for the singular values such that
1513 ! singular values less than the threshold result in zero when inverted.
1514 tref = max(m, n) * epsilon(t) * s(1)
1515 if (present(tol)) then
1516 t = tol
1517 else
1518 t = tref
1519 end if
1520 if (t < smlnum) then
1521 ! ! The supplied tolerance is too small, simply fall back to the
1522 ! ! default, but issue a warning to the user
1523 ! t = tref
1524 ! call report_warning("mtx_rank", "The supplied tolerance was " // &
1525 ! "smaller than a value that would result in an overflow " // &
1526 ! "condition, or is negative; therefore, the tolerance has " // &
1527 ! "been reset to its default value.")
1528 end if
1529
1530 ! Count the singular values that are larger than the tolerance value
1531 do i = 1, mn
1532 if (s(i) < t) exit
1533 rnk = rnk + 1
1534 end do
1535
1536 ! Formatting
1537100 format(i0, a)
1538 end function
1539
1540! ------------------------------------------------------------------------------
1541 module function mtx_rank_cmplx(a, tol, work, olwork, rwork, err) result(rnk)
1542 ! Arguments
1543 complex(real64), intent(inout), dimension(:,:) :: a
1544 real(real64), intent(in), optional :: tol
1545 complex(real64), intent(out), target, optional, dimension(:) :: work
1546 integer(int32), intent(out), optional :: olwork
1547 real(real64), intent(out), target, optional, dimension(:) :: rwork
1548 class(errors), intent(inout), optional, target :: err
1549 integer(int32) :: rnk
1550
1551 ! External Function Interfaces
1552 interface
1553 function dlamch(cmach) result(x)
1554 use, intrinsic :: iso_fortran_env, only : real64
1555 character, intent(in) :: cmach
1556 real(real64) :: x
1557 end function
1558 end interface
1559
1560 ! Local Variables
1561 integer(int32) :: i, m, n, mn, istat, lwork, flag, lrwork
1562 real(real64), pointer, dimension(:) :: s, rwptr, rw
1563 real(real64), allocatable, target, dimension(:) :: rwrk
1564 complex(real64), allocatable, target, dimension(:) :: wrk
1565 complex(real64), pointer, dimension(:) :: wptr
1566 real(real64) :: t, tref, smlnum
1567 real(real64), dimension(1) :: dummy
1568 complex(real64), dimension(1) :: cdummy, temp
1569 class(errors), pointer :: errmgr
1570 type(errors), target :: deferr
1571 character(len = 128) :: errmsg
1572
1573 ! Initialization
1574 m = size(a, 1)
1575 n = size(a, 2)
1576 mn = min(m, n)
1577 lrwork = 6 * mn
1578 smlnum = dlamch('s')
1579 rnk = 0
1580 if (present(err)) then
1581 errmgr => err
1582 else
1583 errmgr => deferr
1584 end if
1585
1586 ! Workspace Query
1587 call zgesvd('N', 'N', m, n, a, m, dummy, cdummy, m, cdummy, n, temp, &
1588 -1, dummy, flag)
1589 lwork = int(temp(1), int32)
1590 if (present(olwork)) then
1591 olwork = lwork
1592 return
1593 end if
1594
1595 ! Local Memory Allocation
1596 if (present(work)) then
1597 if (size(work) < lwork) then
1598 ! ERROR: WORK not sized correctly
1599 call errmgr%report_error("mtx_rank_cmplx", &
1600 "Incorrectly sized input array WORK, argument 5.", &
1601 la_array_size_error)
1602 return
1603 end if
1604 wptr => work(1:lwork)
1605 else
1606 allocate(wrk(lwork), stat = istat)
1607 if (istat /= 0) then
1608 ! ERROR: Out of memory
1609 call errmgr%report_error("mtx_rank_cmplx", &
1610 "Insufficient memory available.", &
1611 la_out_of_memory_error)
1612 return
1613 end if
1614 wptr => wrk
1615 end if
1616
1617 if (present(rwork)) then
1618 if (size(rwork) < lrwork) then
1619 ! ERROR: RWORK not sized correctly
1620 call errmgr%report_error("mtx_rank_cmplx", &
1621 "Incorrectly sized input array RWORK.", &
1622 la_array_size_error)
1623 return
1624 end if
1625 rwptr => rwork(1:lrwork)
1626 else
1627 allocate(rwrk(lrwork), stat = istat)
1628 if (istat /= 0) then
1629 end if
1630 rwptr => rwrk
1631 end if
1632 s => rwptr(1:mn)
1633 rw => rwptr(mn+1:lrwork)
1634
1635 ! Compute the singular values of A
1636 call zgesvd('N', 'N', m, n, a, m, s, cdummy, m, cdummy, n, wptr, &
1637 lwork - mn, rw, flag)
1638 if (flag > 0) then
1639 write(errmsg, 100) flag, " superdiagonals could not " // &
1640 "converge to zero as part of the QR iteration process."
1641 call errmgr%report_warning("mtx_rank_cmplx", errmsg, la_convergence_error)
1642 end if
1643
1644 ! Determine the threshold tolerance for the singular values such that
1645 ! singular values less than the threshold result in zero when inverted.
1646 tref = max(m, n) * epsilon(t) * s(1)
1647 if (present(tol)) then
1648 t = tol
1649 else
1650 t = tref
1651 end if
1652 if (t < smlnum) then
1653 ! ! The supplied tolerance is too small, simply fall back to the
1654 ! ! default, but issue a warning to the user
1655 ! t = tref
1656 ! call report_warning("mtx_rank", "The supplied tolerance was " // &
1657 ! "smaller than a value that would result in an overflow " // &
1658 ! "condition, or is negative; therefore, the tolerance has " // &
1659 ! "been reset to its default value.")
1660 end if
1661
1662 ! Count the singular values that are larger than the tolerance value
1663 do i = 1, mn
1664 if (s(i) < t) exit
1665 rnk = rnk + 1
1666 end do
1667
1668 ! Formatting
1669100 format(i0, a)
1670 end function
1671
1672! ------------------------------------------------------------------------------
1673 module function det_dbl(a, iwork, err) result(x)
1674 ! Arguments
1675 real(real64), intent(inout), dimension(:,:) :: a
1676 integer(int32), intent(out), target, optional, dimension(:) :: iwork
1677 class(errors), intent(inout), optional, target :: err
1678 real(real64) :: x
1679
1680 ! Parameters
1681 real(real64), parameter :: zero = 0.0d0
1682 real(real64), parameter :: one = 1.0d0
1683 real(real64), parameter :: ten = 1.0d1
1684 real(real64), parameter :: p1 = 1.0d-1
1685
1686 ! Local Variables
1687 integer(int32) :: i, ep, n, istat, flag
1688 integer(int32), pointer, dimension(:) :: ipvt
1689 integer(int32), allocatable, target, dimension(:) :: iwrk
1690 real(real64) :: temp
1691 class(errors), pointer :: errmgr
1692 type(errors), target :: deferr
1693
1694 ! Initialization
1695 n = size(a, 1)
1696 x = zero
1697 if (present(err)) then
1698 errmgr => err
1699 else
1700 errmgr => deferr
1701 end if
1702
1703 ! Input Check
1704 if (size(a, 2) /= n) then
1705 call errmgr%report_error("det", &
1706 "The supplied matrix must be square.", la_array_size_error)
1707 return
1708 end if
1709
1710 ! Local Memory Allocation
1711 if (present(iwork)) then
1712 if (size(iwork) < n) then
1713 ! ERROR: WORK not sized correctly
1714 call errmgr%report_error("det", &
1715 "Incorrectly sized input array IWORK, argument 2.", &
1716 la_array_size_error)
1717 return
1718 end if
1719 ipvt => iwork(1:n)
1720 else
1721 allocate(iwrk(n), stat = istat)
1722 if (istat /= 0) then
1723 ! ERROR: Out of memory
1724 call errmgr%report_error("det", &
1725 "Insufficient memory available.", &
1726 la_out_of_memory_error)
1727 return
1728 end if
1729 ipvt => iwrk
1730 end if
1731
1732 ! Compute the LU factorization of A
1733 call dgetrf(n, n, a, n, ipvt, flag)
1734 if (flag > 0) then
1735 ! A singular matrix has a determinant of zero
1736 x = zero
1737 return
1738 end if
1739
1740 ! Compute the product of the diagonal of A
1741 temp = one
1742 ep = 0
1743 do i = 1, n
1744 if (ipvt(i) /= i) temp = -temp
1745
1746 temp = a(i,i) * temp
1747 if (temp == zero) then
1748 x = zero
1749 exit
1750 end if
1751
1752 do while (abs(temp) < one)
1753 temp = ten * temp
1754 ep = ep - 1
1755 end do
1756
1757 do while (abs(temp) > ten)
1758 temp = p1 * temp
1759 ep = ep + 1
1760 end do
1761 end do
1762 x = temp * ten**ep
1763 end function
1764
1765! ------------------------------------------------------------------------------
1766 module function det_cmplx(a, iwork, err) result(x)
1767 ! Arguments
1768 complex(real64), intent(inout), dimension(:,:) :: a
1769 integer(int32), intent(out), target, optional, dimension(:) :: iwork
1770 class(errors), intent(inout), optional, target :: err
1771 complex(real64) :: x
1772
1773 ! Parameters
1774 complex(real64), parameter :: zero = (0.0d0, 0.0d0)
1775 complex(real64), parameter :: one = (1.0d0, 0.0d0)
1776 complex(real64), parameter :: ten = (1.0d1, 0.0d0)
1777 complex(real64), parameter :: p1 = (1.0d-1, 0.0d0)
1778 real(real64), parameter :: real_one = 1.0d0
1779 real(real64), parameter :: real_ten = 1.0d1
1780
1781 ! Local Variables
1782 integer(int32) :: i, ep, n, istat, flag
1783 integer(int32), pointer, dimension(:) :: ipvt
1784 integer(int32), allocatable, target, dimension(:) :: iwrk
1785 complex(real64) :: temp
1786 class(errors), pointer :: errmgr
1787 type(errors), target :: deferr
1788
1789 ! Initialization
1790 n = size(a, 1)
1791 x = zero
1792 if (present(err)) then
1793 errmgr => err
1794 else
1795 errmgr => deferr
1796 end if
1797
1798 ! Input Check
1799 if (size(a, 2) /= n) then
1800 call errmgr%report_error("det_cmplx", &
1801 "The supplied matrix must be square.", la_array_size_error)
1802 return
1803 end if
1804
1805 ! Local Memory Allocation
1806 if (present(iwork)) then
1807 if (size(iwork) < n) then
1808 ! ERROR: WORK not sized correctly
1809 call errmgr%report_error("det_cmplx", &
1810 "Incorrectly sized input array IWORK, argument 2.", &
1811 la_array_size_error)
1812 return
1813 end if
1814 ipvt => iwork(1:n)
1815 else
1816 allocate(iwrk(n), stat = istat)
1817 if (istat /= 0) then
1818 ! ERROR: Out of memory
1819 call errmgr%report_error("det_cmplx", &
1820 "Insufficient memory available.", &
1821 la_out_of_memory_error)
1822 return
1823 end if
1824 ipvt => iwrk
1825 end if
1826
1827 ! Compute the LU factorization of A
1828 call zgetrf(n, n, a, n, ipvt, flag)
1829 if (flag > 0) then
1830 ! A singular matrix has a determinant of zero
1831 x = zero
1832 return
1833 end if
1834
1835 ! Compute the product of the diagonal of A
1836 temp = one
1837 ep = 0
1838 do i = 1, n
1839 if (ipvt(i) /= i) temp = -temp
1840
1841 temp = a(i,i) * temp
1842 if (temp == zero) then
1843 x = zero
1844 exit
1845 end if
1846
1847 do while (abs(temp) < real_one)
1848 temp = ten * temp
1849 ep = ep - 1
1850 end do
1851
1852 do while (abs(temp) > real_ten)
1853 temp = p1 * temp
1854 ep = ep + 1
1855 end do
1856 end do
1857 x = temp * ten**ep
1858 end function
1859
1860! ******************************************************************************
1861! ARRAY SWAPPING ROUTINE
1862! ------------------------------------------------------------------------------
1863 module subroutine swap_dbl(x, y, err)
1864 ! Arguments
1865 real(real64), intent(inout), dimension(:) :: x, y
1866 class(errors), intent(inout), optional, target :: err
1867
1868 ! Local Variables
1869 integer(int32) :: i, n
1870 real(real64) :: temp
1871 class(errors), pointer :: errmgr
1872 type(errors), target :: deferr
1873
1874 ! Initialization
1875 n = size(x)
1876 if (present(err)) then
1877 errmgr => err
1878 else
1879 errmgr => deferr
1880 end if
1881
1882 ! Input Check
1883 if (size(y) /= n) then
1884 call errmgr%report_error("swap", &
1885 "The input arrays are not the same size.", &
1886 la_array_size_error)
1887 return
1888 end if
1889
1890 ! Process
1891 do i = 1, n
1892 temp = x(i)
1893 x(i) = y(i)
1894 y(i) = temp
1895 end do
1896 end subroutine
1897
1898! ------------------------------------------------------------------------------
1899 module subroutine swap_cmplx(x, y, err)
1900 ! Arguments
1901 complex(real64), intent(inout), dimension(:) :: x, y
1902 class(errors), intent(inout), optional, target :: err
1903
1904 ! Local Variables
1905 integer(int32) :: i, n
1906 complex(real64) :: temp
1907 class(errors), pointer :: errmgr
1908 type(errors), target :: deferr
1909
1910 ! Initialization
1911 n = size(x)
1912 if (present(err)) then
1913 errmgr => err
1914 else
1915 errmgr => deferr
1916 end if
1917
1918 ! Input Check
1919 if (size(y) /= n) then
1920 call errmgr%report_error("swap_cmplx", &
1921 "The input arrays are not the same size.", &
1922 la_array_size_error)
1923 return
1924 end if
1925
1926 ! Process
1927 do i = 1, n
1928 temp = x(i)
1929 x(i) = y(i)
1930 y(i) = temp
1931 end do
1932 end subroutine
1933
1934! ******************************************************************************
1935! ARRAY MULTIPLICIATION ROUTINES
1936! ------------------------------------------------------------------------------
1937 module subroutine recip_mult_array_dbl(a, x)
1938 ! Arguments
1939 real(real64), intent(in) :: a
1940 real(real64), intent(inout), dimension(:) :: x
1941
1942 ! External Function Interfaces
1943 interface
1944 function dlamch(cmach) result(x)
1945 use, intrinsic :: iso_fortran_env, only : real64
1946 character, intent(in) :: cmach
1947 real(real64) :: x
1948 end function
1949 end interface
1950
1951 ! Parameters
1952 real(real64), parameter :: zero = 0.0d0
1953 real(real64), parameter :: one = 1.0d0
1954 real(real64), parameter :: twotho = 2.0d3
1955
1956 ! Local Variables
1957 logical :: done
1958 real(real64) :: bignum, cden, cden1, cnum, cnum1, mul, smlnum
1959
1960 ! Initialization
1961 smlnum = dlamch('s')
1962 bignum = one / smlnum
1963 if (log10(bignum) > twotho) then
1964 smlnum = sqrt(smlnum)
1965 bignum = sqrt(bignum)
1966 end if
1967
1968 ! Initialize the denominator to A, and the numerator to ONE
1969 cden = a
1970 cnum = one
1971
1972 ! Process
1973 do
1974 cden1 = cden * smlnum
1975 cnum1 = cnum / bignum
1976 if (abs(cden1) > abs(cnum) .and. cnum /= zero) then
1977 mul = smlnum
1978 done = .false.
1979 cden = cden1
1980 else if (abs(cnum1) > abs(cden)) then
1981 mul = bignum
1982 done = .false.
1983 cnum = cnum1
1984 else
1985 mul = cnum / cden
1986 done = .true.
1987 end if
1988
1989 ! Scale the vector X by MUL
1990 x = mul * x
1991
1992 ! Exit if done
1993 if (done) exit
1994 end do
1995 end subroutine
1996
1997! ******************************************************************************
1998! TRIANGULAR MATRIX MULTIPLICATION ROUTINES
1999! ------------------------------------------------------------------------------
2000 module subroutine tri_mtx_mult_dbl(upper, alpha, a, beta, b, err)
2001 ! Arguments
2002 logical, intent(in) :: upper
2003 real(real64), intent(in) :: alpha, beta
2004 real(real64), intent(in), dimension(:,:) :: a
2005 real(real64), intent(inout), dimension(:,:) :: b
2006 class(errors), intent(inout), optional, target :: err
2007
2008 ! Parameters
2009 real(real64), parameter :: zero = 0.0d0
2010
2011 ! Local Variables
2012 integer(int32) :: i, j, k, n, d1, d2, flag
2013 real(real64) :: temp
2014 class(errors), pointer :: errmgr
2015 type(errors), target :: deferr
2016 character(len = 128) :: errmsg
2017
2018 ! Initialization
2019 n = size(a, 1)
2020 d1 = n
2021 d2 = n
2022 if (present(err)) then
2023 errmgr => err
2024 else
2025 errmgr => deferr
2026 end if
2027
2028 ! Input Check
2029 flag = 0
2030 if (size(a, 2) /= n) then
2031 flag = 3
2032 d2 = size(a, 2)
2033 else if (size(b, 1) /= n .or. size(b, 2) /= n) then
2034 flag = 5
2035 d1 = size(b, 1)
2036 d2 = size(b, 2)
2037 end if
2038 if (flag /= 0) then
2039 ! ERROR: Incorrectly sized matrix
2040 write(errmsg, 100) "The matrix at input ", flag, &
2041 " was not sized appropriately. A matrix of ", n, "-by-", n, &
2042 "was expected, but a matrix of ", d1, "-by-", d2, " was found."
2043 call errmgr%report_error("tri_mtx_mult_dbl", trim(errmsg), &
2044 la_array_size_error)
2045 return
2046 end if
2047
2048 ! Process
2049 if (upper) then
2050 ! Form: B = alpha * A**T * A + beta * B
2051 if (beta == zero) then
2052 do j = 1, n
2053 do i = 1, j
2054 temp = zero
2055 do k = 1, j
2056 temp = temp + a(k,i) * a(k,j)
2057 end do
2058 temp = alpha * temp
2059 b(i,j) = temp
2060 if (i /= j) b(j,i) = temp
2061 end do
2062 end do
2063 else
2064 do j = 1, n
2065 do i = 1, j
2066 temp = zero
2067 do k = 1, j
2068 temp = temp + a(k,i) * a(k,j)
2069 end do
2070 temp = alpha * temp
2071 b(i,j) = temp + beta * b(i,j)
2072 if (i /= j) b(j,i) = temp + beta * b(j,i)
2073 end do
2074 end do
2075 end if
2076 else
2077 ! Form: B = alpha * A * A**T + beta * B
2078 if (beta == zero) then
2079 do j = 1, n
2080 do i = j, n
2081 temp = zero
2082 do k = 1, j
2083 temp = temp + a(i,k) * a(j,k)
2084 end do
2085 temp = alpha * temp
2086 b(i,j) = temp
2087 if (i /= j) b(j,i) = temp
2088 end do
2089 end do
2090 else
2091 do j = 1, n
2092 do i = j, n
2093 temp = zero
2094 do k = 1, j
2095 temp = temp + a(i,k) * a(j,k)
2096 end do
2097 temp = alpha * temp
2098 b(i,j) = temp + beta * b(i,j)
2099 if (i /= j) b(j,i) = temp + beta * b(j,i)
2100 end do
2101 end do
2102 end if
2103 end if
2104
2105 ! Formatting
2106100 format(a, i0, a, i0, a, i0, a, i0, a, i0, a)
2107 end subroutine
2108
2109! ------------------------------------------------------------------------------
2110 module subroutine tri_mtx_mult_cmplx(upper, alpha, a, beta, b, err)
2111 ! Arguments
2112 logical, intent(in) :: upper
2113 complex(real64), intent(in) :: alpha, beta
2114 complex(real64), intent(in), dimension(:,:) :: a
2115 complex(real64), intent(inout), dimension(:,:) :: b
2116 class(errors), intent(inout), optional, target :: err
2117
2118 ! Parameters
2119 complex(real64), parameter :: zero = (0.0d0, 0.0d0)
2120
2121 ! Local Variables
2122 integer(int32) :: i, j, k, n, d1, d2, flag
2123 complex(real64) :: temp
2124 class(errors), pointer :: errmgr
2125 type(errors), target :: deferr
2126 character(len = 128) :: errmsg
2127
2128 ! Initialization
2129 n = size(a, 1)
2130 d1 = n
2131 d2 = n
2132 if (present(err)) then
2133 errmgr => err
2134 else
2135 errmgr => deferr
2136 end if
2137
2138 ! Input Check
2139 flag = 0
2140 if (size(a, 2) /= n) then
2141 flag = 3
2142 d2 = size(a, 2)
2143 else if (size(b, 1) /= n .or. size(b, 2) /= n) then
2144 flag = 5
2145 d1 = size(b, 1)
2146 d2 = size(b, 2)
2147 end if
2148 if (flag /= 0) then
2149 ! ERROR: Incorrectly sized matrix
2150 write(errmsg, 100) "The matrix at input ", flag, &
2151 " was not sized appropriately. A matrix of ", n, "-by-", n, &
2152 "was expected, but a matrix of ", d1, "-by-", d2, " was found."
2153 call errmgr%report_error("tri_mtx_mult_cmplx", trim(errmsg), &
2154 la_array_size_error)
2155 return
2156 end if
2157
2158 ! Process
2159 if (upper) then
2160 ! Form: B = alpha * A**T * A + beta * B
2161 if (beta == zero) then
2162 do j = 1, n
2163 do i = 1, j
2164 temp = zero
2165 do k = 1, j
2166 temp = temp + a(k,i) * a(k,j)
2167 end do
2168 temp = alpha * temp
2169 b(i,j) = temp
2170 if (i /= j) b(j,i) = temp
2171 end do
2172 end do
2173 else
2174 do j = 1, n
2175 do i = 1, j
2176 temp = zero
2177 do k = 1, j
2178 temp = temp + a(k,i) * a(k,j)
2179 end do
2180 temp = alpha * temp
2181 b(i,j) = temp + beta * b(i,j)
2182 if (i /= j) b(j,i) = temp + beta * b(j,i)
2183 end do
2184 end do
2185 end if
2186 else
2187 ! Form: B = alpha * A * A**T + beta * B
2188 if (beta == zero) then
2189 do j = 1, n
2190 do i = j, n
2191 temp = zero
2192 do k = 1, j
2193 temp = temp + a(i,k) * a(j,k)
2194 end do
2195 temp = alpha * temp
2196 b(i,j) = temp
2197 if (i /= j) b(j,i) = temp
2198 end do
2199 end do
2200 else
2201 do j = 1, n
2202 do i = j, n
2203 temp = zero
2204 do k = 1, j
2205 temp = temp + a(i,k) * a(j,k)
2206 end do
2207 temp = alpha * temp
2208 b(i,j) = temp + beta * b(i,j)
2209 if (i /= j) b(j,i) = temp + beta * b(j,i)
2210 end do
2211 end do
2212 end if
2213 end if
2214
2215 ! Formatting
2216100 format(a, i0, a, i0, a, i0, a, i0, a, i0, a)
2217 end subroutine
2218
2219! ------------------------------------------------------------------------------
2220end submodule
Provides a set of common linear algebra routines.
Definition: linalg.f90:145