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