Annotation of rpl/lapack/lapack/dgttrs.f, revision 1.9
1.9 ! bertrand 1: *> \brief \b DGTTRS
! 2: *
! 3: * =========== DOCUMENTATION ===========
! 4: *
! 5: * Online html documentation available at
! 6: * http://www.netlib.org/lapack/explore-html/
! 7: *
! 8: *> \htmlonly
! 9: *> Download DGTTRS + dependencies
! 10: *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dgttrs.f">
! 11: *> [TGZ]</a>
! 12: *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dgttrs.f">
! 13: *> [ZIP]</a>
! 14: *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dgttrs.f">
! 15: *> [TXT]</a>
! 16: *> \endhtmlonly
! 17: *
! 18: * Definition:
! 19: * ===========
! 20: *
! 21: * SUBROUTINE DGTTRS( TRANS, N, NRHS, DL, D, DU, DU2, IPIV, B, LDB,
! 22: * INFO )
! 23: *
! 24: * .. Scalar Arguments ..
! 25: * CHARACTER TRANS
! 26: * INTEGER INFO, LDB, N, NRHS
! 27: * ..
! 28: * .. Array Arguments ..
! 29: * INTEGER IPIV( * )
! 30: * DOUBLE PRECISION B( LDB, * ), D( * ), DL( * ), DU( * ), DU2( * )
! 31: * ..
! 32: *
! 33: *
! 34: *> \par Purpose:
! 35: * =============
! 36: *>
! 37: *> \verbatim
! 38: *>
! 39: *> DGTTRS solves one of the systems of equations
! 40: *> A*X = B or A**T*X = B,
! 41: *> with a tridiagonal matrix A using the LU factorization computed
! 42: *> by DGTTRF.
! 43: *> \endverbatim
! 44: *
! 45: * Arguments:
! 46: * ==========
! 47: *
! 48: *> \param[in] TRANS
! 49: *> \verbatim
! 50: *> TRANS is CHARACTER*1
! 51: *> Specifies the form of the system of equations.
! 52: *> = 'N': A * X = B (No transpose)
! 53: *> = 'T': A**T* X = B (Transpose)
! 54: *> = 'C': A**T* X = B (Conjugate transpose = Transpose)
! 55: *> \endverbatim
! 56: *>
! 57: *> \param[in] N
! 58: *> \verbatim
! 59: *> N is INTEGER
! 60: *> The order of the matrix A.
! 61: *> \endverbatim
! 62: *>
! 63: *> \param[in] NRHS
! 64: *> \verbatim
! 65: *> NRHS is INTEGER
! 66: *> The number of right hand sides, i.e., the number of columns
! 67: *> of the matrix B. NRHS >= 0.
! 68: *> \endverbatim
! 69: *>
! 70: *> \param[in] DL
! 71: *> \verbatim
! 72: *> DL is DOUBLE PRECISION array, dimension (N-1)
! 73: *> The (n-1) multipliers that define the matrix L from the
! 74: *> LU factorization of A.
! 75: *> \endverbatim
! 76: *>
! 77: *> \param[in] D
! 78: *> \verbatim
! 79: *> D is DOUBLE PRECISION array, dimension (N)
! 80: *> The n diagonal elements of the upper triangular matrix U from
! 81: *> the LU factorization of A.
! 82: *> \endverbatim
! 83: *>
! 84: *> \param[in] DU
! 85: *> \verbatim
! 86: *> DU is DOUBLE PRECISION array, dimension (N-1)
! 87: *> The (n-1) elements of the first super-diagonal of U.
! 88: *> \endverbatim
! 89: *>
! 90: *> \param[in] DU2
! 91: *> \verbatim
! 92: *> DU2 is DOUBLE PRECISION array, dimension (N-2)
! 93: *> The (n-2) elements of the second super-diagonal of U.
! 94: *> \endverbatim
! 95: *>
! 96: *> \param[in] IPIV
! 97: *> \verbatim
! 98: *> IPIV is INTEGER array, dimension (N)
! 99: *> The pivot indices; for 1 <= i <= n, row i of the matrix was
! 100: *> interchanged with row IPIV(i). IPIV(i) will always be either
! 101: *> i or i+1; IPIV(i) = i indicates a row interchange was not
! 102: *> required.
! 103: *> \endverbatim
! 104: *>
! 105: *> \param[in,out] B
! 106: *> \verbatim
! 107: *> B is DOUBLE PRECISION array, dimension (LDB,NRHS)
! 108: *> On entry, the matrix of right hand side vectors B.
! 109: *> On exit, B is overwritten by the solution vectors X.
! 110: *> \endverbatim
! 111: *>
! 112: *> \param[in] LDB
! 113: *> \verbatim
! 114: *> LDB is INTEGER
! 115: *> The leading dimension of the array B. LDB >= max(1,N).
! 116: *> \endverbatim
! 117: *>
! 118: *> \param[out] INFO
! 119: *> \verbatim
! 120: *> INFO is INTEGER
! 121: *> = 0: successful exit
! 122: *> < 0: if INFO = -i, the i-th argument had an illegal value
! 123: *> \endverbatim
! 124: *
! 125: * Authors:
! 126: * ========
! 127: *
! 128: *> \author Univ. of Tennessee
! 129: *> \author Univ. of California Berkeley
! 130: *> \author Univ. of Colorado Denver
! 131: *> \author NAG Ltd.
! 132: *
! 133: *> \date November 2011
! 134: *
! 135: *> \ingroup doubleOTHERcomputational
! 136: *
! 137: * =====================================================================
1.1 bertrand 138: SUBROUTINE DGTTRS( TRANS, N, NRHS, DL, D, DU, DU2, IPIV, B, LDB,
139: $ INFO )
140: *
1.9 ! bertrand 141: * -- LAPACK computational routine (version 3.4.0) --
1.1 bertrand 142: * -- LAPACK is a software package provided by Univ. of Tennessee, --
143: * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
1.9 ! bertrand 144: * November 2011
1.1 bertrand 145: *
146: * .. Scalar Arguments ..
147: CHARACTER TRANS
148: INTEGER INFO, LDB, N, NRHS
149: * ..
150: * .. Array Arguments ..
151: INTEGER IPIV( * )
152: DOUBLE PRECISION B( LDB, * ), D( * ), DL( * ), DU( * ), DU2( * )
153: * ..
154: *
155: * =====================================================================
156: *
157: * .. Local Scalars ..
158: LOGICAL NOTRAN
159: INTEGER ITRANS, J, JB, NB
160: * ..
161: * .. External Functions ..
162: INTEGER ILAENV
163: EXTERNAL ILAENV
164: * ..
165: * .. External Subroutines ..
166: EXTERNAL DGTTS2, XERBLA
167: * ..
168: * .. Intrinsic Functions ..
169: INTRINSIC MAX, MIN
170: * ..
171: * .. Executable Statements ..
172: *
173: INFO = 0
174: NOTRAN = ( TRANS.EQ.'N' .OR. TRANS.EQ.'n' )
175: IF( .NOT.NOTRAN .AND. .NOT.( TRANS.EQ.'T' .OR. TRANS.EQ.
176: $ 't' ) .AND. .NOT.( TRANS.EQ.'C' .OR. TRANS.EQ.'c' ) ) THEN
177: INFO = -1
178: ELSE IF( N.LT.0 ) THEN
179: INFO = -2
180: ELSE IF( NRHS.LT.0 ) THEN
181: INFO = -3
182: ELSE IF( LDB.LT.MAX( N, 1 ) ) THEN
183: INFO = -10
184: END IF
185: IF( INFO.NE.0 ) THEN
186: CALL XERBLA( 'DGTTRS', -INFO )
187: RETURN
188: END IF
189: *
190: * Quick return if possible
191: *
192: IF( N.EQ.0 .OR. NRHS.EQ.0 )
193: $ RETURN
194: *
195: * Decode TRANS
196: *
197: IF( NOTRAN ) THEN
198: ITRANS = 0
199: ELSE
200: ITRANS = 1
201: END IF
202: *
203: * Determine the number of right-hand sides to solve at a time.
204: *
205: IF( NRHS.EQ.1 ) THEN
206: NB = 1
207: ELSE
208: NB = MAX( 1, ILAENV( 1, 'DGTTRS', TRANS, N, NRHS, -1, -1 ) )
209: END IF
210: *
211: IF( NB.GE.NRHS ) THEN
212: CALL DGTTS2( ITRANS, N, NRHS, DL, D, DU, DU2, IPIV, B, LDB )
213: ELSE
214: DO 10 J = 1, NRHS, NB
215: JB = MIN( NRHS-J+1, NB )
216: CALL DGTTS2( ITRANS, N, JB, DL, D, DU, DU2, IPIV, B( 1, J ),
217: $ LDB )
218: 10 CONTINUE
219: END IF
220: *
221: * End of DGTTRS
222: *
223: END
CVSweb interface <joel.bertrand@systella.fr>