Skip to content

Commit b1dec60

Browse files
committed
cholesky submodule
1 parent 3d2a6ec commit b1dec60

File tree

2 files changed

+35
-17
lines changed

2 files changed

+35
-17
lines changed

src/stdlib_linalg.fypp

Lines changed: 31 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -16,6 +16,7 @@ module stdlib_linalg
1616
implicit none
1717
private
1818

19+
public :: cholesky
1920
public :: det
2021
public :: operator(.det.)
2122
public :: diag
@@ -42,6 +43,36 @@ module stdlib_linalg
4243
! Export linalg error handling
4344
public :: linalg_state_type, linalg_error_handling
4445

46+
interface cholesky
47+
#:for rk,rt,ri in RC_KINDS_TYPES
48+
#:if rk!="xdp"
49+
pure module subroutine stdlib_linalg_${ri}$_cholesky_inplace(a,lower,other_zeroed,err)
50+
!> Input matrix a[m,n]
51+
${rt}$, intent(inout), target :: a(:,:)
52+
!> [optional] is the lower or upper triangular factor required? Default = lower
53+
logical(lk), optional, intent(in) :: lower
54+
!> [optional] should the unused half of the return matrix be zeroed out? Default: yes
55+
logical(lk), optional, intent(in) :: other_zeroed
56+
!> [optional] state return flag. On error if not requested, the code will stop
57+
type(linalg_state_type), optional, intent(out) :: err
58+
end subroutine stdlib_linalg_${ri}$_cholesky_inplace
59+
60+
pure module subroutine stdlib_linalg_${ri}$_cholesky(a,c,lower,other_zeroed,err)
61+
!> Input matrix a[n,n]
62+
${rt}$, intent(in) :: a(:,:)
63+
!> Output matrix with Cholesky factors c[n,n]
64+
${rt}$, intent(out) :: c(:,:)
65+
!> [optional] is the lower or upper triangular factor required? Default = lower
66+
logical(lk), optional, intent(in) :: lower
67+
!> [optional] should the unused half of the return matrix be zeroed out? Default: yes
68+
logical(lk), optional, intent(in) :: other_zeroed
69+
!> [optional] state return flag. On error if not requested, the code will stop
70+
type(linalg_state_type), optional, intent(out) :: err
71+
end subroutine stdlib_linalg_${ri}$_cholesky
72+
#:endif
73+
#:endfor
74+
end interface cholesky
75+
4576
interface diag
4677
!! version: experimental
4778
!!

src/stdlib_linalg_cholesky.fypp

Lines changed: 4 additions & 17 deletions
Original file line numberDiff line numberDiff line change
@@ -1,16 +1,12 @@
11
#:include "common.fypp"
22
#:set RC_KINDS_TYPES = REAL_KINDS_TYPES + CMPLX_KINDS_TYPES
33
! Cholesky factorization of a matrix, based on LAPACK *POTRF functions
4-
module stdlib_linalg_cholesky
4+
submodule (stdlib_linalg) stdlib_linalg_cholesky
55
use stdlib_linalg_constants
66
use stdlib_linalg_lapack, only: potrf
77
use stdlib_linalg_state, only: linalg_state_type, linalg_error_handling, LINALG_ERROR, &
88
LINALG_INTERNAL_ERROR, LINALG_VALUE_ERROR
99
implicit none(type,external)
10-
private
11-
12-
!> Cholesky factorization of a matrix
13-
public :: cholesky
1410

1511
! Returns Lower Cholesky factor
1612
! NumPy: L = cholesky(a)
@@ -22,15 +18,6 @@ module stdlib_linalg_cholesky
2218
! SciPy: [C, lower] = cho_factor(a, lower=False, overwrite_a=False, check_finite=True)
2319
! SciPy: cho_solve(c_and_lower, b, overwrite_b=False, check_finite=True)
2420

25-
interface cholesky
26-
#:for rk,rt,ri in RC_KINDS_TYPES
27-
#:if rk!="xdp"
28-
module procedure stdlib_linalg_${ri}$_cholesky_inplace
29-
module procedure stdlib_linalg_${ri}$_cholesky
30-
#:endif
31-
#:endfor
32-
end interface cholesky
33-
3421
character(*), parameter :: this = 'cholesky'
3522

3623
contains
@@ -65,7 +52,7 @@ module stdlib_linalg_cholesky
6552

6653
! Compute the Cholesky factorization of a symmetric / Hermitian matrix, A = L*L^T = U^T*U.
6754
! The factorization is returned in-place, overwriting matrix A
68-
pure subroutine stdlib_linalg_${ri}$_cholesky_inplace(a,lower,other_zeroed,err)
55+
pure module subroutine stdlib_linalg_${ri}$_cholesky_inplace(a,lower,other_zeroed,err)
6956
!> Input matrix a[m,n]
7057
${rt}$, intent(inout), target :: a(:,:)
7158
!> [optional] is the lower or upper triangular factor required? Default = lower
@@ -127,7 +114,7 @@ module stdlib_linalg_cholesky
127114

128115
! Compute the Cholesky factorization of a symmetric / Hermitian matrix, A = L*L^T = U^T*U.
129116
! The factorization is returned as a separate matrix
130-
pure subroutine stdlib_linalg_${ri}$_cholesky(a,c,lower,other_zeroed,err)
117+
pure module subroutine stdlib_linalg_${ri}$_cholesky(a,c,lower,other_zeroed,err)
131118
!> Input matrix a[n,n]
132119
${rt}$, intent(in) :: a(:,:)
133120
!> Output matrix with Cholesky factors c[n,n]
@@ -172,4 +159,4 @@ module stdlib_linalg_cholesky
172159
#:endif
173160
#:endfor
174161

175-
end module stdlib_linalg_cholesky
162+
end submodule stdlib_linalg_cholesky

0 commit comments

Comments
 (0)