LAPACK 3.12.1
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 165 of file dla_gbrcond.f.

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