15 use,
intrinsic :: iso_fortran_env, only : int32, real64
19 public :: mat_rank1_update
26 public :: mat_qr_rank1_update
28 public :: mat_cholesky
29 public :: mat_cholesky_rank1_update
30 public :: mat_cholesky_rank1_downdate
32 public :: mat_pinverse
49 module procedure :: mat_mult_diag_1
50 module procedure :: mat_mult_diag_2
51 module procedure :: mat_mult_diag_3
52 module procedure :: mat_mult_diag_1_cmplx
53 module procedure :: mat_mult_diag_2_cmplx
54 module procedure :: mat_mult_diag_3_cmplx
61 module procedure :: mat_mult_upper_tri_1
62 module procedure :: mat_mult_upper_tri_2
63 module procedure :: mat_mult_upper_tri_1_cmplx
64 module procedure :: mat_mult_upper_tri_2_cmplx
71 module procedure :: mat_mult_lower_tri_1
72 module procedure :: mat_mult_lower_tri_2
73 module procedure :: mat_mult_lower_tri_1_cmplx
74 module procedure :: mat_mult_lower_tri_2_cmplx
81 module procedure :: mat_solve_upper_tri_1
82 module procedure :: mat_solve_upper_tri_2
83 module procedure :: mat_solve_upper_tri_1_cmplx
84 module procedure :: mat_solve_upper_tri_2_cmplx
91 module procedure :: mat_solve_lower_tri_1
92 module procedure :: mat_solve_lower_tri_2
93 module procedure :: mat_solve_lower_tri_1_cmplx
94 module procedure :: mat_solve_lower_tri_2_cmplx
101 module procedure :: mat_lu_dbl
102 module procedure :: mat_lu_cmplx
109 module procedure :: mat_eigen_1
110 module procedure :: mat_eigen_2
117 real(real64),
allocatable,
dimension(:,:) :: l
119 real(real64),
allocatable,
dimension(:,:) :: u
121 real(real64),
allocatable,
dimension(:,:) :: p
128 complex(real64),
allocatable,
dimension(:,:) :: l
130 complex(real64),
allocatable,
dimension(:,:) :: u
132 real(real64),
allocatable,
dimension(:,:) :: p
139 real(real64),
allocatable,
dimension(:,:) :: q
141 real(real64),
allocatable,
dimension(:,:) :: r
144 real(real64),
allocatable,
dimension(:,:) :: p
151 complex(real64),
allocatable,
dimension(:,:) :: q
153 complex(real64),
allocatable,
dimension(:,:) :: r
156 complex(real64),
allocatable,
dimension(:,:) :: p
164 real(real64),
allocatable,
dimension(:,:) :: u
166 real(real64),
allocatable,
dimension(:,:) :: s
168 real(real64),
allocatable,
dimension(:,:) :: vt
176 complex(real64),
allocatable,
dimension(:,:) :: u
178 real(real64),
allocatable,
dimension(:,:) :: s
180 complex(real64),
allocatable,
dimension(:,:) :: vt
188 complex(real64),
allocatable,
dimension(:) :: values
191 complex(real64),
allocatable,
dimension(:,:) :: vectors
204 function mat_rank1_update(a, x, y)
result(b)
206 real(real64),
intent(in),
dimension(:,:) :: a
207 real(real64),
intent(in),
dimension(:) :: x, y
208 real(real64),
dimension(size(a, 1), size(a, 2)) :: b
223 function mat_mult_diag_1(a, b)
result(c)
225 real(real64),
intent(in),
dimension(:) :: a
226 real(real64),
intent(in),
dimension(:,:) :: b
227 real(real64),
dimension(size(a), size(b, 2)) :: c
230 if (
size(b, 1) >
size(a))
then
231 call diag_mtx_mult(.true., .false., 1.0d0, a, b(1:
size(a),:), &
246 function mat_mult_diag_2(a, b)
result(c)
248 real(real64),
intent(in),
dimension(:) :: a, b
249 real(real64),
dimension(size(a)) :: c
252 real(real64),
dimension(size(a), 1) :: bc, cc
255 bc(:,1) = b(1:min(
size(a),
size(b)))
268 function mat_mult_diag_3(a, b)
result(c)
270 real(real64),
intent(in),
dimension(:,:) :: a
271 real(real64),
intent(in),
dimension(:) :: b
272 real(real64),
dimension(size(a, 1), size(b)) :: c
275 if (
size(a, 2) >
size(b))
then
276 call diag_mtx_mult(.false., .false., 1.0d0, b, a(:,1:
size(b)), &
291 function mat_mult_diag_1_cmplx(a, b)
result(c)
293 complex(real64),
intent(in),
dimension(:) :: a
294 complex(real64),
intent(in),
dimension(:,:) :: b
295 complex(real64),
dimension(size(a), size(b, 2)) :: c
298 complex(real64),
parameter :: zero = (0.0d0, 0.0d0)
299 complex(real64),
parameter :: one = (1.0d0, 0.0d0)
302 if (
size(b, 1) >
size(a))
then
303 call diag_mtx_mult(.true., la_no_operation, one, a, b(1:
size(a),:), &
306 call diag_mtx_mult(.true., la_no_operation, one, a, b, zero, c)
318 function mat_mult_diag_2_cmplx(a, b)
result(c)
320 complex(real64),
intent(in),
dimension(:) :: a, b
321 complex(real64),
dimension(size(a)) :: c
324 complex(real64),
parameter :: zero = (0.0d0, 0.0d0)
325 complex(real64),
parameter :: one = (1.0d0, 0.0d0)
328 complex(real64),
dimension(size(a), 1) :: bc, cc
331 bc(:,1) = b(1:min(
size(a),
size(b)))
332 call diag_mtx_mult(.true., la_no_operation, one, a, bc, zero, cc)
344 function mat_mult_diag_3_cmplx(a, b)
result(c)
346 complex(real64),
intent(in),
dimension(:,:) :: a
347 complex(real64),
intent(in),
dimension(:) :: b
348 complex(real64),
dimension(size(a, 1), size(b)) :: c
351 if (
size(a, 2) >
size(b))
then
352 call diag_mtx_mult(.false., la_no_operation, 1.0d0, b, a(:,1:
size(b)), &
355 call diag_mtx_mult(.false., la_no_operation, 1.0d0, b, a, 0.0d0, c)
366 function mat_mult_upper_tri_1(a, b)
result(c)
368 real(real64),
intent(in),
dimension(:,:) :: a, b
369 real(real64),
dimension(size(a, 1), size(b, 2)) :: c
373 call dtrmm(
'L',
'U',
'N',
'N',
size(b, 1),
size(b, 2), 1.0d0, &
374 a,
size(a, 1), c,
size(c, 1))
384 function mat_mult_upper_tri_2(a, b)
result(c)
386 real(real64),
intent(in),
dimension(:,:) :: a
387 real(real64),
intent(in),
dimension(:) :: b
388 real(real64),
dimension(size(a, 1)) :: c
392 call dtrmv(
'U',
'N',
'N',
size(a, 1), a,
size(a, 1), c, 1)
402 function mat_mult_lower_tri_1(a, b)
result(c)
404 real(real64),
intent(in),
dimension(:,:) :: a, b
405 real(real64),
dimension(size(a, 1), size(b, 2)) :: c
409 call dtrmm(
'L',
'L',
'N',
'N',
size(b, 1),
size(b, 2), 1.0d0, &
410 a,
size(a, 1), c,
size(c, 1))
420 function mat_mult_lower_tri_2(a, b)
result(c)
422 real(real64),
intent(in),
dimension(:,:) :: a
423 real(real64),
intent(in),
dimension(:) :: b
424 real(real64),
dimension(size(a, 1)) :: c
428 call dtrmv(
'L',
'N',
'N',
size(a, 1), a,
size(a, 1), c, 1)
438 function mat_mult_upper_tri_1_cmplx(a, b)
result(c)
440 complex(real64),
intent(in),
dimension(:,:) :: a, b
441 complex(real64),
dimension(size(a, 1), size(b, 2)) :: c
445 call ztrmm(
'L',
'U',
'N',
'N',
size(b, 1),
size(b, 2), 1.0d0, &
446 a,
size(a, 1), c,
size(c, 1))
456 function mat_mult_upper_tri_2_cmplx(a, b)
result(c)
458 complex(real64),
intent(in),
dimension(:,:) :: a
459 complex(real64),
intent(in),
dimension(:) :: b
460 complex(real64),
dimension(size(a, 1)) :: c
464 call ztrmv(
'U',
'N',
'N',
size(a, 1), a,
size(a, 1), c, 1)
474 function mat_mult_lower_tri_1_cmplx(a, b)
result(c)
476 complex(real64),
intent(in),
dimension(:,:) :: a, b
477 complex(real64),
dimension(size(a, 1), size(b, 2)) :: c
481 call ztrmm(
'L',
'L',
'N',
'N',
size(b, 1),
size(b, 2), 1.0d0, &
482 a,
size(a, 1), c,
size(c, 1))
492 function mat_mult_lower_tri_2_cmplx(a, b)
result(c)
494 complex(real64),
intent(in),
dimension(:,:) :: a
495 complex(real64),
intent(in),
dimension(:) :: b
496 complex(real64),
dimension(size(a, 1)) :: c
500 call ztrmv(
'L',
'N',
'N',
size(a, 1), a,
size(a, 1), c, 1)
508 function mat_det(a)
result(x)
510 real(real64),
intent(in),
dimension(:,:) :: a
514 real(real64),
dimension(size(a, 1), size(a, 2)) :: b
527 function mat_lu_dbl(a)
result(x)
529 real(real64),
intent(in),
dimension(:,:) :: a
534 integer(int32),
allocatable,
dimension(:) :: ipvt
548 call form_lu(x%l, ipvt, x%u, x%p)
557 function mat_lu_cmplx(a)
result(x)
559 complex(real64),
intent(in),
dimension(:,:) :: a
564 integer(int32),
allocatable,
dimension(:) :: ipvt
578 call form_lu(x%l, ipvt, x%u, x%p)
591 function mat_qr(a, pvt)
result(x)
593 real(real64),
intent(in),
dimension(:,:) :: a
594 logical,
intent(in),
optional :: pvt
599 integer(int32) :: m, n, mn
600 integer(int32),
allocatable,
dimension(:) :: jpvt
601 real(real64),
allocatable,
dimension(:) :: tau
605 if (
present(pvt)) use_pivot = pvt
620 call form_qr(x%r, tau, jpvt, x%q, x%p)
637 function mat_qr_rank1_update(q, r, x, y)
result(rst)
639 real(real64),
intent(in),
dimension(:,:) :: q, r
640 real(real64),
intent(in),
dimension(:) :: x, y
644 integer(int32) :: i, m, n
645 real(real64),
allocatable,
dimension(:) :: xc, yc
673 function mat_svd(a)
result(x)
675 real(real64),
intent(in),
dimension(:,:) :: a
679 integer(int32) :: i, m, n, mn
680 real(real64),
allocatable,
dimension(:) :: s
681 real(real64),
allocatable,
dimension(:,:) :: ac
695 call svd(ac, s, x%u, x%vt)
714 function mat_cholesky(a, upper)
result(r)
716 real(real64),
intent(in),
dimension(:,:) :: a
717 logical,
intent(in),
optional :: upper
718 real(real64),
dimension(size(a, 1), size(a, 2)) :: r
721 logical :: compute_upper
724 compute_upper = .true.
725 if (
present(upper)) compute_upper = upper
737 function mat_cholesky_rank1_update(a, x)
result(r)
739 real(real64),
intent(in),
dimension(:,:) :: a
740 real(real64),
intent(in),
dimension(:) :: x
741 real(real64),
dimension(size(a, 1), size(a, 2)) :: r
744 real(real64),
dimension(size(x)) :: xc
759 function mat_cholesky_rank1_downdate(a, x)
result(r)
761 real(real64),
intent(in),
dimension(:,:) :: a
762 real(real64),
intent(in),
dimension(:) :: x
763 real(real64),
dimension(size(a, 1), size(a, 2)) :: r
766 real(real64),
dimension(size(x)) :: xc
779 function mat_inverse(a)
result(x)
781 real(real64),
intent(in),
dimension(:,:) :: a
782 real(real64),
dimension(size(a, 2), size(a, 1)) :: x
794 function mat_pinverse(a)
result(x)
796 real(real64),
intent(in),
dimension(:,:) :: a
797 real(real64),
dimension(size(a, 2), size(a, 1)) :: x
800 real(real64),
dimension(size(a, 1), size(a, 2)) :: ac
814 function mat_solve_upper_tri_1(a, b)
result(x)
816 real(real64),
intent(in),
dimension(:,:) :: a, b
817 real(real64),
dimension(size(b, 1), size(b, 2)) :: x
832 function mat_solve_upper_tri_2(a, b)
result(x)
834 real(real64),
intent(in),
dimension(:,:) :: a
835 real(real64),
intent(in),
dimension(:) :: b
836 real(real64),
dimension(size(b)) :: x
850 function mat_solve_upper_tri_1_cmplx(a, b)
result(x)
852 complex(real64),
intent(in),
dimension(:,:) :: a, b
853 complex(real64),
dimension(size(b, 1), size(b, 2)) :: x
858 (1.0d0, 0.0d0), a, x)
868 function mat_solve_upper_tri_2_cmplx(a, b)
result(x)
870 complex(real64),
intent(in),
dimension(:,:) :: a
871 complex(real64),
intent(in),
dimension(:) :: b
872 complex(real64),
dimension(size(b)) :: x
886 function mat_solve_lower_tri_1(a, b)
result(x)
888 real(real64),
intent(in),
dimension(:,:) :: a, b
889 real(real64),
dimension(size(b, 1), size(b, 2)) :: x
904 function mat_solve_lower_tri_2(a, b)
result(x)
906 real(real64),
intent(in),
dimension(:,:) :: a
907 real(real64),
intent(in),
dimension(:) :: b
908 real(real64),
dimension(size(b)) :: x
922 function mat_solve_lower_tri_1_cmplx(a, b)
result(x)
924 complex(real64),
intent(in),
dimension(:,:) :: a, b
925 complex(real64),
dimension(size(b, 1), size(b, 2)) :: x
930 (1.0d0, 0.0d0), a, x)
940 function mat_solve_lower_tri_2_cmplx(a, b)
result(x)
942 complex(real64),
intent(in),
dimension(:,:) :: a
943 complex(real64),
intent(in),
dimension(:) :: b
944 complex(real64),
dimension(size(b)) :: x
958 function mat_eigen_1(a)
result(x)
960 real(real64),
intent(in),
dimension(:,:) :: a
965 real(real64),
dimension(size(a, 1), size(a, 2)) :: ac
969 allocate(x%values(n))
970 allocate(x%vectors(n,n))
974 call eigen(ac, x%values, x%vectors)
977 call sort(x%values, x%vectors, .true.)
988 function mat_eigen_2(a, b)
result(x)
990 real(real64),
intent(in),
dimension(:,:) :: a, b
994 integer(int32) :: i, j, n
995 real(real64),
dimension(size(a, 1), size(a, 2)) :: ac
996 real(real64),
dimension(size(b, 1), size(b, 2)) :: bc
1000 allocate(x%values(n))
1001 allocate(x%vectors(n,n))
1010 call eigen(ac, bc, x%values, vecs = x%vectors)
1013 call sort(x%values, x%vectors, .true.)
1021 pure function identity(n)
result(x)
1022 integer(int32),
intent(in) :: n
1023 real(real64),
dimension(n, n) :: x
Computes the Cholesky factorization of a symmetric, positive definite matrix.
Computes the rank 1 downdate to a Cholesky factored matrix (upper triangular).
Computes the rank 1 update to a Cholesky factored matrix (upper triangular).
Computes the determinant of a square matrix.
Multiplies a diagonal matrix with another matrix or array.
Computes the eigenvalues, and optionally the eigenvectors, of a matrix.
Computes the LU factorization of an M-by-N matrix.
Computes the inverse of a square matrix.
Computes the Moore-Penrose pseudo-inverse of a M-by-N matrix using the singular value decomposition o...
Computes the QR factorization of an M-by-N matrix.
Computes the rank 1 update to an M-by-N QR factored matrix A (M >= N) where , and such that .
Performs the rank-1 update to matrix A such that: , where is an M-by-N matrix, is a scalar,...
Solves a triangular system of equations.
Computes the singular value decomposition of a matrix A. The SVD is defined as: , where is an M-by-M...
Computes the eigenvalues and eigenvectors (right) of a general N-by-N matrix.
Computes the LU factorization of a square matrix. Notice, partial row pivoting is utilized.
Computes the matrix operation: C = A * B, where A is a diagonal matrix.
Computes the matrix operation C = A * B, where A is a lower triangular matrix.
Computes the matrix operation C = A * B, where A is an upper triangular matrix.
Solves the lower triangular system A X = B, where A is a lower triangular matrix.
Solves the upper triangular system A X = B, where A is an upper triangular matrix.
Provides an immutable interface to many of the core linear algebra routines in this library....
Provides a set of common linear algebra routines.
Defines a container for the output of an Eigen analysis of a square matrix.
Defines a container for the output of an LU factorization.
Defines a container for the output of an LU factorization.
Defines a container for the output of a QR factorization.
Defines a container for the output of a QR factorization.
Defines a container for the output of a singular value decomposition of a matrix.
Defines a container for the output of a singular value decomposition of a matrix.