LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Modules Pages

◆ dpbsvx()

subroutine dpbsvx ( character fact,
character uplo,
integer n,
integer kd,
integer nrhs,
double precision, dimension( ldab, * ) ab,
integer ldab,
double precision, dimension( ldafb, * ) afb,
integer ldafb,
character equed,
double precision, dimension( * ) s,
double precision, dimension( ldb, * ) b,
integer ldb,
double precision, dimension( ldx, * ) x,
integer ldx,
double precision rcond,
double precision, dimension( * ) ferr,
double precision, dimension( * ) berr,
double precision, dimension( * ) work,
integer, dimension( * ) iwork,
integer info )

DPBSVX computes the solution to system of linear equations A * X = B for OTHER matrices

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

Purpose:
!> !> DPBSVX uses the Cholesky factorization A = U**T*U or A = L*L**T to !> compute the solution to a real system of linear equations !> A * X = B, !> where A is an N-by-N symmetric positive definite band matrix and X !> and B are N-by-NRHS matrices. !> !> Error bounds on the solution and a condition estimate are also !> provided. !>
Description:
!> !> The following steps are performed: !> !> 1. If FACT = 'E', real scaling factors are computed to equilibrate !> the system: !> diag(S) * A * diag(S) * inv(diag(S)) * X = diag(S) * B !> Whether or not the system will be equilibrated depends on the !> scaling of the matrix A, but if equilibration is used, A is !> overwritten by diag(S)*A*diag(S) and B by diag(S)*B. !> !> 2. If FACT = 'N' or 'E', the Cholesky decomposition is used to !> factor the matrix A (after equilibration if FACT = 'E') as !> A = U**T * U, if UPLO = 'U', or !> A = L * L**T, if UPLO = 'L', !> where U is an upper triangular band matrix, and L is a lower !> triangular band matrix. !> !> 3. If the leading principal minor of order i is not positive, !> then the routine returns with INFO = i. Otherwise, the factored !> form of A is used to estimate the condition number of the matrix !> A. If the reciprocal of the condition number is less than machine !> precision, INFO = N+1 is returned as a warning, but the routine !> still goes on to solve for X and compute error bounds as !> described below. !> !> 4. The system of equations is solved for X using the factored form !> of A. !> !> 5. Iterative refinement is applied to improve the computed solution !> matrix and calculate error bounds and backward error estimates !> for it. !> !> 6. If equilibration was used, the matrix X is premultiplied by !> diag(S) so that it solves the original system before !> equilibration. !>
Parameters
[in]FACT
!> FACT is CHARACTER*1 !> Specifies whether or not the factored form of the matrix A is !> supplied on entry, and if not, whether the matrix A should be !> equilibrated before it is factored. !> = 'F': On entry, AFB contains the factored form of A. !> If EQUED = 'Y', the matrix A has been equilibrated !> with scaling factors given by S. AB and AFB will not !> be modified. !> = 'N': The matrix A will be copied to AFB and factored. !> = 'E': The matrix A will be equilibrated if necessary, then !> copied to AFB and factored. !>
[in]UPLO
!> UPLO is CHARACTER*1 !> = 'U': Upper triangle of A is stored; !> = 'L': Lower triangle of A is stored. !>
[in]N
!> N is INTEGER !> The number of linear equations, i.e., the order of the !> matrix A. N >= 0. !>
[in]KD
!> KD is INTEGER !> The number of superdiagonals of the matrix A if UPLO = 'U', !> or the number of subdiagonals if UPLO = 'L'. KD >= 0. !>
[in]NRHS
!> NRHS is INTEGER !> The number of right-hand sides, i.e., the number of columns !> of the matrices B and X. NRHS >= 0. !>
[in,out]AB
!> AB is DOUBLE PRECISION array, dimension (LDAB,N) !> On entry, the upper or lower triangle of the symmetric band !> matrix A, stored in the first KD+1 rows of the array, except !> if FACT = 'F' and EQUED = 'Y', then A must contain the !> equilibrated matrix diag(S)*A*diag(S). The j-th column of A !> is stored in the j-th column of the array AB as follows: !> if UPLO = 'U', AB(KD+1+i-j,j) = A(i,j) for max(1,j-KD)<=i<=j; !> if UPLO = 'L', AB(1+i-j,j) = A(i,j) for j<=i<=min(N,j+KD). !> See below for further details. !> !> On exit, if FACT = 'E' and EQUED = 'Y', A is overwritten by !> diag(S)*A*diag(S). !>
[in]LDAB
!> LDAB is INTEGER !> The leading dimension of the array A. LDAB >= KD+1. !>
[in,out]AFB
!> AFB is DOUBLE PRECISION array, dimension (LDAFB,N) !> If FACT = 'F', then AFB is an input argument and on entry !> contains the triangular factor U or L from the Cholesky !> factorization A = U**T*U or A = L*L**T of the band matrix !> A, in the same storage format as A (see AB). If EQUED = 'Y', !> then AFB is the factored form of the equilibrated matrix A. !> !> If FACT = 'N', then AFB is an output argument and on exit !> returns the triangular factor U or L from the Cholesky !> factorization A = U**T*U or A = L*L**T. !> !> If FACT = 'E', then AFB is an output argument and on exit !> returns the triangular factor U or L from the Cholesky !> factorization A = U**T*U or A = L*L**T of the equilibrated !> matrix A (see the description of A for the form of the !> equilibrated matrix). !>
[in]LDAFB
!> LDAFB is INTEGER !> The leading dimension of the array AFB. LDAFB >= KD+1. !>
[in,out]EQUED
!> EQUED is CHARACTER*1 !> Specifies the form of equilibration that was done. !> = 'N': No equilibration (always true if FACT = 'N'). !> = 'Y': Equilibration was done, i.e., A has been replaced by !> diag(S) * A * diag(S). !> EQUED is an input argument if FACT = 'F'; otherwise, it is an !> output argument. !>
[in,out]S
!> S is DOUBLE PRECISION array, dimension (N) !> The scale factors for A; not accessed if EQUED = 'N'. S is !> an input argument if FACT = 'F'; otherwise, S is an output !> argument. If FACT = 'F' and EQUED = 'Y', each element of S !> must be positive. !>
[in,out]B
!> B is DOUBLE PRECISION array, dimension (LDB,NRHS) !> On entry, the N-by-NRHS right hand side matrix B. !> On exit, if EQUED = 'N', B is not modified; if EQUED = 'Y', !> B is overwritten by diag(S) * B. !>
[in]LDB
!> LDB is INTEGER !> The leading dimension of the array B. LDB >= max(1,N). !>
[out]X
!> X is DOUBLE PRECISION array, dimension (LDX,NRHS) !> If INFO = 0 or INFO = N+1, the N-by-NRHS solution matrix X to !> the original system of equations. Note that if EQUED = 'Y', !> A and B are modified on exit, and the solution to the !> equilibrated system is inv(diag(S))*X. !>
[in]LDX
!> LDX is INTEGER !> The leading dimension of the array X. LDX >= max(1,N). !>
[out]RCOND
!> RCOND is DOUBLE PRECISION !> The estimate of the reciprocal condition number of the matrix !> A after equilibration (if done). If RCOND is less than the !> machine precision (in particular, if RCOND = 0), the matrix !> is singular to working precision. This condition is !> indicated by a return code of INFO > 0. !>
[out]FERR
!> FERR is DOUBLE PRECISION array, dimension (NRHS) !> The estimated forward error bound for each solution vector !> X(j) (the j-th column of the solution matrix X). !> If XTRUE is the true solution corresponding to X(j), FERR(j) !> is an estimated upper bound for the magnitude of the largest !> element in (X(j) - XTRUE) divided by the magnitude of the !> largest element in X(j). The estimate is as reliable as !> the estimate for RCOND, and is almost always a slight !> overestimate of the true error. !>
[out]BERR
!> BERR is DOUBLE PRECISION array, dimension (NRHS) !> The componentwise relative backward error of each solution !> vector X(j) (i.e., the smallest relative change in !> any element of A or B that makes X(j) an exact solution). !>
[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 !> > 0: if INFO = i, and i is !> <= N: the leading principal minor of order i of A !> is not positive, so the factorization could not !> be completed, and the solution has not been !> computed. RCOND = 0 is returned. !> = N+1: U is nonsingular, but RCOND is less than machine !> precision, meaning that the matrix is singular !> to working precision. Nevertheless, the !> solution and error bounds are computed because !> there are a number of situations where the !> computed solution can be more accurate than the !> value of RCOND would suggest. !>
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Further Details:
!> !> The band storage scheme is illustrated by the following example, when !> N = 6, KD = 2, and UPLO = 'U': !> !> Two-dimensional storage of the symmetric matrix A: !> !> a11 a12 a13 !> a22 a23 a24 !> a33 a34 a35 !> a44 a45 a46 !> a55 a56 !> (aij=conjg(aji)) a66 !> !> Band storage of the upper triangle of A: !> !> * * a13 a24 a35 a46 !> * a12 a23 a34 a45 a56 !> a11 a22 a33 a44 a55 a66 !> !> Similarly, if UPLO = 'L' the format of A is as follows: !> !> a11 a22 a33 a44 a55 a66 !> a21 a32 a43 a54 a65 * !> a31 a42 a53 a64 * * !> !> Array elements marked * are not used by the routine. !>

Definition at line 338 of file dpbsvx.f.

342*
343* -- LAPACK driver routine --
344* -- LAPACK is a software package provided by Univ. of Tennessee, --
345* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
346*
347* .. Scalar Arguments ..
348 CHARACTER EQUED, FACT, UPLO
349 INTEGER INFO, KD, LDAB, LDAFB, LDB, LDX, N, NRHS
350 DOUBLE PRECISION RCOND
351* ..
352* .. Array Arguments ..
353 INTEGER IWORK( * )
354 DOUBLE PRECISION AB( LDAB, * ), AFB( LDAFB, * ), B( LDB, * ),
355 $ BERR( * ), FERR( * ), S( * ), WORK( * ),
356 $ X( LDX, * )
357* ..
358*
359* =====================================================================
360*
361* .. Parameters ..
362 DOUBLE PRECISION ZERO, ONE
363 parameter( zero = 0.0d+0, one = 1.0d+0 )
364* ..
365* .. Local Scalars ..
366 LOGICAL EQUIL, NOFACT, RCEQU, UPPER
367 INTEGER I, INFEQU, J, J1, J2
368 DOUBLE PRECISION AMAX, ANORM, BIGNUM, SCOND, SMAX, SMIN, SMLNUM
369* ..
370* .. External Functions ..
371 LOGICAL LSAME
372 DOUBLE PRECISION DLAMCH, DLANSB
373 EXTERNAL lsame, dlamch, dlansb
374* ..
375* .. External Subroutines ..
376 EXTERNAL dcopy, dlacpy, dlaqsb, dpbcon, dpbequ,
377 $ dpbrfs,
379* ..
380* .. Intrinsic Functions ..
381 INTRINSIC max, min
382* ..
383* .. Executable Statements ..
384*
385 info = 0
386 nofact = lsame( fact, 'N' )
387 equil = lsame( fact, 'E' )
388 upper = lsame( uplo, 'U' )
389 IF( nofact .OR. equil ) THEN
390 equed = 'N'
391 rcequ = .false.
392 ELSE
393 rcequ = lsame( equed, 'Y' )
394 smlnum = dlamch( 'Safe minimum' )
395 bignum = one / smlnum
396 END IF
397*
398* Test the input parameters.
399*
400 IF( .NOT.nofact .AND.
401 $ .NOT.equil .AND.
402 $ .NOT.lsame( fact, 'F' ) )
403 $ THEN
404 info = -1
405 ELSE IF( .NOT.upper .AND. .NOT.lsame( uplo, 'L' ) ) THEN
406 info = -2
407 ELSE IF( n.LT.0 ) THEN
408 info = -3
409 ELSE IF( kd.LT.0 ) THEN
410 info = -4
411 ELSE IF( nrhs.LT.0 ) THEN
412 info = -5
413 ELSE IF( ldab.LT.kd+1 ) THEN
414 info = -7
415 ELSE IF( ldafb.LT.kd+1 ) THEN
416 info = -9
417 ELSE IF( lsame( fact, 'F' ) .AND. .NOT.
418 $ ( rcequ .OR. lsame( equed, 'N' ) ) ) THEN
419 info = -10
420 ELSE
421 IF( rcequ ) THEN
422 smin = bignum
423 smax = zero
424 DO 10 j = 1, n
425 smin = min( smin, s( j ) )
426 smax = max( smax, s( j ) )
427 10 CONTINUE
428 IF( smin.LE.zero ) THEN
429 info = -11
430 ELSE IF( n.GT.0 ) THEN
431 scond = max( smin, smlnum ) / min( smax, bignum )
432 ELSE
433 scond = one
434 END IF
435 END IF
436 IF( info.EQ.0 ) THEN
437 IF( ldb.LT.max( 1, n ) ) THEN
438 info = -13
439 ELSE IF( ldx.LT.max( 1, n ) ) THEN
440 info = -15
441 END IF
442 END IF
443 END IF
444*
445 IF( info.NE.0 ) THEN
446 CALL xerbla( 'DPBSVX', -info )
447 RETURN
448 END IF
449*
450 IF( equil ) THEN
451*
452* Compute row and column scalings to equilibrate the matrix A.
453*
454 CALL dpbequ( uplo, n, kd, ab, ldab, s, scond, amax, infequ )
455 IF( infequ.EQ.0 ) THEN
456*
457* Equilibrate the matrix.
458*
459 CALL dlaqsb( uplo, n, kd, ab, ldab, s, scond, amax,
460 $ equed )
461 rcequ = lsame( equed, 'Y' )
462 END IF
463 END IF
464*
465* Scale the right-hand side.
466*
467 IF( rcequ ) THEN
468 DO 30 j = 1, nrhs
469 DO 20 i = 1, n
470 b( i, j ) = s( i )*b( i, j )
471 20 CONTINUE
472 30 CONTINUE
473 END IF
474*
475 IF( nofact .OR. equil ) THEN
476*
477* Compute the Cholesky factorization A = U**T *U or A = L*L**T.
478*
479 IF( upper ) THEN
480 DO 40 j = 1, n
481 j1 = max( j-kd, 1 )
482 CALL dcopy( j-j1+1, ab( kd+1-j+j1, j ), 1,
483 $ afb( kd+1-j+j1, j ), 1 )
484 40 CONTINUE
485 ELSE
486 DO 50 j = 1, n
487 j2 = min( j+kd, n )
488 CALL dcopy( j2-j+1, ab( 1, j ), 1, afb( 1, j ), 1 )
489 50 CONTINUE
490 END IF
491*
492 CALL dpbtrf( uplo, n, kd, afb, ldafb, info )
493*
494* Return if INFO is non-zero.
495*
496 IF( info.GT.0 )THEN
497 rcond = zero
498 RETURN
499 END IF
500 END IF
501*
502* Compute the norm of the matrix A.
503*
504 anorm = dlansb( '1', uplo, n, kd, ab, ldab, work )
505*
506* Compute the reciprocal of the condition number of A.
507*
508 CALL dpbcon( uplo, n, kd, afb, ldafb, anorm, rcond, work,
509 $ iwork,
510 $ info )
511*
512* Compute the solution matrix X.
513*
514 CALL dlacpy( 'Full', n, nrhs, b, ldb, x, ldx )
515 CALL dpbtrs( uplo, n, kd, nrhs, afb, ldafb, x, ldx, info )
516*
517* Use iterative refinement to improve the computed solution and
518* compute error bounds and backward error estimates for it.
519*
520 CALL dpbrfs( uplo, n, kd, nrhs, ab, ldab, afb, ldafb, b, ldb,
521 $ x,
522 $ ldx, ferr, berr, work, iwork, info )
523*
524* Transform the solution matrix X to a solution of the original
525* system.
526*
527 IF( rcequ ) THEN
528 DO 70 j = 1, nrhs
529 DO 60 i = 1, n
530 x( i, j ) = s( i )*x( i, j )
531 60 CONTINUE
532 70 CONTINUE
533 DO 80 j = 1, nrhs
534 ferr( j ) = ferr( j ) / scond
535 80 CONTINUE
536 END IF
537*
538* Set INFO = N+1 if the matrix is singular to working precision.
539*
540 IF( rcond.LT.dlamch( 'Epsilon' ) )
541 $ info = n + 1
542*
543 RETURN
544*
545* End of DPBSVX
546*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine dcopy(n, dx, incx, dy, incy)
DCOPY
Definition dcopy.f:82
subroutine dlacpy(uplo, m, n, a, lda, b, ldb)
DLACPY copies all or part of one two-dimensional array to another.
Definition dlacpy.f:101
double precision function dlamch(cmach)
DLAMCH
Definition dlamch.f:69
double precision function dlansb(norm, uplo, n, k, ab, ldab, work)
DLANSB returns the value of the 1-norm, or the Frobenius norm, or the infinity norm,...
Definition dlansb.f:127
subroutine dlaqsb(uplo, n, kd, ab, ldab, s, scond, amax, equed)
DLAQSB scales a symmetric/Hermitian band matrix, using scaling factors computed by spbequ.
Definition dlaqsb.f:139
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
subroutine dpbcon(uplo, n, kd, ab, ldab, anorm, rcond, work, iwork, info)
DPBCON
Definition dpbcon.f:130
subroutine dpbequ(uplo, n, kd, ab, ldab, s, scond, amax, info)
DPBEQU
Definition dpbequ.f:128
subroutine dpbrfs(uplo, n, kd, nrhs, ab, ldab, afb, ldafb, b, ldb, x, ldx, ferr, berr, work, iwork, info)
DPBRFS
Definition dpbrfs.f:187
subroutine dpbtrf(uplo, n, kd, ab, ldab, info)
DPBTRF
Definition dpbtrf.f:140
subroutine dpbtrs(uplo, n, kd, nrhs, ab, ldab, b, ldb, info)
DPBTRS
Definition dpbtrs.f:119
Here is the call graph for this function:
Here is the caller graph for this function: