1
1
#:include "common.fypp"
2
2
#:set RC_KINDS_TYPES = REAL_KINDS_TYPES + CMPLX_KINDS_TYPES
3
- module stdlib_linalg_pseudoinverse
3
+ submodule(stdlib_linalg) stdlib_linalg_pseudoinverse
4
4
!! Compute the (Moore-Penrose) pseudo-inverse of a matrix.
5
5
use stdlib_linalg_constants
6
6
use stdlib_linalg_blas
@@ -9,35 +9,6 @@ module stdlib_linalg_pseudoinverse
9
9
use stdlib_linalg, only: svd
10
10
use iso_fortran_env,only:real32,real64,real128,int8,int16,int32,int64,stderr => error_unit
11
11
implicit none(type,external)
12
- private
13
-
14
- !> Pseudo-inverse: Function interface
15
- public :: pinv
16
- !> Pseudo-inverse: Subroutine interface (pre-allocated)
17
- public :: pseudoinvert
18
- !> Operator interface: .pinv.A returns the pseudo-inverse of A
19
- public :: operator(.pinv.)
20
-
21
- ! Function interface
22
- interface pinv
23
- #:for rk,rt,ri in RC_KINDS_TYPES
24
- module procedure stdlib_linalg_pseudoinverse_${ri}$
25
- #:endfor
26
- end interface pinv
27
-
28
- ! Subroutine interface
29
- interface pseudoinvert
30
- #:for rk,rt,ri in RC_KINDS_TYPES
31
- module procedure stdlib_linalg_pseudoinvert_${ri}$
32
- #:endfor
33
- end interface pseudoinvert
34
-
35
- ! Operator interface
36
- interface operator(.pinv.)
37
- #:for rk,rt,ri in RC_KINDS_TYPES
38
- module procedure stdlib_linalg_pinv_${ri}$_operator
39
- #:endfor
40
- end interface operator(.pinv.)
41
12
42
13
character(*), parameter :: this = 'pseudo-inverse'
43
14
@@ -46,7 +17,7 @@ module stdlib_linalg_pseudoinverse
46
17
#:for rk,rt,ri in RC_KINDS_TYPES
47
18
48
19
! Compute the in-place pseudo-inverse of matrix a
49
- subroutine stdlib_linalg_pseudoinvert_${ri}$(a,pinva,rtol,err)
20
+ module subroutine stdlib_linalg_pseudoinvert_${ri}$(a,pinva,rtol,err)
50
21
!> Input matrix a[m,n]
51
22
${rt}$, intent(inout) :: a(:,:)
52
23
!> Output pseudo-inverse matrix
@@ -115,7 +86,7 @@ module stdlib_linalg_pseudoinverse
115
86
end subroutine stdlib_linalg_pseudoinvert_${ri}$
116
87
117
88
! Function interface
118
- function stdlib_linalg_pseudoinverse_${ri}$(a,rtol,err) result(pinva)
89
+ module function stdlib_linalg_pseudoinverse_${ri}$(a,rtol,err) result(pinva)
119
90
!> Input matrix a[m,n]
120
91
${rt}$, intent(in), target :: a(:,:)
121
92
!> [optional] ....
@@ -134,7 +105,7 @@ module stdlib_linalg_pseudoinverse
134
105
end function stdlib_linalg_pseudoinverse_${ri}$
135
106
136
107
! Inverse matrix operator
137
- function stdlib_linalg_pinv_${ri}$_operator(a) result(pinva)
108
+ module function stdlib_linalg_pinv_${ri}$_operator(a) result(pinva)
138
109
!> Input matrix a[m,n]
139
110
${rt}$, intent(in), target :: a(:,:)
140
111
!> Result matrix
@@ -150,4 +121,4 @@ module stdlib_linalg_pseudoinverse
150
121
151
122
#:endfor
152
123
153
- end module stdlib_linalg_pseudoinverse
124
+ end submodule stdlib_linalg_pseudoinverse
0 commit comments