-
Notifications
You must be signed in to change notification settings - Fork 191
linalg: determinant #798
New issue
Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.
By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.
Already on GitHub? Sign in to your account
linalg: determinant #798
Changes from all commits
4a4b899
d54cfa3
c987fa1
b9a3b0e
4a29219
3d05cc3
4beab1f
15023e1
5923a54
31cdc84
45a606f
cacb585
5c16ff8
ab030c5
13bd98a
7bf7141
504d90d
e80b508
c6076ea
eaebe5c
5d52d48
cff995d
2fe2428
45447a8
3ee20ad
1e50115
fe933ad
59dae89
3c72b06
File filter
Filter by extension
Conversations
Jump to
Diff view
Diff view
There are no files selected for viewing
Original file line number | Diff line number | Diff line change |
---|---|---|
@@ -0,0 +1,14 @@ | ||
program example_determinant | ||
use stdlib_kinds, only: dp | ||
use stdlib_linalg, only: det, linalg_state_type | ||
implicit none | ||
type(linalg_state_type) :: err | ||
|
||
real(dp) :: d | ||
|
||
! Compute determinate of a real matrix | ||
d = det(reshape([real(dp)::1,2,3,4],[2,2])) | ||
|
||
print *, d ! a*d-b*c = -2.0 | ||
|
||
end program example_determinant |
Original file line number | Diff line number | Diff line change |
---|---|---|
@@ -0,0 +1,13 @@ | ||
program example_determinant2 | ||
use stdlib_kinds, only: dp | ||
use stdlib_linalg, only: operator(.det.) | ||
implicit none | ||
|
||
real(dp) :: d | ||
|
||
! Compute determinate of a real matrix | ||
d = .det.reshape([real(dp)::1,2,3,4],[2,2]) | ||
|
||
print *, d ! a*d-b*c = -2.0 | ||
|
||
end program example_determinant2 |
Original file line number | Diff line number | Diff line change |
---|---|---|
@@ -1,15 +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 | ||
implicit none | ||
private | ||
|
||
public :: det | ||
public :: operator(.det.) | ||
public :: diag | ||
public :: eye | ||
public :: trace | ||
|
@@ -23,6 +27,9 @@ module stdlib_linalg | |
public :: is_hermitian | ||
public :: is_triangular | ||
public :: is_hessenberg | ||
|
||
! Export linalg error handling | ||
public :: linalg_state_type, linalg_error_handling | ||
Comment on lines
+30
to
+32
There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. This let me think that we should review the structure of the specs. Currently we have 1 page per module. However, if these become available through |
||
|
||
interface diag | ||
!! version: experimental | ||
|
@@ -214,6 +221,109 @@ module stdlib_linalg | |
#:endfor | ||
end interface is_hessenberg | ||
|
||
|
||
interface det | ||
!! version: experimental | ||
!! | ||
!! Computes the determinant of a square matrix | ||
!! ([Specification](../page/specs/stdlib_linalg.html#det-computes-the-determinant-of-a-square-matrix)) | ||
!! | ||
!!### Summary | ||
jvdp1 marked this conversation as resolved.
Show resolved
Hide resolved
|
||
!! 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 only. | ||
!!@note BLAS/LAPACK backends do not currently support extended precision (``xdp``). | ||
!! | ||
!!### Example | ||
!! | ||
!!```fortran | ||
perazz marked this conversation as resolved.
Show resolved
Hide resolved
|
||
!! | ||
!! 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.) | ||
!! version: experimental | ||
!! | ||
!! Determinant operator of a square matrix | ||
!! ([Specification](../page/specs/stdlib_linalg.html#det-determinant-operator-of-a-square-matrix)) | ||
!! | ||
!!### Summary | ||
jvdp1 marked this conversation as resolved.
Show resolved
Hide resolved
|
||
!! 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. | ||
!!@note BLAS/LAPACK backends do not currently support extended precision (``xdp``). | ||
!! | ||
!!### 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 | ||
jvdp1 marked this conversation as resolved.
Show resolved
Hide resolved
|
||
!> Matrix determinant | ||
${rt}$ :: det | ||
end function stdlib_linalg_${rt[0]}$${rk}$determinant | ||
pure 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 | ||
|
||
|
||
|
Uh oh!
There was an error while loading. Please reload this page.