.TH DGTSV l "15 June 2000" "LAPACK version 3.0" ")"
.SH NAME
DGTSV - solve the equation  A*X = B,
.SH SYNOPSIS
.TP 18
SUBROUTINE DGTSV(
N, NRHS, DL, D, DU, B, LDB, INFO )
.TP 18
.ti +4
INTEGER
INFO, LDB, N, NRHS
.TP 18
.ti +4
DOUBLE
PRECISION B( LDB, * ), D( * ), DL( * ), DU( * )
.SH PURPOSE
DGTSV solves the equation A*X = B, 
where A is an n by n tridiagonal matrix, by Gaussian elimination with
partial pivoting.
.br

Note that the equation  A'*X = B  may be solved by interchanging the
order of the arguments DU and DL.
.br

.SH ARGUMENTS
.TP 8
N       (input) INTEGER
The order of the matrix A.  N >= 0.
.TP 8
NRHS    (input) INTEGER
The number of right hand sides, i.e., the number of columns
of the matrix B.  NRHS >= 0.
.TP 8
DL      (input/output) DOUBLE PRECISION 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-2) elements of the
second super-diagonal of the upper triangular matrix U from
the LU factorization of A, in DL(1), ..., DL(n-2).
.TP 8
D       (input/output) DOUBLE PRECISION array, dimension (N)
On entry, D must contain the diagonal elements of A.

On exit, D is overwritten by the n diagonal elements of U.
.TP 8
DU      (input/output) DOUBLE PRECISION 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.
.TP 8
B       (input/output) DOUBLE PRECISION array, dimension (LDB,NRHS)
On entry, the N by NRHS matrix of right hand side matrix B.
On exit, if INFO = 0, the N by NRHS solution matrix X.
.TP 8
LDB     (input) INTEGER
The leading dimension of the array B.  LDB >= max(1,N).
.TP 8
INFO    (output) INTEGER
= 0: successful exit
.br
< 0: if INFO = -i, the i-th argument had an illegal value
.br
> 0: if INFO = i, U(i,i) is exactly zero, and the solution
has not been computed.  The factorization has not been
completed unless i = N.
