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

◆ zla_gbrcond_x()

double precision function zla_gbrcond_x ( character  trans,
integer  n,
integer  kl,
integer  ku,
complex*16, dimension( ldab, * )  ab,
integer  ldab,
complex*16, dimension( ldafb, * )  afb,
integer  ldafb,
integer, dimension( * )  ipiv,
complex*16, dimension( * )  x,
integer  info,
complex*16, dimension( * )  work,
double precision, dimension( * )  rwork 
)

ZLA_GBRCOND_X computes the infinity norm condition number of op(A)*diag(x) for general banded matrices.

Download ZLA_GBRCOND_X + dependencies [TGZ] [ZIP] [TXT]

Purpose:
    ZLA_GBRCOND_X Computes the infinity norm condition number of
    op(A) * diag(X) where X is a COMPLEX*16 vector.
Parameters
[in]TRANS
          TRANS is CHARACTER*1
     Specifies the form of the system of equations:
       = 'N':  A * X = B     (No transpose)
       = 'T':  A**T * X = B  (Transpose)
       = 'C':  A**H * X = B  (Conjugate Transpose = Transpose)
[in]N
          N is INTEGER
     The number of linear equations, i.e., the order 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]AB
          AB is COMPLEX*16 array, dimension (LDAB,N)
     On entry, the matrix A in band storage, in rows 1 to KL+KU+1.
     The j-th column of A is stored in the j-th column of the
     array AB as follows:
     AB(KU+1+i-j,j) = A(i,j) for max(1,j-KU)<=i<=min(N,j+kl)
[in]LDAB
          LDAB is INTEGER
     The leading dimension of the array AB.  LDAB >= KL+KU+1.
[in]AFB
          AFB is COMPLEX*16 array, dimension (LDAFB,N)
     Details of the LU factorization of the band matrix A, 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.
[in]LDAFB
          LDAFB is INTEGER
     The leading dimension of the array AFB.  LDAFB >= 2*KL+KU+1.
[in]IPIV
          IPIV is INTEGER array, dimension (N)
     The pivot indices from the factorization A = P*L*U
     as computed by ZGBTRF; row i of the matrix was interchanged
     with row IPIV(i).
[in]X
          X is COMPLEX*16 array, dimension (N)
     The vector X in the formula op(A) * diag(X).
[out]INFO
          INFO is INTEGER
       = 0:  Successful exit.
     i > 0:  The ith argument is invalid.
[out]WORK
          WORK is COMPLEX*16 array, dimension (2*N).
     Workspace.
[out]RWORK
          RWORK is DOUBLE PRECISION array, dimension (N).
     Workspace.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.

Definition at line 152 of file zla_gbrcond_x.f.

155*
156* -- LAPACK computational routine --
157* -- LAPACK is a software package provided by Univ. of Tennessee, --
158* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
159*
160* .. Scalar Arguments ..
161 CHARACTER TRANS
162 INTEGER N, KL, KU, KD, KE, LDAB, LDAFB, INFO
163* ..
164* .. Array Arguments ..
165 INTEGER IPIV( * )
166 COMPLEX*16 AB( LDAB, * ), AFB( LDAFB, * ), WORK( * ),
167 $ X( * )
168 DOUBLE PRECISION RWORK( * )
169*
170*
171* =====================================================================
172*
173* .. Local Scalars ..
174 LOGICAL NOTRANS
175 INTEGER KASE, I, J
176 DOUBLE PRECISION AINVNM, ANORM, TMP
177 COMPLEX*16 ZDUM
178* ..
179* .. Local Arrays ..
180 INTEGER ISAVE( 3 )
181* ..
182* .. External Functions ..
183 LOGICAL LSAME
184 EXTERNAL lsame
185* ..
186* .. External Subroutines ..
187 EXTERNAL zlacn2, zgbtrs, xerbla
188* ..
189* .. Intrinsic Functions ..
190 INTRINSIC abs, max
191* ..
192* .. Statement Functions ..
193 DOUBLE PRECISION CABS1
194* ..
195* .. Statement Function Definitions ..
196 cabs1( zdum ) = abs( dble( zdum ) ) + abs( dimag( zdum ) )
197* ..
198* .. Executable Statements ..
199*
200 zla_gbrcond_x = 0.0d+0
201*
202 info = 0
203 notrans = lsame( trans, 'N' )
204 IF ( .NOT. notrans .AND. .NOT. lsame(trans, 'T') .AND. .NOT.
205 $ lsame( trans, 'C' ) ) THEN
206 info = -1
207 ELSE IF( n.LT.0 ) THEN
208 info = -2
209 ELSE IF( kl.LT.0 .OR. kl.GT.n-1 ) THEN
210 info = -3
211 ELSE IF( ku.LT.0 .OR. ku.GT.n-1 ) THEN
212 info = -4
213 ELSE IF( ldab.LT.kl+ku+1 ) THEN
214 info = -6
215 ELSE IF( ldafb.LT.2*kl+ku+1 ) THEN
216 info = -8
217 END IF
218 IF( info.NE.0 ) THEN
219 CALL xerbla( 'ZLA_GBRCOND_X', -info )
220 RETURN
221 END IF
222*
223* Compute norm of op(A)*op2(C).
224*
225 kd = ku + 1
226 ke = kl + 1
227 anorm = 0.0d+0
228 IF ( notrans ) THEN
229 DO i = 1, n
230 tmp = 0.0d+0
231 DO j = max( i-kl, 1 ), min( i+ku, n )
232 tmp = tmp + cabs1( ab( kd+i-j, j) * x( j ) )
233 END DO
234 rwork( i ) = tmp
235 anorm = max( anorm, tmp )
236 END DO
237 ELSE
238 DO i = 1, n
239 tmp = 0.0d+0
240 DO j = max( i-kl, 1 ), min( i+ku, n )
241 tmp = tmp + cabs1( ab( ke-i+j, i ) * x( j ) )
242 END DO
243 rwork( i ) = tmp
244 anorm = max( anorm, tmp )
245 END DO
246 END IF
247*
248* Quick return if possible.
249*
250 IF( n.EQ.0 ) THEN
251 zla_gbrcond_x = 1.0d+0
252 RETURN
253 ELSE IF( anorm .EQ. 0.0d+0 ) THEN
254 RETURN
255 END IF
256*
257* Estimate the norm of inv(op(A)).
258*
259 ainvnm = 0.0d+0
260*
261 kase = 0
262 10 CONTINUE
263 CALL zlacn2( n, work( n+1 ), work, ainvnm, kase, isave )
264 IF( kase.NE.0 ) THEN
265 IF( kase.EQ.2 ) THEN
266*
267* Multiply by R.
268*
269 DO i = 1, n
270 work( i ) = work( i ) * rwork( i )
271 END DO
272*
273 IF ( notrans ) THEN
274 CALL zgbtrs( 'No transpose', n, kl, ku, 1, afb, ldafb,
275 $ ipiv, work, n, info )
276 ELSE
277 CALL zgbtrs( 'Conjugate transpose', n, kl, ku, 1, afb,
278 $ ldafb, ipiv, work, n, info )
279 ENDIF
280*
281* Multiply by inv(X).
282*
283 DO i = 1, n
284 work( i ) = work( i ) / x( i )
285 END DO
286 ELSE
287*
288* Multiply by inv(X**H).
289*
290 DO i = 1, n
291 work( i ) = work( i ) / x( i )
292 END DO
293*
294 IF ( notrans ) THEN
295 CALL zgbtrs( 'Conjugate transpose', n, kl, ku, 1, afb,
296 $ ldafb, ipiv, work, n, info )
297 ELSE
298 CALL zgbtrs( 'No transpose', n, kl, ku, 1, afb, ldafb,
299 $ ipiv, work, n, info )
300 END IF
301*
302* Multiply by R.
303*
304 DO i = 1, n
305 work( i ) = work( i ) * rwork( i )
306 END DO
307 END IF
308 GO TO 10
309 END IF
310*
311* Compute the estimate of the reciprocal condition number.
312*
313 IF( ainvnm .NE. 0.0d+0 )
314 $ zla_gbrcond_x = 1.0d+0 / ainvnm
315*
316 RETURN
317*
318* End of ZLA_GBRCOND_X
319*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine zgbtrs(trans, n, kl, ku, nrhs, ab, ldab, ipiv, b, ldb, info)
ZGBTRS
Definition zgbtrs.f:138
double precision function zla_gbrcond_x(trans, n, kl, ku, ab, ldab, afb, ldafb, ipiv, x, info, work, rwork)
ZLA_GBRCOND_X computes the infinity norm condition number of op(A)*diag(x) for general banded matrice...
subroutine zlacn2(n, v, x, est, kase, isave)
ZLACN2 estimates the 1-norm of a square matrix, using reverse communication for evaluating matrix-vec...
Definition zlacn2.f:133
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
Here is the call graph for this function:
Here is the caller graph for this function: