Annotation of rpl/lapack/lapack/dgbsv.f, revision 1.1
1.1 ! bertrand 1: SUBROUTINE DGBSV( N, KL, KU, NRHS, AB, LDAB, IPIV, B, LDB, INFO )
! 2: *
! 3: * -- LAPACK driver routine (version 3.2) --
! 4: * -- LAPACK is a software package provided by Univ. of Tennessee, --
! 5: * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
! 6: * November 2006
! 7: *
! 8: * .. Scalar Arguments ..
! 9: INTEGER INFO, KL, KU, LDAB, LDB, N, NRHS
! 10: * ..
! 11: * .. Array Arguments ..
! 12: INTEGER IPIV( * )
! 13: DOUBLE PRECISION AB( LDAB, * ), B( LDB, * )
! 14: * ..
! 15: *
! 16: * Purpose
! 17: * =======
! 18: *
! 19: * DGBSV computes the solution to a real system of linear equations
! 20: * A * X = B, where A is a band matrix of order N with KL subdiagonals
! 21: * and KU superdiagonals, and X and B are N-by-NRHS matrices.
! 22: *
! 23: * The LU decomposition with partial pivoting and row interchanges is
! 24: * used to factor A as A = L * U, where L is a product of permutation
! 25: * and unit lower triangular matrices with KL subdiagonals, and U is
! 26: * upper triangular with KL+KU superdiagonals. The factored form of A
! 27: * is then used to solve the system of equations A * X = B.
! 28: *
! 29: * Arguments
! 30: * =========
! 31: *
! 32: * N (input) INTEGER
! 33: * The number of linear equations, i.e., the order of the
! 34: * matrix A. N >= 0.
! 35: *
! 36: * KL (input) INTEGER
! 37: * The number of subdiagonals within the band of A. KL >= 0.
! 38: *
! 39: * KU (input) INTEGER
! 40: * The number of superdiagonals within the band of A. KU >= 0.
! 41: *
! 42: * NRHS (input) INTEGER
! 43: * The number of right hand sides, i.e., the number of columns
! 44: * of the matrix B. NRHS >= 0.
! 45: *
! 46: * AB (input/output) DOUBLE PRECISION array, dimension (LDAB,N)
! 47: * On entry, the matrix A in band storage, in rows KL+1 to
! 48: * 2*KL+KU+1; rows 1 to KL of the array need not be set.
! 49: * The j-th column of A is stored in the j-th column of the
! 50: * array AB as follows:
! 51: * AB(KL+KU+1+i-j,j) = A(i,j) for max(1,j-KU)<=i<=min(N,j+KL)
! 52: * On exit, details of the factorization: U is stored as an
! 53: * upper triangular band matrix with KL+KU superdiagonals in
! 54: * rows 1 to KL+KU+1, and the multipliers used during the
! 55: * factorization are stored in rows KL+KU+2 to 2*KL+KU+1.
! 56: * See below for further details.
! 57: *
! 58: * LDAB (input) INTEGER
! 59: * The leading dimension of the array AB. LDAB >= 2*KL+KU+1.
! 60: *
! 61: * IPIV (output) INTEGER array, dimension (N)
! 62: * The pivot indices that define the permutation matrix P;
! 63: * row i of the matrix was interchanged with row IPIV(i).
! 64: *
! 65: * B (input/output) DOUBLE PRECISION array, dimension (LDB,NRHS)
! 66: * On entry, the N-by-NRHS right hand side matrix B.
! 67: * On exit, if INFO = 0, the N-by-NRHS solution matrix X.
! 68: *
! 69: * LDB (input) INTEGER
! 70: * The leading dimension of the array B. LDB >= max(1,N).
! 71: *
! 72: * INFO (output) INTEGER
! 73: * = 0: successful exit
! 74: * < 0: if INFO = -i, the i-th argument had an illegal value
! 75: * > 0: if INFO = i, U(i,i) is exactly zero. The factorization
! 76: * has been completed, but the factor U is exactly
! 77: * singular, and the solution has not been computed.
! 78: *
! 79: * Further Details
! 80: * ===============
! 81: *
! 82: * The band storage scheme is illustrated by the following example, when
! 83: * M = N = 6, KL = 2, KU = 1:
! 84: *
! 85: * On entry: On exit:
! 86: *
! 87: * * * * + + + * * * u14 u25 u36
! 88: * * * + + + + * * u13 u24 u35 u46
! 89: * * a12 a23 a34 a45 a56 * u12 u23 u34 u45 u56
! 90: * a11 a22 a33 a44 a55 a66 u11 u22 u33 u44 u55 u66
! 91: * a21 a32 a43 a54 a65 * m21 m32 m43 m54 m65 *
! 92: * a31 a42 a53 a64 * * m31 m42 m53 m64 * *
! 93: *
! 94: * Array elements marked * are not used by the routine; elements marked
! 95: * + need not be set on entry, but are required by the routine to store
! 96: * elements of U because of fill-in resulting from the row interchanges.
! 97: *
! 98: * =====================================================================
! 99: *
! 100: * .. External Subroutines ..
! 101: EXTERNAL DGBTRF, DGBTRS, XERBLA
! 102: * ..
! 103: * .. Intrinsic Functions ..
! 104: INTRINSIC MAX
! 105: * ..
! 106: * .. Executable Statements ..
! 107: *
! 108: * Test the input parameters.
! 109: *
! 110: INFO = 0
! 111: IF( N.LT.0 ) THEN
! 112: INFO = -1
! 113: ELSE IF( KL.LT.0 ) THEN
! 114: INFO = -2
! 115: ELSE IF( KU.LT.0 ) THEN
! 116: INFO = -3
! 117: ELSE IF( NRHS.LT.0 ) THEN
! 118: INFO = -4
! 119: ELSE IF( LDAB.LT.2*KL+KU+1 ) THEN
! 120: INFO = -6
! 121: ELSE IF( LDB.LT.MAX( N, 1 ) ) THEN
! 122: INFO = -9
! 123: END IF
! 124: IF( INFO.NE.0 ) THEN
! 125: CALL XERBLA( 'DGBSV ', -INFO )
! 126: RETURN
! 127: END IF
! 128: *
! 129: * Compute the LU factorization of the band matrix A.
! 130: *
! 131: CALL DGBTRF( N, N, KL, KU, AB, LDAB, IPIV, INFO )
! 132: IF( INFO.EQ.0 ) THEN
! 133: *
! 134: * Solve the system A*X = B, overwriting B with X.
! 135: *
! 136: CALL DGBTRS( 'No transpose', N, KL, KU, NRHS, AB, LDAB, IPIV,
! 137: $ B, LDB, INFO )
! 138: END IF
! 139: RETURN
! 140: *
! 141: * End of DGBSV
! 142: *
! 143: END
CVSweb interface <joel.bertrand@systella.fr>