1
2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 |
SUBROUTINE SGTTRF( N, DL, D, DU, DU2, IPIV, INFO )
* * -- LAPACK routine (version 3.2) -- * -- LAPACK is a software package provided by Univ. of Tennessee, -- * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..-- * November 2006 * * .. Scalar Arguments .. INTEGER INFO, N * .. * .. Array Arguments .. INTEGER IPIV( * ) REAL D( * ), DL( * ), DU( * ), DU2( * ) * .. * * Purpose * ======= * * SGTTRF computes an LU factorization of a real tridiagonal matrix A * using elimination with partial pivoting and row interchanges. * * The factorization has the form * A = L * U * where L is a product of permutation and unit lower bidiagonal * matrices and U is upper triangular with nonzeros in only the main * diagonal and first two superdiagonals. * * Arguments * ========= * * N (input) INTEGER * The order of the matrix A. * * DL (input/output) REAL array, dimension (N-1) * On entry, DL must contain the (n-1) sub-diagonal elements of * A. * * On exit, DL is overwritten by the (n-1) multipliers that * define the matrix L from the LU factorization of A. * * D (input/output) REAL array, dimension (N) * On entry, D must contain the diagonal elements of A. * * On exit, D is overwritten by the n diagonal elements of the * upper triangular matrix U from the LU factorization of A. * * DU (input/output) REAL array, dimension (N-1) * On entry, DU must contain the (n-1) super-diagonal elements * of A. * * On exit, DU is overwritten by the (n-1) elements of the first * super-diagonal of U. * * DU2 (output) REAL array, dimension (N-2) * On exit, DU2 is overwritten by the (n-2) elements of the * second super-diagonal of U. * * IPIV (output) INTEGER array, dimension (N) * The pivot indices; for 1 <= i <= n, row i of the matrix was * interchanged with row IPIV(i). IPIV(i) will always be either * i or i+1; IPIV(i) = i indicates a row interchange was not * required. * * INFO (output) INTEGER * = 0: successful exit * < 0: if INFO = -k, the k-th argument had an illegal value * > 0: if INFO = k, U(k,k) is exactly zero. The factorization * has been completed, but the factor U is exactly * singular, and division by zero will occur if it is used * to solve a system of equations. * * ===================================================================== * * .. Parameters .. REAL ZERO PARAMETER ( ZERO = 0.0E+0 ) * .. * .. Local Scalars .. INTEGER I REAL FACT, TEMP * .. * .. Intrinsic Functions .. INTRINSIC ABS * .. * .. External Subroutines .. EXTERNAL XERBLA * .. * .. Executable Statements .. * INFO = 0 IF( N.LT.0 ) THEN INFO = -1 CALL XERBLA( 'SGTTRF', -INFO ) RETURN END IF * * Quick return if possible * IF( N.EQ.0 ) $ RETURN * * Initialize IPIV(i) = i and DU2(I) = 0 * DO 10 I = 1, N IPIV( I ) = I 10 CONTINUE DO 20 I = 1, N - 2 DU2( I ) = ZERO 20 CONTINUE * DO 30 I = 1, N - 2 IF( ABS( D( I ) ).GE.ABS( DL( I ) ) ) THEN * * No row interchange required, eliminate DL(I) * IF( D( I ).NE.ZERO ) THEN FACT = DL( I ) / D( I ) DL( I ) = FACT D( I+1 ) = D( I+1 ) - FACT*DU( I ) END IF ELSE * * Interchange rows I and I+1, eliminate DL(I) * FACT = D( I ) / DL( I ) D( I ) = DL( I ) DL( I ) = FACT TEMP = DU( I ) DU( I ) = D( I+1 ) D( I+1 ) = TEMP - FACT*D( I+1 ) DU2( I ) = DU( I+1 ) DU( I+1 ) = -FACT*DU( I+1 ) IPIV( I ) = I + 1 END IF 30 CONTINUE IF( N.GT.1 ) THEN I = N - 1 IF( ABS( D( I ) ).GE.ABS( DL( I ) ) ) THEN IF( D( I ).NE.ZERO ) THEN FACT = DL( I ) / D( I ) DL( I ) = FACT D( I+1 ) = D( I+1 ) - FACT*DU( I ) END IF ELSE FACT = D( I ) / DL( I ) D( I ) = DL( I ) DL( I ) = FACT TEMP = DU( I ) DU( I ) = D( I+1 ) D( I+1 ) = TEMP - FACT*D( I+1 ) IPIV( I ) = I + 1 END IF END IF * * Check for a zero on the diagonal of U. * DO 40 I = 1, N IF( D( I ).EQ.ZERO ) THEN INFO = I GO TO 50 END IF 40 CONTINUE 50 CONTINUE * RETURN * * End of SGTTRF * END |