Skip to content

Commit 29a6e04

Browse files
committed
Improved error handling.
1 parent 9d37af7 commit 29a6e04

File tree

2 files changed

+112
-7
lines changed

2 files changed

+112
-7
lines changed

src/stdlib_specialmatrices.fypp

Lines changed: 27 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -10,6 +10,8 @@ module stdlib_specialmatrices
1010
!! applications. ([Specifications](../page/specs/stdlib_specialmatrices.html))
1111
use stdlib_linalg_constants
1212
use stdlib_constants
13+
use stdlib_linalg_state, only: linalg_state_type, linalg_error_handling, LINALG_ERROR, &
14+
LINALG_INTERNAL_ERROR, LINALG_VALUE_ERROR
1315
implicit none
1416
private
1517
public :: tridiagonal
@@ -80,7 +82,7 @@ module stdlib_specialmatrices
8082
!! A = Tridiagonal(a, b, c, n)
8183
!! ```
8284
#:for k1, t1, s1 in (KINDS_TYPES)
83-
pure module function initialize_tridiagonal_${s1}$(dl, dv, du) result(A)
85+
pure module function initialize_tridiagonal_pure_${s1}$(dl, dv, du) result(A)
8486
!! Construct a `tridiagonal` matrix from the rank-1 arrays
8587
!! `dl`, `dv` and `du`.
8688
${t1}$, intent(in) :: dl(:), dv(:), du(:)
@@ -89,7 +91,7 @@ module stdlib_specialmatrices
8991
!! Corresponding Tridiagonal matrix.
9092
end function
9193

92-
pure module function initialize_constant_tridiagonal_${s1}$(dl, dv, du, n) result(A)
94+
pure module function initialize_constant_tridiagonal_pure_${s1}$(dl, dv, du, n) result(A)
9395
!! Construct a `tridiagonal` matrix with constant elements.
9496
${t1}$, intent(in) :: dl, dv, du
9597
!! Tridiagonal matrix elements.
@@ -98,6 +100,29 @@ module stdlib_specialmatrices
98100
type(tridiagonal_${s1}$_type) :: A
99101
!! Corresponding Tridiagonal matrix.
100102
end function
103+
104+
module function initialize_tridiagonal_impure_${s1}$(dl, dv, du, err) result(A)
105+
!! Construct a `tridiagonal` matrix from the rank-1 arrays
106+
!! `dl`, `dv` and `du`.
107+
${t1}$, intent(in) :: dl(:), dv(:), du(:)
108+
!! Tridiagonal matrix elements.
109+
type(linalg_state_type), intent(out) :: err
110+
!! Error handling.
111+
type(tridiagonal_${s1}$_type) :: A
112+
!! Corresponding Tridiagonal matrix.
113+
end function
114+
115+
module function initialize_constant_tridiagonal_impure_${s1}$(dl, dv, du, n, err) result(A)
116+
!! Construct a `tridiagonal` matrix with constant elements.
117+
${t1}$, intent(in) :: dl, dv, du
118+
!! Tridiagonal matrix elements.
119+
integer(ilp), intent(in) :: n
120+
!! Matrix dimension.
121+
type(linalg_state_type), intent(out) :: err
122+
!! Error handling.
123+
type(tridiagonal_${s1}$_type) :: A
124+
!! Corresponding Tridiagonal matrix.
125+
end function
101126
#:endfor
102127
end interface
103128

src/stdlib_specialmatrices_tridiagonal.fypp

Lines changed: 85 additions & 5 deletions
Original file line numberDiff line numberDiff line change
@@ -5,6 +5,8 @@
55
#:set KINDS_TYPES = R_KINDS_TYPES+C_KINDS_TYPES
66
submodule (stdlib_specialmatrices) tridiagonal_matrices
77
use stdlib_linalg_lapack, only: lagtm
8+
9+
character(len=*), parameter :: this = "tridiagonal matrices"
810
contains
911

1012
!--------------------------------
@@ -14,42 +16,120 @@ submodule (stdlib_specialmatrices) tridiagonal_matrices
1416
!--------------------------------
1517

1618
#:for k1, t1, s1 in (KINDS_TYPES)
17-
pure module function initialize_tridiagonal_${s1}$(dl, dv, du) result(A)
19+
pure module function initialize_tridiagonal_pure_${s1}$(dl, dv, du) result(A)
20+
!! Construct a `tridiagonal` matrix from the rank-1 arrays
21+
!! `dl`, `dv` and `du`.
22+
${t1}$, intent(in) :: dl(:), dv(:), du(:)
23+
!! tridiagonal matrix elements.
24+
type(tridiagonal_${s1}$_type) :: A
25+
!! Corresponding tridiagonal matrix.
26+
27+
! Internal variables.
28+
integer(ilp) :: n
29+
type(linalg_state_type) :: err0
30+
31+
! Sanity check.
32+
n = size(dv, kind=ilp)
33+
if (n <= 0) then
34+
err0 = linalg_state_type(this, LINALG_VALUE_ERROR, "Matrix size needs to be positive, n = ", n, ".")
35+
call linalg_error_handling(err0)
36+
endif
37+
if (size(dl, kind=ilp) /= n-1) then
38+
err0 = linalg_state_type(this, LINALG_VALUE_ERROR, "Vector dl does not have the correct length.")
39+
call linalg_error_handling(err0)
40+
endif
41+
if (size(du, kind=ilp) /= n-1) then
42+
err0 = linalg_state_type(this, LINALG_VALUE_ERROR, "Vector du does not have the correct length.")
43+
call linalg_error_handling(err0)
44+
endif
45+
46+
! Description of the matrix.
47+
A%n = n
48+
! Matrix elements.
49+
A%dl = dl ; A%dv = dv ; A%du = du
50+
end function
51+
52+
pure module function initialize_constant_tridiagonal_pure_${s1}$(dl, dv, du, n) result(A)
53+
!! Construct a `tridiagonal` matrix with constant elements.
54+
${t1}$, intent(in) :: dl, dv, du
55+
!! tridiagonal matrix elements.
56+
integer(ilp), intent(in) :: n
57+
!! Matrix dimension.
58+
type(tridiagonal_${s1}$_type) :: A
59+
!! Corresponding tridiagonal matrix.
60+
61+
! Internal variables.
62+
integer(ilp) :: i
63+
type(linalg_state_type) :: err0
64+
65+
! Description of the matrix.
66+
A%n = n
67+
if (n <= 0) then
68+
err0 = linalg_state_type(this, LINALG_VALUE_ERROR, "Matrix size needs to be positive, n = ", n, ".")
69+
call linalg_error_handling(err0)
70+
endif
71+
! Matrix elements.
72+
A%dl = [(dl, i = 1, n-1)]
73+
A%dv = [(dv, i = 1, n)]
74+
A%du = [(du, i = 1, n-1)]
75+
end function
76+
77+
module function initialize_tridiagonal_impure_${s1}$(dl, dv, du, err) result(A)
1878
!! Construct a `tridiagonal` matrix from the rank-1 arrays
1979
!! `dl`, `dv` and `du`.
2080
${t1}$, intent(in) :: dl(:), dv(:), du(:)
2181
!! tridiagonal matrix elements.
82+
type(linalg_state_type), intent(out) :: err
83+
!! Error handling.
2284
type(tridiagonal_${s1}$_type) :: A
2385
!! Corresponding tridiagonal matrix.
2486

2587
! Internal variables.
2688
integer(ilp) :: n
89+
type(linalg_state_type) :: err0
2790

2891
! Sanity check.
29-
n = size(dv)
30-
if (size(dl) /= n-1) error stop "Vector dl does not have the correct length."
31-
if (size(du) /= n-1) error stop "Vector du does not have the correct length."
92+
n = size(dv, kind=ilp)
93+
if (n <= 0) then
94+
err0 = linalg_state_type(this, LINALG_VALUE_ERROR, "Matrix size needs to be positive, n = ", n, ".")
95+
call linalg_error_handling(err0, err)
96+
endif
97+
if (size(dl, kind=ilp) /= n-1) then
98+
err0 = linalg_state_type(this, LINALG_VALUE_ERROR, "Vector dl does not have the correct length.")
99+
call linalg_error_handling(err0, err)
100+
endif
101+
if (size(du, kind=ilp) /= n-1) then
102+
err0 = linalg_state_type(this, LINALG_VALUE_ERROR, "Vector du does not have the correct length.")
103+
call linalg_error_handling(err0, err)
104+
endif
32105

33106
! Description of the matrix.
34107
A%n = n
35108
! Matrix elements.
36109
A%dl = dl ; A%dv = dv ; A%du = du
37110
end function
38111

39-
pure module function initialize_constant_tridiagonal_${s1}$(dl, dv, du, n) result(A)
112+
module function initialize_constant_tridiagonal_impure_${s1}$(dl, dv, du, n, err) result(A)
40113
!! Construct a `tridiagonal` matrix with constant elements.
41114
${t1}$, intent(in) :: dl, dv, du
42115
!! tridiagonal matrix elements.
43116
integer(ilp), intent(in) :: n
44117
!! Matrix dimension.
118+
type(linalg_state_type), intent(out) :: err
119+
!! Error handling
45120
type(tridiagonal_${s1}$_type) :: A
46121
!! Corresponding tridiagonal matrix.
47122

48123
! Internal variables.
49124
integer(ilp) :: i
125+
type(linalg_state_type) :: err0
50126

51127
! Description of the matrix.
52128
A%n = n
129+
if (n <= 0) then
130+
err0 = linalg_state_type(this, LINALG_VALUE_ERROR, "Matrix size needs to be positive, n = ", n, ".")
131+
call linalg_error_handling(err0, err)
132+
endif
53133
! Matrix elements.
54134
A%dl = [(dl, i = 1, n-1)]
55135
A%dv = [(dv, i = 1, n)]

0 commit comments

Comments
 (0)