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

◆ dla_gbrcond()

double precision function dla_gbrcond ( character  trans,
integer  n,
integer  kl,
integer  ku,
double precision, dimension( ldab, * )  ab,
integer  ldab,
double precision, dimension( ldafb, * )  afb,
integer  ldafb,
integer, dimension( * )  ipiv,
integer  cmode,
double precision, dimension( * )  c,
integer  info,
double precision, dimension( * )  work,
integer, dimension( * )  iwork 
)

DLA_GBRCOND estimates the Skeel condition number for a general banded matrix.

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

Purpose:
    DLA_GBRCOND Estimates the Skeel condition number of  op(A) * op2(C)
    where op2 is determined by CMODE as follows
    CMODE =  1    op2(C) = C
    CMODE =  0    op2(C) = I
    CMODE = -1    op2(C) = inv(C)
    The Skeel condition number  cond(A) = norminf( |inv(A)||A| )
    is computed by computing scaling factors R such that
    diag(R)*A*op2(C) is row equilibrated and computing the standard
    infinity-norm condition number.
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 DOUBLE PRECISION 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 DOUBLE PRECISION array, dimension (LDAFB,N)
     Details of the LU factorization of the band matrix A, as
     computed by DGBTRF.  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 DGBTRF; row i of the matrix was interchanged
     with row IPIV(i).
[in]CMODE
          CMODE is INTEGER
     Determines op2(C) in the formula op(A) * op2(C) as follows:
     CMODE =  1    op2(C) = C
     CMODE =  0    op2(C) = I
     CMODE = -1    op2(C) = inv(C)
[in]C
          C is DOUBLE PRECISION array, dimension (N)
     The vector C in the formula op(A) * op2(C).
[out]INFO
          INFO is INTEGER
       = 0:  Successful exit.
     i > 0:  The ith argument is invalid.
[out]WORK
          WORK is DOUBLE PRECISION array, dimension (5*N).
     Workspace.
[out]IWORK
          IWORK is INTEGER array, dimension (N).
     Workspace.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.

Definition at line 167 of file dla_gbrcond.f.

170*
171* -- LAPACK computational routine --
172* -- LAPACK is a software package provided by Univ. of Tennessee, --
173* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
174*
175* .. Scalar Arguments ..
176 CHARACTER TRANS
177 INTEGER N, LDAB, LDAFB, INFO, KL, KU, CMODE
178* ..
179* .. Array Arguments ..
180 INTEGER IWORK( * ), IPIV( * )
181 DOUBLE PRECISION AB( LDAB, * ), AFB( LDAFB, * ), WORK( * ),
182 $ C( * )
183* ..
184*
185* =====================================================================
186*
187* .. Local Scalars ..
188 LOGICAL NOTRANS
189 INTEGER KASE, I, J, KD, KE
190 DOUBLE PRECISION AINVNM, TMP
191* ..
192* .. Local Arrays ..
193 INTEGER ISAVE( 3 )
194* ..
195* .. External Functions ..
196 LOGICAL LSAME
197 EXTERNAL lsame
198* ..
199* .. External Subroutines ..
200 EXTERNAL dlacn2, dgbtrs, xerbla
201* ..
202* .. Intrinsic Functions ..
203 INTRINSIC abs, max
204* ..
205* .. Executable Statements ..
206*
207 dla_gbrcond = 0.0d+0
208*
209 info = 0
210 notrans = lsame( trans, 'N' )
211 IF ( .NOT. notrans .AND. .NOT. lsame(trans, 'T')
212 $ .AND. .NOT. lsame(trans, 'C') ) THEN
213 info = -1
214 ELSE IF( n.LT.0 ) THEN
215 info = -2
216 ELSE IF( kl.LT.0 .OR. kl.GT.n-1 ) THEN
217 info = -3
218 ELSE IF( ku.LT.0 .OR. ku.GT.n-1 ) THEN
219 info = -4
220 ELSE IF( ldab.LT.kl+ku+1 ) THEN
221 info = -6
222 ELSE IF( ldafb.LT.2*kl+ku+1 ) THEN
223 info = -8
224 END IF
225 IF( info.NE.0 ) THEN
226 CALL xerbla( 'DLA_GBRCOND', -info )
227 RETURN
228 END IF
229 IF( n.EQ.0 ) THEN
230 dla_gbrcond = 1.0d+0
231 RETURN
232 END IF
233*
234* Compute the equilibration matrix R such that
235* inv(R)*A*C has unit 1-norm.
236*
237 kd = ku + 1
238 ke = kl + 1
239 IF ( notrans ) THEN
240 DO i = 1, n
241 tmp = 0.0d+0
242 IF ( cmode .EQ. 1 ) THEN
243 DO j = max( i-kl, 1 ), min( i+ku, n )
244 tmp = tmp + abs( ab( kd+i-j, j ) * c( j ) )
245 END DO
246 ELSE IF ( cmode .EQ. 0 ) THEN
247 DO j = max( i-kl, 1 ), min( i+ku, n )
248 tmp = tmp + abs( ab( kd+i-j, j ) )
249 END DO
250 ELSE
251 DO j = max( i-kl, 1 ), min( i+ku, n )
252 tmp = tmp + abs( ab( kd+i-j, j ) / c( j ) )
253 END DO
254 END IF
255 work( 2*n+i ) = tmp
256 END DO
257 ELSE
258 DO i = 1, n
259 tmp = 0.0d+0
260 IF ( cmode .EQ. 1 ) THEN
261 DO j = max( i-kl, 1 ), min( i+ku, n )
262 tmp = tmp + abs( ab( ke-i+j, i ) * c( j ) )
263 END DO
264 ELSE IF ( cmode .EQ. 0 ) THEN
265 DO j = max( i-kl, 1 ), min( i+ku, n )
266 tmp = tmp + abs( ab( ke-i+j, i ) )
267 END DO
268 ELSE
269 DO j = max( i-kl, 1 ), min( i+ku, n )
270 tmp = tmp + abs( ab( ke-i+j, i ) / c( j ) )
271 END DO
272 END IF
273 work( 2*n+i ) = tmp
274 END DO
275 END IF
276*
277* Estimate the norm of inv(op(A)).
278*
279 ainvnm = 0.0d+0
280
281 kase = 0
282 10 CONTINUE
283 CALL dlacn2( n, work( n+1 ), work, iwork, ainvnm, kase, isave )
284 IF( kase.NE.0 ) THEN
285 IF( kase.EQ.2 ) THEN
286*
287* Multiply by R.
288*
289 DO i = 1, n
290 work( i ) = work( i ) * work( 2*n+i )
291 END DO
292
293 IF ( notrans ) THEN
294 CALL dgbtrs( 'No transpose', n, kl, ku, 1, afb, ldafb,
295 $ ipiv, work, n, info )
296 ELSE
297 CALL dgbtrs( 'Transpose', n, kl, ku, 1, afb, ldafb, ipiv,
298 $ work, n, info )
299 END IF
300*
301* Multiply by inv(C).
302*
303 IF ( cmode .EQ. 1 ) THEN
304 DO i = 1, n
305 work( i ) = work( i ) / c( i )
306 END DO
307 ELSE IF ( cmode .EQ. -1 ) THEN
308 DO i = 1, n
309 work( i ) = work( i ) * c( i )
310 END DO
311 END IF
312 ELSE
313*
314* Multiply by inv(C**T).
315*
316 IF ( cmode .EQ. 1 ) THEN
317 DO i = 1, n
318 work( i ) = work( i ) / c( i )
319 END DO
320 ELSE IF ( cmode .EQ. -1 ) THEN
321 DO i = 1, n
322 work( i ) = work( i ) * c( i )
323 END DO
324 END IF
325
326 IF ( notrans ) THEN
327 CALL dgbtrs( 'Transpose', n, kl, ku, 1, afb, ldafb, ipiv,
328 $ work, n, info )
329 ELSE
330 CALL dgbtrs( 'No transpose', n, kl, ku, 1, afb, ldafb,
331 $ ipiv, work, n, info )
332 END IF
333*
334* Multiply by R.
335*
336 DO i = 1, n
337 work( i ) = work( i ) * work( 2*n+i )
338 END DO
339 END IF
340 GO TO 10
341 END IF
342*
343* Compute the estimate of the reciprocal condition number.
344*
345 IF( ainvnm .NE. 0.0d+0 )
346 $ dla_gbrcond = ( 1.0d+0 / ainvnm )
347*
348 RETURN
349*
350* End of DLA_GBRCOND
351*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine dgbtrs(trans, n, kl, ku, nrhs, ab, ldab, ipiv, b, ldb, info)
DGBTRS
Definition dgbtrs.f:138
double precision function dla_gbrcond(trans, n, kl, ku, ab, ldab, afb, ldafb, ipiv, cmode, c, info, work, iwork)
DLA_GBRCOND estimates the Skeel condition number for a general banded matrix.
subroutine dlacn2(n, v, x, isgn, est, kase, isave)
DLACN2 estimates the 1-norm of a square matrix, using reverse communication for evaluating matrix-vec...
Definition dlacn2.f:136
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: