LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches

◆ zgbt01()

subroutine zgbt01 ( integer  m,
integer  n,
integer  kl,
integer  ku,
complex*16, dimension( lda, * )  a,
integer  lda,
complex*16, dimension( ldafac, * )  afac,
integer  ldafac,
integer, dimension( * )  ipiv,
complex*16, dimension( * )  work,
double precision  resid 
)

ZGBT01

Purpose:
 ZGBT01 reconstructs a band matrix A from its L*U factorization and
 computes the residual:
    norm(L*U - A) / ( N * norm(A) * EPS ),
 where EPS is the machine epsilon.

 The expression L*U - A is computed one column at a time, so A and
 AFAC are not modified.
Parameters
[in]M
          M is INTEGER
          The number of rows of the matrix A.  M >= 0.
[in]N
          N is INTEGER
          The number of columns of the matrix A.  N >= 0.
[in]KL
          KL is INTEGER
          The number of subdiagonals within the band of A.  KL >= 0.
[in]KU
          KU is INTEGER
          The number of superdiagonals within the band of A.  KU >= 0.
[in,out]A
          A is COMPLEX*16 array, dimension (LDA,N)
          The original matrix A in band storage, stored in rows 1 to
          KL+KU+1.
[in]LDA
          LDA is INTEGER.
          The leading dimension of the array A.  LDA >= max(1,KL+KU+1).
[in]AFAC
          AFAC is COMPLEX*16 array, dimension (LDAFAC,N)
          The factored form of the matrix A.  AFAC contains the banded
          factors L and U from the L*U factorization, as computed by
          ZGBTRF.  U is stored as an upper triangular band matrix with
          KL+KU superdiagonals in rows 1 to KL+KU+1, and the
          multipliers used during the factorization are stored in rows
          KL+KU+2 to 2*KL+KU+1.  See ZGBTRF for further details.
[in]LDAFAC
          LDAFAC is INTEGER
          The leading dimension of the array AFAC.
          LDAFAC >= max(1,2*KL*KU+1).
[in]IPIV
          IPIV is INTEGER array, dimension (min(M,N))
          The pivot indices from ZGBTRF.
[out]WORK
          WORK is COMPLEX*16 array, dimension (2*KL+KU+1)
[out]RESID
          RESID is DOUBLE PRECISION
          norm(L*U - A) / ( N * norm(A) * EPS )
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.

Definition at line 124 of file zgbt01.f.

126*
127* -- LAPACK test routine --
128* -- LAPACK is a software package provided by Univ. of Tennessee, --
129* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
130*
131* .. Scalar Arguments ..
132 INTEGER KL, KU, LDA, LDAFAC, M, N
133 DOUBLE PRECISION RESID
134* ..
135* .. Array Arguments ..
136 INTEGER IPIV( * )
137 COMPLEX*16 A( LDA, * ), AFAC( LDAFAC, * ), WORK( * )
138* ..
139*
140* =====================================================================
141*
142* .. Parameters ..
143 DOUBLE PRECISION ZERO, ONE
144 parameter( zero = 0.0d+0, one = 1.0d+0 )
145* ..
146* .. Local Scalars ..
147 INTEGER I, I1, I2, IL, IP, IW, J, JL, JU, JUA, KD, LENJ
148 DOUBLE PRECISION ANORM, EPS
149 COMPLEX*16 T
150* ..
151* .. External Functions ..
152 DOUBLE PRECISION DLAMCH, DZASUM
153 EXTERNAL dlamch, dzasum
154* ..
155* .. External Subroutines ..
156 EXTERNAL zaxpy, zcopy
157* ..
158* .. Intrinsic Functions ..
159 INTRINSIC dble, dcmplx, max, min
160* ..
161* .. Executable Statements ..
162*
163* Quick exit if M = 0 or N = 0.
164*
165 resid = zero
166 IF( m.LE.0 .OR. n.LE.0 )
167 $ RETURN
168*
169* Determine EPS and the norm of A.
170*
171 eps = dlamch( 'Epsilon' )
172 kd = ku + 1
173 anorm = zero
174 DO 10 j = 1, n
175 i1 = max( kd+1-j, 1 )
176 i2 = min( kd+m-j, kl+kd )
177 IF( i2.GE.i1 )
178 $ anorm = max( anorm, dzasum( i2-i1+1, a( i1, j ), 1 ) )
179 10 CONTINUE
180*
181* Compute one column at a time of L*U - A.
182*
183 kd = kl + ku + 1
184 DO 40 j = 1, n
185*
186* Copy the J-th column of U to WORK.
187*
188 ju = min( kl+ku, j-1 )
189 jl = min( kl, m-j )
190 lenj = min( m, j ) - j + ju + 1
191 IF( lenj.GT.0 ) THEN
192 CALL zcopy( lenj, afac( kd-ju, j ), 1, work, 1 )
193 DO 20 i = lenj + 1, ju + jl + 1
194 work( i ) = zero
195 20 CONTINUE
196*
197* Multiply by the unit lower triangular matrix L. Note that L
198* is stored as a product of transformations and permutations.
199*
200 DO 30 i = min( m-1, j ), j - ju, -1
201 il = min( kl, m-i )
202 IF( il.GT.0 ) THEN
203 iw = i - j + ju + 1
204 t = work( iw )
205 CALL zaxpy( il, t, afac( kd+1, i ), 1, work( iw+1 ),
206 $ 1 )
207 ip = ipiv( i )
208 IF( i.NE.ip ) THEN
209 ip = ip - j + ju + 1
210 work( iw ) = work( ip )
211 work( ip ) = t
212 END IF
213 END IF
214 30 CONTINUE
215*
216* Subtract the corresponding column of A.
217*
218 jua = min( ju, ku )
219 IF( jua+jl+1.GT.0 )
220 $ CALL zaxpy( jua+jl+1, -dcmplx( one ), a( ku+1-jua, j ),
221 $ 1, work( ju+1-jua ), 1 )
222*
223* Compute the 1-norm of the column.
224*
225 resid = max( resid, dzasum( ju+jl+1, work, 1 ) )
226 END IF
227 40 CONTINUE
228*
229* Compute norm(L*U - A) / ( N * norm(A) * EPS )
230*
231 IF( anorm.LE.zero ) THEN
232 IF( resid.NE.zero )
233 $ resid = one / eps
234 ELSE
235 resid = ( ( resid / dble( n ) ) / anorm ) / eps
236 END IF
237*
238 RETURN
239*
240* End of ZGBT01
241*
double precision function dzasum(n, zx, incx)
DZASUM
Definition dzasum.f:72
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
double precision function dlamch(cmach)
DLAMCH
Definition dlamch.f:69
Here is the call graph for this function:
Here is the caller graph for this function: