LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
ztbt02.f
Go to the documentation of this file.
1*> \brief \b ZTBT02
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8* Definition:
9* ===========
10*
11* SUBROUTINE ZTBT02( UPLO, TRANS, DIAG, N, KD, NRHS, AB, LDAB, X,
12* LDX, B, LDB, WORK, RWORK, RESID )
13*
14* .. Scalar Arguments ..
15* CHARACTER DIAG, TRANS, UPLO
16* INTEGER KD, LDAB, LDB, LDX, N, NRHS
17* DOUBLE PRECISION RESID
18* ..
19* .. Array Arguments ..
20* DOUBLE PRECISION RWORK( * )
21* COMPLEX*16 AB( LDAB, * ), B( LDB, * ), WORK( * ),
22* $ X( LDX, * )
23* ..
24*
25*
26*> \par Purpose:
27* =============
28*>
29*> \verbatim
30*>
31*> ZTBT02 computes the residual for the computed solution to a
32*> triangular system of linear equations op(A)*X = B, when A is a
33*> triangular band matrix. The test ratio is the maximum over
34*> norm(b - op(A)*x) / ( ||op(A)||_1 * norm(x) * EPS ),
35*> where op(A) = A, A**T, or A**H, b is the column of B, x is the
36*> solution vector, and EPS is the machine epsilon.
37*> \endverbatim
38*
39* Arguments:
40* ==========
41*
42*> \param[in] UPLO
43*> \verbatim
44*> UPLO is CHARACTER*1
45*> Specifies whether the matrix A is upper or lower triangular.
46*> = 'U': Upper triangular
47*> = 'L': Lower triangular
48*> \endverbatim
49*>
50*> \param[in] TRANS
51*> \verbatim
52*> TRANS is CHARACTER*1
53*> Specifies the operation applied to A.
54*> = 'N': A * X = B (No transpose)
55*> = 'T': A**T * X = B (Transpose)
56*> = 'C': A**H * X = B (Conjugate transpose)
57*> \endverbatim
58*>
59*> \param[in] DIAG
60*> \verbatim
61*> DIAG is CHARACTER*1
62*> Specifies whether or not the matrix A is unit triangular.
63*> = 'N': Non-unit triangular
64*> = 'U': Unit triangular
65*> \endverbatim
66*>
67*> \param[in] N
68*> \verbatim
69*> N is INTEGER
70*> The order of the matrix A. N >= 0.
71*> \endverbatim
72*>
73*> \param[in] KD
74*> \verbatim
75*> KD is INTEGER
76*> The number of superdiagonals or subdiagonals of the
77*> triangular band matrix A. KD >= 0.
78*> \endverbatim
79*>
80*> \param[in] NRHS
81*> \verbatim
82*> NRHS is INTEGER
83*> The number of right hand sides, i.e., the number of columns
84*> of the matrices X and B. NRHS >= 0.
85*> \endverbatim
86*>
87*> \param[in] AB
88*> \verbatim
89*> AB is COMPLEX*16 array, dimension (LDA,N)
90*> The upper or lower triangular band matrix A, stored in the
91*> first kd+1 rows of the array. The j-th column of A is stored
92*> in the j-th column of the array AB as follows:
93*> if UPLO = 'U', AB(kd+1+i-j,j) = A(i,j) for max(1,j-kd)<=i<=j;
94*> if UPLO = 'L', AB(1+i-j,j) = A(i,j) for j<=i<=min(n,j+kd).
95*> \endverbatim
96*>
97*> \param[in] LDAB
98*> \verbatim
99*> LDAB is INTEGER
100*> The leading dimension of the array AB. LDAB >= max(1,KD+1).
101*> \endverbatim
102*>
103*> \param[in] X
104*> \verbatim
105*> X is COMPLEX*16 array, dimension (LDX,NRHS)
106*> The computed solution vectors for the system of linear
107*> equations.
108*> \endverbatim
109*>
110*> \param[in] LDX
111*> \verbatim
112*> LDX is INTEGER
113*> The leading dimension of the array X. LDX >= max(1,N).
114*> \endverbatim
115*>
116*> \param[in] B
117*> \verbatim
118*> B is COMPLEX*16 array, dimension (LDB,NRHS)
119*> The right hand side vectors for the system of linear
120*> equations.
121*> \endverbatim
122*>
123*> \param[in] LDB
124*> \verbatim
125*> LDB is INTEGER
126*> The leading dimension of the array B. LDB >= max(1,N).
127*> \endverbatim
128*>
129*> \param[out] WORK
130*> \verbatim
131*> WORK is COMPLEX*16 array, dimension (N)
132*> \endverbatim
133*>
134*> \param[out] RWORK
135*> \verbatim
136*> RWORK is DOUBLE PRECISION array, dimension (N)
137*> \endverbatim
138*>
139*> \param[out] RESID
140*> \verbatim
141*> RESID is DOUBLE PRECISION
142*> The maximum over the number of right hand sides of
143*> norm(op(A)*x - b) / ( norm(op(A)) * norm(x) * EPS ).
144*> \endverbatim
145*
146* Authors:
147* ========
148*
149*> \author Univ. of Tennessee
150*> \author Univ. of California Berkeley
151*> \author Univ. of Colorado Denver
152*> \author NAG Ltd.
153*
154*> \ingroup complex16_lin
155*
156* =====================================================================
157 SUBROUTINE ztbt02( UPLO, TRANS, DIAG, N, KD, NRHS, AB, LDAB, X,
158 $ LDX, B, LDB, WORK, RWORK, RESID )
159*
160* -- LAPACK test routine --
161* -- LAPACK is a software package provided by Univ. of Tennessee, --
162* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
163*
164* .. Scalar Arguments ..
165 CHARACTER DIAG, TRANS, UPLO
166 INTEGER KD, LDAB, LDB, LDX, N, NRHS
167 DOUBLE PRECISION RESID
168* ..
169* .. Array Arguments ..
170 DOUBLE PRECISION RWORK( * )
171 COMPLEX*16 AB( LDAB, * ), B( LDB, * ), WORK( * ),
172 $ x( ldx, * )
173* ..
174*
175* =====================================================================
176*
177* .. Parameters ..
178 DOUBLE PRECISION ZERO, ONE
179 parameter( zero = 0.0d+0, one = 1.0d+0 )
180* ..
181* .. Local Scalars ..
182 INTEGER J
183 DOUBLE PRECISION ANORM, BNORM, EPS, XNORM
184* ..
185* .. External Functions ..
186 LOGICAL LSAME
187 DOUBLE PRECISION DLAMCH, DZASUM, ZLANTB
188 EXTERNAL lsame, dlamch, dzasum, zlantb
189* ..
190* .. External Subroutines ..
191 EXTERNAL zaxpy, zcopy, ztbmv
192* ..
193* .. Intrinsic Functions ..
194 INTRINSIC dcmplx, max
195* ..
196* .. Executable Statements ..
197*
198* Quick exit if N = 0 or NRHS = 0
199*
200 IF( n.LE.0 .OR. nrhs.LE.0 ) THEN
201 resid = zero
202 RETURN
203 END IF
204*
205* Compute the 1-norm of op(A).
206*
207 IF( lsame( trans, 'N' ) ) THEN
208 anorm = zlantb( '1', uplo, diag, n, kd, ab, ldab, rwork )
209 ELSE
210 anorm = zlantb( 'I', uplo, diag, n, kd, ab, ldab, rwork )
211 END IF
212*
213* Exit with RESID = 1/EPS if ANORM = 0.
214*
215 eps = dlamch( 'Epsilon' )
216 IF( anorm.LE.zero ) THEN
217 resid = one / eps
218 RETURN
219 END IF
220*
221* Compute the maximum over the number of right hand sides of
222* norm(op(A)*x - b) / ( norm(op(A)) * norm(x) * EPS ).
223*
224 resid = zero
225 DO 10 j = 1, nrhs
226 CALL zcopy( n, x( 1, j ), 1, work, 1 )
227 CALL ztbmv( uplo, trans, diag, n, kd, ab, ldab, work, 1 )
228 CALL zaxpy( n, dcmplx( -one ), b( 1, j ), 1, work, 1 )
229 bnorm = dzasum( n, work, 1 )
230 xnorm = dzasum( n, x( 1, j ), 1 )
231 IF( xnorm.LE.zero ) THEN
232 resid = one / eps
233 ELSE
234 resid = max( resid, ( ( bnorm / anorm ) / xnorm ) / eps )
235 END IF
236 10 CONTINUE
237*
238 RETURN
239*
240* End of ZTBT02
241*
242 END
subroutine zaxpy(n, za, zx, incx, zy, incy)
ZAXPY
Definition zaxpy.f:88
subroutine zcopy(n, zx, incx, zy, incy)
ZCOPY
Definition zcopy.f:81
subroutine ztbmv(uplo, trans, diag, n, k, a, lda, x, incx)
ZTBMV
Definition ztbmv.f:186
subroutine ztbt02(uplo, trans, diag, n, kd, nrhs, ab, ldab, x, ldx, b, ldb, work, rwork, resid)
ZTBT02
Definition ztbt02.f:159