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

◆ sgbcon()

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

SGBCON

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

Purpose:
!>
!> SGBCON 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 SGBTRF.
!>
!> 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 REAL array, dimension (LDAB,N)
!>          Details of the LU factorization of the band matrix A, as
!>          computed by SGBTRF.  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 REAL
!>          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 REAL
!>          The reciprocal of the condition number of the matrix A,
!>          computed as RCOND = 1/(norm(A) * norm(inv(A))).
!> 
[out]WORK
!>          WORK is REAL 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 sgbcon.f.

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