Skip to content

Commit 07c244d

Browse files
committed
Moved handle_gesdd_info and signature.
1 parent d0181ca commit 07c244d

File tree

2 files changed

+36
-35
lines changed

2 files changed

+36
-35
lines changed

src/lapack/stdlib_linalg_lapack_aux.fypp

Lines changed: 33 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -4,7 +4,7 @@ module stdlib_linalg_lapack_aux
44
use stdlib_linalg_constants
55
use stdlib_linalg_blas
66
use stdlib_linalg_state, only: linalg_state_type, linalg_error_handling, LINALG_ERROR, &
7-
LINALG_INTERNAL_ERROR, LINALG_VALUE_ERROR
7+
LINALG_INTERNAL_ERROR, LINALG_VALUE_ERROR, LINALG_SUCCESS
88
use ieee_arithmetic, only: ieee_support_inf, ieee_support_nan
99
implicit none
1010
private
@@ -42,6 +42,7 @@ module stdlib_linalg_lapack_aux
4242
#:endfor
4343
public :: handle_potrf_info
4444
public :: handle_getri_info
45+
public :: handle_gesdd_info
4546

4647
! SELCTG is a LOGICAL FUNCTION of three DOUBLE PRECISION arguments
4748
! used to select eigenvalues to sort to the top left of the Schur form.
@@ -1331,4 +1332,35 @@ module stdlib_linalg_lapack_aux
13311332
end select
13321333
end subroutine handle_getri_info
13331334

1335+
elemental subroutine handle_gesdd_info(this,err,info,m,n)
1336+
character(len=*), intent(in) :: this
1337+
!> Error handler
1338+
type(linalg_state_type), intent(inout) :: err
1339+
!> GESDD return flag
1340+
integer(ilp), intent(in) :: info
1341+
!> Input matrix size
1342+
integer(ilp), intent(in) :: m,n
1343+
1344+
select case (info)
1345+
case (0)
1346+
! Success!
1347+
err%state = LINALG_SUCCESS
1348+
case (-1)
1349+
err = linalg_state_type(this,LINALG_INTERNAL_ERROR,'Invalid task ID on input to GESDD.')
1350+
case (-5,-3:-2)
1351+
err = linalg_state_type(this,LINALG_VALUE_ERROR,'invalid matrix size: a=',[m,n])
1352+
case (-8)
1353+
err = linalg_state_type(this,LINALG_VALUE_ERROR,'invalid matrix U size, with a=',[m,n])
1354+
case (-10)
1355+
err = linalg_state_type(this,LINALG_VALUE_ERROR,'invalid matrix V size, with a=',[m,n])
1356+
case (-4)
1357+
err = linalg_state_type(this,LINALG_VALUE_ERROR,'A contains invalid/NaN values.')
1358+
case (1:)
1359+
err = linalg_state_type(this,LINALG_ERROR,'SVD computation did not converge.')
1360+
case default
1361+
err = linalg_state_type(this,LINALG_INTERNAL_ERROR,'Unknown error returned by GESDD.')
1362+
end select
1363+
1364+
end subroutine handle_gesdd_info
1365+
13341366
end module stdlib_linalg_lapack_aux

src/stdlib_linalg_svd.fypp

Lines changed: 3 additions & 34 deletions
Original file line numberDiff line numberDiff line change
@@ -4,6 +4,7 @@ submodule(stdlib_linalg) stdlib_linalg_svd
44
!! Singular-Value Decomposition
55
use stdlib_linalg_constants
66
use stdlib_linalg_lapack, only: gesdd
7+
use stdlib_linalg_lapack_aux, only: handle_gesdd_info
78
use stdlib_linalg_state, only: linalg_state_type, linalg_error_handling, LINALG_ERROR, &
89
LINALG_INTERNAL_ERROR, LINALG_VALUE_ERROR, LINALG_SUCCESS
910
implicit none
@@ -26,38 +27,6 @@ submodule(stdlib_linalg) stdlib_linalg_svd
2627

2728
contains
2829

29-
!> Process GESDD output flag
30-
elemental subroutine handle_gesdd_info(err,info,m,n)
31-
!> Error handler
32-
type(linalg_state_type), intent(inout) :: err
33-
!> GESDD return flag
34-
integer(ilp), intent(in) :: info
35-
!> Input matrix size
36-
integer(ilp), intent(in) :: m,n
37-
38-
select case (info)
39-
case (0)
40-
! Success!
41-
err%state = LINALG_SUCCESS
42-
case (-1)
43-
err = linalg_state_type(this,LINALG_INTERNAL_ERROR,'Invalid task ID on input to GESDD.')
44-
case (-5,-3:-2)
45-
err = linalg_state_type(this,LINALG_VALUE_ERROR,'invalid matrix size: a=',[m,n])
46-
case (-8)
47-
err = linalg_state_type(this,LINALG_VALUE_ERROR,'invalid matrix U size, with a=',[m,n])
48-
case (-10)
49-
err = linalg_state_type(this,LINALG_VALUE_ERROR,'invalid matrix V size, with a=',[m,n])
50-
case (-4)
51-
err = linalg_state_type(this,LINALG_VALUE_ERROR,'A contains invalid/NaN values.')
52-
case (1:)
53-
err = linalg_state_type(this,LINALG_ERROR,'SVD computation did not converge.')
54-
case default
55-
err = linalg_state_type(this,LINALG_INTERNAL_ERROR,'Unknown error returned by GESDD.')
56-
end select
57-
58-
end subroutine handle_gesdd_info
59-
60-
6130
#:for rk,rt,ri in RC_KINDS_TYPES
6231

6332
!> Singular values of matrix A
@@ -265,7 +234,7 @@ submodule(stdlib_linalg) stdlib_linalg_svd
265234
lwork = -1_ilp
266235
call gesdd(task,m,n,amat,lda,s,umat,ldu,vtmat,ldvt,&
267236
work_dummy,lwork,#{if rt.startswith('complex')}#rwork,#{endif}#iwork,info)
268-
call handle_gesdd_info(err0,info,m,n)
237+
call handle_gesdd_info(this,err0,info,m,n)
269238

270239
! Compute SVD
271240
if (info==0) then
@@ -281,7 +250,7 @@ submodule(stdlib_linalg) stdlib_linalg_svd
281250
!> Compute SVD
282251
call gesdd(task,m,n,amat,lda,s,umat,ldu,vtmat,ldvt,&
283252
work,lwork,#{if rt.startswith('complex')}#rwork,#{endif}#iwork,info)
284-
call handle_gesdd_info(err0,info,m,n)
253+
call handle_gesdd_info(this,err0,info,m,n)
285254

286255
endif
287256

0 commit comments

Comments
 (0)