Skip to content

Commit 4f64bfa

Browse files
committed
submodule version
1 parent 4beab1f commit 4f64bfa

File tree

2 files changed

+94
-74
lines changed

2 files changed

+94
-74
lines changed

src/stdlib_linalg.fypp

+90-3
Original file line numberDiff line numberDiff line change
@@ -1,18 +1,19 @@
11
#:include "common.fypp"
2-
#:set RCI_KINDS_TYPES = REAL_KINDS_TYPES + CMPLX_KINDS_TYPES + INT_KINDS_TYPES
2+
#:set RC_KINDS_TYPES = REAL_KINDS_TYPES + CMPLX_KINDS_TYPES
3+
#:set RCI_KINDS_TYPES = RC_KINDS_TYPES + INT_KINDS_TYPES
34
module stdlib_linalg
45
!!Provides a support for various linear algebra procedures
56
!! ([Specification](../page/specs/stdlib_linalg.html))
6-
use stdlib_kinds, only: sp, dp, xdp, qp, &
7+
use stdlib_kinds, only: sp, dp, xdp, qp, lk, &
78
int8, int16, int32, int64
89
use stdlib_error, only: error_stop
910
use stdlib_optval, only: optval
1011
use stdlib_linalg_state, only: linalg_state_type, linalg_error_handling
11-
use stdlib_linalg_determinant, only: det
1212
implicit none
1313
private
1414

1515
public :: det
16+
public :: operator(.det.)
1617
public :: diag
1718
public :: eye
1819
public :: trace
@@ -220,6 +221,92 @@ module stdlib_linalg
220221
#:endfor
221222
end interface is_hessenberg
222223

224+
225+
interface det
226+
!!### Summary
227+
!! Interface for computing matrix determinant.
228+
!!
229+
!!### Description
230+
!!
231+
!! This interface provides methods for computing the determinant of a matrix.
232+
!! Supported data types include real and complex.
233+
!!
234+
!!@note The provided functions are intended for square matrices.
235+
!!
236+
!!### Example
237+
!!
238+
!!```fortran
239+
!!
240+
!! real(sp) :: a(3,3), d
241+
!! type(linalg_state_type) :: state
242+
!! a = reshape([1, 2, 3, 4, 5, 6, 7, 8, 9], [3, 3])
243+
!!
244+
!! d = det(a,err=state)
245+
!! if (state%ok()) then
246+
!! print *, 'Success! det=',d
247+
!! else
248+
!! print *, state%print()
249+
!! endif
250+
!!
251+
!!```
252+
!!
253+
#:for rk,rt in RC_KINDS_TYPES
254+
#:if rk!="xdp"
255+
module function stdlib_linalg_${rt[0]}$${rk}$determinant(a,overwrite_a,err) result(det)
256+
!> Input matrix a[m,n]
257+
${rt}$, intent(inout), target :: a(:,:)
258+
!> [optional] Can A data be overwritten and destroyed?
259+
logical(lk), optional, intent(in) :: overwrite_a
260+
!> State return flag.
261+
type(linalg_state_type), intent(out) :: err
262+
!> Matrix determinant
263+
${rt}$ :: det
264+
end function stdlib_linalg_${rt[0]}$${rk}$determinant
265+
module function stdlib_linalg_pure_${rt[0]}$${rk}$determinant(a) result(det)
266+
!> Input matrix a[m,n]
267+
${rt}$, intent(in) :: a(:,:)
268+
!> Matrix determinant
269+
${rt}$ :: det
270+
end function stdlib_linalg_pure_${rt[0]}$${rk}$determinant
271+
272+
#:endif
273+
#:endfor
274+
end interface det
275+
276+
interface operator(.det.)
277+
!!### Summary
278+
!! Pure operator interface for computing matrix determinant.
279+
!!
280+
!!### Description
281+
!!
282+
!! This pure operator interface provides a convenient way to compute the determinant of a matrix.
283+
!! Supported data types include real and complex.
284+
!!
285+
!!@note The provided functions are intended for square matrices.
286+
!!
287+
!!### Example
288+
!!
289+
!!```fortran
290+
!!
291+
!! real(sp) :: matrix(3,3), d
292+
!! matrix = reshape([1, 2, 3, 4, 5, 6, 7, 8, 9], [3, 3])
293+
!! d = .det.matrix
294+
!!
295+
!!```
296+
!
297+
#:for rk,rt in RC_KINDS_TYPES
298+
#:if rk!="xdp"
299+
module function stdlib_linalg_pure_${rt[0]}$${rk}$determinant(a) result(det)
300+
!> Input matrix a[m,n]
301+
${rt}$, intent(in) :: a(:,:)
302+
!> Matrix determinant
303+
${rt}$ :: det
304+
end function stdlib_linalg_pure_${rt[0]}$${rk}$determinant
305+
#:endif
306+
#:endfor
307+
end interface operator(.det.)
308+
309+
223310
contains
224311

225312

src/stdlib_linalg_determinant.fypp

+4-71
Original file line numberDiff line numberDiff line change
@@ -1,6 +1,6 @@
11
#:include "common.fypp"
22
#:set RC_KINDS_TYPES = REAL_KINDS_TYPES + CMPLX_KINDS_TYPES
3-
module stdlib_linalg_determinant
3+
submodule (stdlib_linalg) stdlib_linalg_determinant
44
!! Determinant of a rectangular matrix
55
use stdlib_linalg_constants
66
use stdlib_linalg_lapack, only: getrf
@@ -10,80 +10,13 @@ module stdlib_linalg_determinant
1010
private
1111

1212
! Function interface
13-
public :: det
14-
public :: operator(.det.)
15-
1613
character(*), parameter :: this = 'determinant'
1714

18-
interface det
19-
!!### Summary
20-
!! Interface for computing matrix determinant.
21-
!!
22-
!!### Description
23-
!!
24-
!! This interface provides methods for computing the determinant of a matrix.
25-
!! Supported data types include real and complex.
26-
!!
27-
!!@note The provided functions are intended for square matrices.
28-
!!
29-
!!### Example
30-
!!
31-
!!```fortran
32-
!!
33-
!! real(sp) :: a(3,3), d
34-
!! type(linalg_state_type) :: state
35-
!! a = reshape([1, 2, 3, 4, 5, 6, 7, 8, 9], [3, 3])
36-
!!
37-
!! d = det(a,err=state)
38-
!! if (state%ok()) then
39-
!! print *, 'Success! det=',d
40-
!! else
41-
!! print *, state%print()
42-
!! endif
43-
!!
44-
!!```
45-
!!
46-
#:for rk,rt in RC_KINDS_TYPES
47-
#:if rk!="xdp"
48-
module procedure stdlib_linalg_${rt[0]}$${rk}$determinant
49-
module procedure stdlib_linalg_pure_${rt[0]}$${rk}$determinant
50-
#:endif
51-
#:endfor
52-
end interface det
53-
54-
interface operator(.det.)
55-
!!### Summary
56-
!! Pure operator interface for computing matrix determinant.
57-
!!
58-
!!### Description
59-
!!
60-
!! This pure operator interface provides a convenient way to compute the determinant of a matrix.
61-
!! Supported data types include real and complex.
62-
!!
63-
!!@note The provided functions are intended for square matrices.
64-
!!
65-
!!### Example
66-
!!
67-
!!```fortran
68-
!!
69-
!! real(sp) :: matrix(3,3), d
70-
!! matrix = reshape([1, 2, 3, 4, 5, 6, 7, 8, 9], [3, 3])
71-
!! d = .det.matrix
72-
!!
73-
!!```
74-
!
75-
#:for rk,rt in RC_KINDS_TYPES
76-
#:if rk!="xdp"
77-
module procedure stdlib_linalg_pure_${rt[0]}$${rk}$determinant
78-
#:endif
79-
#:endfor
80-
end interface operator(.det.)
81-
8215
contains
8316

8417
#:for rk,rt in RC_KINDS_TYPES
8518
#:if rk!="xdp"
86-
pure function stdlib_linalg_pure_${rt[0]}$${rk}$determinant(a) result(det)
19+
module function stdlib_linalg_pure_${rt[0]}$${rk}$determinant(a) result(det)
8720
!!### Summary
8821
!! Compute determinant of a real square matrix (pure interface).
8922
!!
@@ -180,7 +113,7 @@ module stdlib_linalg_determinant
180113

181114
end function stdlib_linalg_pure_${rt[0]}$${rk}$determinant
182115

183-
function stdlib_linalg_${rt[0]}$${rk}$determinant(a,overwrite_a,err) result(det)
116+
module function stdlib_linalg_${rt[0]}$${rk}$determinant(a,overwrite_a,err) result(det)
184117
!!### Summary
185118
!! Compute determinant of a square matrix (with error control).
186119
!!
@@ -299,4 +232,4 @@ module stdlib_linalg_determinant
299232
#:endif
300233
#:endfor
301234

302-
end module stdlib_linalg_determinant
235+
end submodule stdlib_linalg_determinant

0 commit comments

Comments
 (0)