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

◆ dgbcon()

subroutine dgbcon ( character norm,
integer n,
integer kl,
integer ku,
double precision, dimension( ldab, * ) ab,
integer ldab,
integer, dimension( * ) ipiv,
double precision anorm,
double precision rcond,
double precision, dimension( * ) work,
integer, dimension( * ) iwork,
integer info )

DGBCON

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

Purpose:
!>
!> DGBCON estimates the reciprocal of the condition number of a real
!> general band matrix A, in either the 1-norm or the infinity-norm,
!> using the LU factorization computed by DGBTRF.
!>
!> An estimate is obtained for norm(inv(A)), and the reciprocal of the
!> condition number is computed as
!>    RCOND = 1 / ( norm(A) * norm(inv(A)) ).
!> 
Parameters
[in]NORM
!>          NORM is CHARACTER*1
!>          Specifies whether the 1-norm condition number or the
!>          infinity-norm condition number is required:
!>          = '1' or 'O':  1-norm;
!>          = 'I':         Infinity-norm.
!> 
[in]N
!>          N is INTEGER
!>          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)
!>          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]LDAB
!>          LDAB is INTEGER
!>          The leading dimension of the array AB.  LDAB >= 2*KL+KU+1.
!> 
[in]IPIV
!>          IPIV is INTEGER array, dimension (N)
!>          The pivot indices; for 1 <= i <= N, row i of the matrix was
!>          interchanged with row IPIV(i).
!> 
[in]ANORM
!>          ANORM is DOUBLE PRECISION
!>          If NORM = '1' or 'O', the 1-norm of the original matrix A.
!>          If NORM = 'I', the infinity-norm of the original matrix A.
!> 
[out]RCOND
!>          RCOND is DOUBLE PRECISION
!>          The reciprocal of the condition number of the matrix A,
!>          computed as RCOND = 1/(norm(A) * norm(inv(A))).
!> 
[out]WORK
!>          WORK is DOUBLE PRECISION array, dimension (3*N)
!> 
[out]IWORK
!>          IWORK is INTEGER array, dimension (N)
!> 
[out]INFO
!>          INFO is INTEGER
!>          = 0:  successful exit
!>          < 0: if INFO = -i, the i-th argument had an illegal value
!> 
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.

Definition at line 142 of file dgbcon.f.

145*
146* -- LAPACK computational routine --
147* -- LAPACK is a software package provided by Univ. of Tennessee, --
148* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
149*
150* .. Scalar Arguments ..
151 CHARACTER NORM
152 INTEGER INFO, KL, KU, LDAB, N
153 DOUBLE PRECISION ANORM, RCOND
154* ..
155* .. Array Arguments ..
156 INTEGER IPIV( * ), IWORK( * )
157 DOUBLE PRECISION AB( LDAB, * ), WORK( * )
158* ..
159*
160* =====================================================================
161*
162* .. Parameters ..
163 DOUBLE PRECISION ONE, ZERO
164 parameter( one = 1.0d+0, zero = 0.0d+0 )
165* ..
166* .. Local Scalars ..
167 LOGICAL LNOTI, ONENRM
168 CHARACTER NORMIN
169 INTEGER IX, J, JP, KASE, KASE1, KD, LM
170 DOUBLE PRECISION AINVNM, SCALE, SMLNUM, T
171* ..
172* .. Local Arrays ..
173 INTEGER ISAVE( 3 )
174* ..
175* .. External Functions ..
176 LOGICAL LSAME
177 INTEGER IDAMAX
178 DOUBLE PRECISION DDOT, DLAMCH
179 EXTERNAL lsame, idamax, ddot, dlamch
180* ..
181* .. External Subroutines ..
182 EXTERNAL daxpy, dlacn2, dlatbs, drscl,
183 $ xerbla
184* ..
185* .. Intrinsic Functions ..
186 INTRINSIC abs, min
187* ..
188* .. Executable Statements ..
189*
190* Test the input parameters.
191*
192 info = 0
193 onenrm = norm.EQ.'1' .OR. lsame( norm, 'O' )
194 IF( .NOT.onenrm .AND. .NOT.lsame( norm, 'I' ) ) THEN
195 info = -1
196 ELSE IF( n.LT.0 ) THEN
197 info = -2
198 ELSE IF( kl.LT.0 ) THEN
199 info = -3
200 ELSE IF( ku.LT.0 ) THEN
201 info = -4
202 ELSE IF( ldab.LT.2*kl+ku+1 ) THEN
203 info = -6
204 ELSE IF( anorm.LT.zero ) THEN
205 info = -8
206 END IF
207 IF( info.NE.0 ) THEN
208 CALL xerbla( 'DGBCON', -info )
209 RETURN
210 END IF
211*
212* Quick return if possible
213*
214 rcond = zero
215 IF( n.EQ.0 ) THEN
216 rcond = one
217 RETURN
218 ELSE IF( anorm.EQ.zero ) THEN
219 RETURN
220 END IF
221*
222 smlnum = dlamch( 'Safe minimum' )
223*
224* Estimate the norm of inv(A).
225*
226 ainvnm = zero
227 normin = 'N'
228 IF( onenrm ) THEN
229 kase1 = 1
230 ELSE
231 kase1 = 2
232 END IF
233 kd = kl + ku + 1
234 lnoti = kl.GT.0
235 kase = 0
236 10 CONTINUE
237 CALL dlacn2( n, work( n+1 ), work, iwork, ainvnm, kase, isave )
238 IF( kase.NE.0 ) THEN
239 IF( kase.EQ.kase1 ) THEN
240*
241* Multiply by inv(L).
242*
243 IF( lnoti ) THEN
244 DO 20 j = 1, n - 1
245 lm = min( kl, n-j )
246 jp = ipiv( j )
247 t = work( jp )
248 IF( jp.NE.j ) THEN
249 work( jp ) = work( j )
250 work( j ) = t
251 END IF
252 CALL daxpy( lm, -t, ab( kd+1, j ), 1, work( j+1 ),
253 $ 1 )
254 20 CONTINUE
255 END IF
256*
257* Multiply by inv(U).
258*
259 CALL dlatbs( 'Upper', 'No transpose', 'Non-unit', normin,
260 $ n,
261 $ kl+ku, ab, ldab, work, scale, work( 2*n+1 ),
262 $ info )
263 ELSE
264*
265* Multiply by inv(U**T).
266*
267 CALL dlatbs( 'Upper', 'Transpose', 'Non-unit', normin, n,
268 $ kl+ku, ab, ldab, work, scale, work( 2*n+1 ),
269 $ info )
270*
271* Multiply by inv(L**T).
272*
273 IF( lnoti ) THEN
274 DO 30 j = n - 1, 1, -1
275 lm = min( kl, n-j )
276 work( j ) = work( j ) - ddot( lm, ab( kd+1, j ), 1,
277 $ work( j+1 ), 1 )
278 jp = ipiv( j )
279 IF( jp.NE.j ) THEN
280 t = work( jp )
281 work( jp ) = work( j )
282 work( j ) = t
283 END IF
284 30 CONTINUE
285 END IF
286 END IF
287*
288* Divide X by 1/SCALE if doing so will not cause overflow.
289*
290 normin = 'Y'
291 IF( scale.NE.one ) THEN
292 ix = idamax( n, work, 1 )
293 IF( scale.LT.abs( work( ix ) )*smlnum .OR. scale.EQ.zero )
294 $ GO TO 40
295 CALL drscl( n, scale, work, 1 )
296 END IF
297 GO TO 10
298 END IF
299*
300* Compute the estimate of the reciprocal condition number.
301*
302 IF( ainvnm.NE.zero )
303 $ rcond = ( one / ainvnm ) / anorm
304*
305 40 CONTINUE
306 RETURN
307*
308* End of DGBCON
309*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine daxpy(n, da, dx, incx, dy, incy)
DAXPY
Definition daxpy.f:89
double precision function ddot(n, dx, incx, dy, incy)
DDOT
Definition ddot.f:82
integer function idamax(n, dx, incx)
IDAMAX
Definition idamax.f:71
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
double precision function dlamch(cmach)
DLAMCH
Definition dlamch.f:69
subroutine dlatbs(uplo, trans, diag, normin, n, kd, ab, ldab, x, scale, cnorm, info)
DLATBS solves a triangular banded system of equations.
Definition dlatbs.f:241
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
subroutine drscl(n, sa, sx, incx)
DRSCL multiplies a vector by the reciprocal of a real scalar.
Definition drscl.f:82
Here is the call graph for this function:
Here is the caller graph for this function: