Skip to content

Commit

Permalink
submodule version
Browse files Browse the repository at this point in the history
  • Loading branch information
perazz committed Apr 14, 2024
1 parent 4beab1f commit 15023e1
Show file tree
Hide file tree
Showing 2 changed files with 96 additions and 75 deletions.
95 changes: 92 additions & 3 deletions src/stdlib_linalg.fypp
Original file line number Diff line number Diff line change
@@ -1,18 +1,19 @@
#:include "common.fypp"
#:set RCI_KINDS_TYPES = REAL_KINDS_TYPES + CMPLX_KINDS_TYPES + INT_KINDS_TYPES
#:set RC_KINDS_TYPES = REAL_KINDS_TYPES + CMPLX_KINDS_TYPES
#:set RCI_KINDS_TYPES = RC_KINDS_TYPES + INT_KINDS_TYPES
module stdlib_linalg
!!Provides a support for various linear algebra procedures
!! ([Specification](../page/specs/stdlib_linalg.html))
use stdlib_kinds, only: sp, dp, xdp, qp, &
use stdlib_kinds, only: sp, dp, xdp, qp, lk, &
int8, int16, int32, int64
use stdlib_error, only: error_stop
use stdlib_optval, only: optval
use stdlib_linalg_state, only: linalg_state_type, linalg_error_handling
use stdlib_linalg_determinant, only: det
implicit none
private

public :: det
public :: operator(.det.)
public :: diag
public :: eye
public :: trace
Expand Down Expand Up @@ -220,6 +221,94 @@ module stdlib_linalg
#:endfor
end interface is_hessenberg


interface det
!!### Summary
!! Interface for computing matrix determinant.
!!
!!### Description
!!
!! This interface provides methods for computing the determinant of a matrix.
!! Supported data types include real and complex.
!!
!!@note The provided functions are intended for square matrices.
!!
!!### Example
!!
!!```fortran
!!
!! real(sp) :: a(3,3), d
!! type(linalg_state_type) :: state
!! a = reshape([1, 2, 3, 4, 5, 6, 7, 8, 9], [3, 3])
!!
!! d = det(a,err=state)
!! if (state%ok()) then
!! print *, 'Success! det=',d
!! else
!! print *, state%print()
!! endif
!!
!!```
!!
#:for rk,rt in RC_KINDS_TYPES
#:if rk!="xdp"
module procedure stdlib_linalg_${rt[0]}$${rk}$determinant
module procedure stdlib_linalg_pure_${rt[0]}$${rk}$determinant
#:endif
#:endfor
end interface det

interface operator(.det.)
!!### Summary
!! Pure operator interface for computing matrix determinant.
!!
!!### Description
!!
!! This pure operator interface provides a convenient way to compute the determinant of a matrix.
!! Supported data types include real and complex.
!!
!!@note The provided functions are intended for square matrices.
!!
!!### Example
!!
!!```fortran
!!
!! real(sp) :: matrix(3,3), d
!! matrix = reshape([1, 2, 3, 4, 5, 6, 7, 8, 9], [3, 3])
!! d = .det.matrix
!!
!!```
!
#:for rk,rt in RC_KINDS_TYPES
#:if rk!="xdp"
module procedure stdlib_linalg_pure_${rt[0]}$${rk}$determinant
#:endif
#:endfor
end interface operator(.det.)

interface
#:for rk,rt in RC_KINDS_TYPES
#:if rk!="xdp"
module function stdlib_linalg_${rt[0]}$${rk}$determinant(a,overwrite_a,err) result(det)
!> Input matrix a[m,n]
${rt}$, intent(inout), target :: a(:,:)
!> [optional] Can A data be overwritten and destroyed?
logical(lk), optional, intent(in) :: overwrite_a
!> State return flag.
type(linalg_state_type), intent(out) :: err
!> Matrix determinant
${rt}$ :: det
end function stdlib_linalg_${rt[0]}$${rk}$determinant
module function stdlib_linalg_pure_${rt[0]}$${rk}$determinant(a) result(det)
!> Input matrix a[m,n]
${rt}$, intent(in) :: a(:,:)
!> Matrix determinant
${rt}$ :: det
end function stdlib_linalg_pure_${rt[0]}$${rk}$determinant
#:endif
#:endfor
end interface

contains


Expand Down
76 changes: 4 additions & 72 deletions src/stdlib_linalg_determinant.fypp
Original file line number Diff line number Diff line change
@@ -1,89 +1,21 @@
#:include "common.fypp"
#:set RC_KINDS_TYPES = REAL_KINDS_TYPES + CMPLX_KINDS_TYPES
module stdlib_linalg_determinant
submodule (stdlib_linalg) stdlib_linalg_determinant
!! Determinant of a rectangular matrix
use stdlib_linalg_constants
use stdlib_linalg_lapack, only: getrf
use stdlib_linalg_state, only: linalg_state_type, linalg_error_handling, LINALG_ERROR, &
LINALG_INTERNAL_ERROR, LINALG_VALUE_ERROR
implicit none(type,external)
private

! Function interface
public :: det
public :: operator(.det.)

character(*), parameter :: this = 'determinant'

interface det
!!### Summary
!! Interface for computing matrix determinant.
!!
!!### Description
!!
!! This interface provides methods for computing the determinant of a matrix.
!! Supported data types include real and complex.
!!
!!@note The provided functions are intended for square matrices.
!!
!!### Example
!!
!!```fortran
!!
!! real(sp) :: a(3,3), d
!! type(linalg_state_type) :: state
!! a = reshape([1, 2, 3, 4, 5, 6, 7, 8, 9], [3, 3])
!!
!! d = det(a,err=state)
!! if (state%ok()) then
!! print *, 'Success! det=',d
!! else
!! print *, state%print()
!! endif
!!
!!```
!!
#:for rk,rt in RC_KINDS_TYPES
#:if rk!="xdp"
module procedure stdlib_linalg_${rt[0]}$${rk}$determinant
module procedure stdlib_linalg_pure_${rt[0]}$${rk}$determinant
#:endif
#:endfor
end interface det

interface operator(.det.)
!!### Summary
!! Pure operator interface for computing matrix determinant.
!!
!!### Description
!!
!! This pure operator interface provides a convenient way to compute the determinant of a matrix.
!! Supported data types include real and complex.
!!
!!@note The provided functions are intended for square matrices.
!!
!!### Example
!!
!!```fortran
!!
!! real(sp) :: matrix(3,3), d
!! matrix = reshape([1, 2, 3, 4, 5, 6, 7, 8, 9], [3, 3])
!! d = .det.matrix
!!
!!```
!
#:for rk,rt in RC_KINDS_TYPES
#:if rk!="xdp"
module procedure stdlib_linalg_pure_${rt[0]}$${rk}$determinant
#:endif
#:endfor
end interface operator(.det.)

contains

#:for rk,rt in RC_KINDS_TYPES
#:if rk!="xdp"
pure function stdlib_linalg_pure_${rt[0]}$${rk}$determinant(a) result(det)
module function stdlib_linalg_pure_${rt[0]}$${rk}$determinant(a) result(det)
!!### Summary
!! Compute determinant of a real square matrix (pure interface).
!!
Expand Down Expand Up @@ -180,7 +112,7 @@ module stdlib_linalg_determinant

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

function stdlib_linalg_${rt[0]}$${rk}$determinant(a,overwrite_a,err) result(det)
module function stdlib_linalg_${rt[0]}$${rk}$determinant(a,overwrite_a,err) result(det)
!!### Summary
!! Compute determinant of a square matrix (with error control).
!!
Expand Down Expand Up @@ -299,4 +231,4 @@ module stdlib_linalg_determinant
#:endif
#:endfor

end module stdlib_linalg_determinant
end submodule stdlib_linalg_determinant

0 comments on commit 15023e1

Please sign in to comment.