- Notifications
You must be signed in to change notification settings - Fork 198
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 1 commit
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
pure interface - Loading branch information
There are no files selected for viewing
| Original file line number | Diff line number | Diff line change |
|---|---|---|
| @@ -1,14 +1,15 @@ | ||
| #:include "common.fypp" | ||
| #:set RC_KINDS_TYPES = REAL_KINDS_TYPES + CMPLX_KINDS_TYPES | ||
| !> Determinant of a rectangular matrix | ||
| module stdlib_linalg_determinant | ||
| ||
| use stdlib_linalg_constants | ||
| use stdlib_linalg_blas | ||
| use stdlib_linalg_lapack | ||
| use stdlib_linalg_state | ||
| use iso_fortran_env,only:real32,real64,real128,int8,int16,int32,int64,stderr => error_unit | ||
| 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 | ||
| | ||
| !> Determinant of a rectangular matrix | ||
| !> Function interface | ||
| public :: det | ||
| | ||
| character(*), parameter :: this = 'determinant' | ||
| | @@ -18,28 +19,106 @@ module stdlib_linalg_determinant | |
| ! IMSL: DET(a) | ||
| | ||
| interface det | ||
| #:for rk,rt,ri in ALL_KINDS_TYPES | ||
| module procedure stdlib_linalg_${ri}$determinant | ||
| #:for rk,rt in RC_KINDS_TYPES | ||
| module procedure stdlib_linalg_${rt[0]}$${rk}$determinant | ||
| module procedure stdlib_linalg_pure_${rt[0]}$${rk}$determinant | ||
| #:endfor | ||
| end interface det | ||
| | ||
| | ||
| contains | ||
| | ||
| #:for rk,rt,ri in ALL_KINDS_TYPES | ||
| ! Compute determinant of a square matrix A | ||
| function stdlib_linalg_${ri}$determinant(a,overwrite_a,err) result(det) | ||
| #:for rk,rt in RC_KINDS_TYPES | ||
| ! Compute determinant of a square matrix A: pure interface | ||
| function stdlib_linalg_pure_${rt[0]}$${rk}$determinant(a) result(det) | ||
| !> Input matrix a[m,n] | ||
| ${rt}$, intent(in), target :: a(:,:) | ||
| !> Result: matrix determinant | ||
| ${rt}$ :: det | ||
| | ||
| !> Local variables | ||
| type(linalg_state_type) :: err0 | ||
| integer(ilp) :: m,n,info,perm,k | ||
| integer(ilp), allocatable :: ipiv(:) | ||
| ${rt}$, allocatable :: amat(:,:) | ||
| | ||
| !> Matrix determinant size | ||
| m = size(a,1,kind=ilp) | ||
| n = size(a,2,kind=ilp) | ||
| | ||
| if (m/=n .or. .not.min(m,n)>=0) then | ||
| err0 = linalg_state_type(this,LINALG_VALUE_ERROR,'invalid or non-square matrix: a=[',m,',',n,']') | ||
| det = 0.0_${rk}$ | ||
| goto 1 | ||
perazz marked this conversation as resolved. Outdated Show resolved Hide resolved | ||
| end if | ||
| | ||
| select case (m) | ||
| case (0) | ||
| | ||
| ! Empty array has determinant 1 because math | ||
| det = 1.0_${rk}$ | ||
| | ||
| case (1) | ||
| | ||
| ! Scalar input | ||
| det = a(1,1) | ||
| | ||
| case default | ||
| | ||
| ! Find determinant from LU decomposition | ||
| | ||
| ! Initialize a matrix temporary | ||
| allocate(amat(m,n),source=a) | ||
| | ||
| ! Pivot indices | ||
| allocate(ipiv(n)) | ||
| | ||
| ! Compute determinant from LU factorization, then calculate the | ||
| ! product of all diagonal entries of the U factor. | ||
| call getrf(m,n,amat,m,ipiv,info) | ||
| | ||
| select case (info) | ||
| case (0) | ||
| ! Success: compute determinant | ||
| | ||
| ! Start with real 1.0 | ||
| det = 1.0_${rk}$ | ||
| perm = 0 | ||
| do k=1,n | ||
| if (ipiv(k)/=k) perm = perm+1 | ||
| det = det*amat(k,k) | ||
| end do | ||
| if (mod(perm,2)/=0) det = -det | ||
| | ||
| case (:-1) | ||
| err0 = linalg_state_type(this,LINALG_ERROR,'invalid matrix size a=[',m,',',n,']') | ||
| case (1:) | ||
| err0 = linalg_state_type(this,LINALG_ERROR,'singular matrix') | ||
| case default | ||
| err0 = linalg_state_type(this,LINALG_INTERNAL_ERROR,'catastrophic error') | ||
| end select | ||
| | ||
| deallocate(amat) | ||
| | ||
| end select | ||
| | ||
| ! Process output and return | ||
| 1 call linalg_error_handling(err0) | ||
| | ||
| end function stdlib_linalg_pure_${rt[0]}$${rk}$determinant | ||
| | ||
| ! Compute determinant of a square matrix A, with error control | ||
| 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 | ||
| !> [optional] state return flag. On error if not requested, the code will stop | ||
| type(linalg_state), optional, intent(out) :: err | ||
| type(linalg_state_type), intent(out) :: err | ||
| !> Result: matrix determinant | ||
| ${rt}$ :: det | ||
| | ||
| !> Local variables | ||
| type(linalg_state) :: err0 | ||
| type(linalg_state_type) :: err0 | ||
| integer(ilp) :: m,n,info,perm,k | ||
| integer(ilp), allocatable :: ipiv(:) | ||
| logical(lk) :: copy_a | ||
| | @@ -50,7 +129,7 @@ module stdlib_linalg_determinant | |
| n = size(a,2,kind=ilp) | ||
| | ||
| if (m/=n .or. .not.min(m,n)>=0) then | ||
| err0 = linalg_state(this,LINALG_VALUE_ERROR,'invalid or non-square matrix: a=[',m,',',n,']') | ||
| err0 = linalg_state_type(this,LINALG_VALUE_ERROR,'invalid or non-square matrix: a=[',m,',',n,']') | ||
| det = 0.0_${rk}$ | ||
| goto 1 | ||
perazz marked this conversation as resolved. Outdated Show resolved Hide resolved | ||
| end if | ||
| | @@ -64,11 +143,13 @@ module stdlib_linalg_determinant | |
| | ||
| select case (m) | ||
| case (0) | ||
| | ||
| ! Empty array has determinant 1 because math | ||
| det = 1.0_${rk}$ | ||
| | ||
| case (1) | ||
| ! Scalar | ||
| | ||
| ! Scalar input | ||
| det = a(1,1) | ||
| | ||
| case default | ||
| | @@ -85,8 +166,8 @@ module stdlib_linalg_determinant | |
| ! Pivot indices | ||
| allocate(ipiv(n)) | ||
| | ||
| ! Compute determinant from LU factorization, then calculate the product of | ||
| ! all diagonal entries of the U factor. | ||
| ! Compute determinant from LU factorization, then calculate the | ||
| ! product of all diagonal entries of the U factor. | ||
| call getrf(m,n,amat,m,ipiv,info) | ||
| | ||
| select case (info) | ||
| | @@ -103,11 +184,11 @@ module stdlib_linalg_determinant | |
| if (mod(perm,2)/=0) det = -det | ||
| | ||
| case (:-1) | ||
| err0 = linalg_state(this,LINALG_ERROR,'invalid matrix size a=[',m,',',n,']') | ||
| err0 = linalg_state_type(this,LINALG_ERROR,'invalid matrix size a=[',m,',',n,']') | ||
| case (1:) | ||
| err0 = linalg_state(this,LINALG_ERROR,'singular matrix') | ||
| err0 = linalg_state_type(this,LINALG_ERROR,'singular matrix') | ||
| case default | ||
| err0 = linalg_state(this,LINALG_INTERNAL_ERROR,'catastrophic error') | ||
| err0 = linalg_state_type(this,LINALG_INTERNAL_ERROR,'catastrophic error') | ||
| end select | ||
| | ||
| if (.not.copy_a) deallocate(amat) | ||
| | @@ -117,7 +198,7 @@ module stdlib_linalg_determinant | |
| ! Process output and return | ||
| 1 call linalg_error_handling(err0,err) | ||
| | ||
| end function stdlib_linalg_${ri}$determinant | ||
| end function stdlib_linalg_${rt[0]}$${rk}$determinant | ||
| | ||
| #:endfor | ||
| | ||
| | ||
There was a problem hiding this comment.
Choose a reason for hiding this comment
The 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
stdlib_linalg, the specs should be modified accordingly (probably in another PR; I can start to work on a proposition when this is ready).