LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine dptsvx ( character  FACT,
integer  N,
integer  NRHS,
double precision, dimension( * )  D,
double precision, dimension( * )  E,
double precision, dimension( * )  DF,
double precision, dimension( * )  EF,
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  INFO 
)

DPTSVX computes the solution to system of linear equations A * X = B for PT matrices

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

Purpose:
 DPTSVX uses the factorization A = L*D*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 tridiagonal 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 = 'N', the matrix A is factored as A = L*D*L**T, where L
    is a unit lower bidiagonal matrix and D is diagonal.  The
    factorization can also be regarded as having the form
    A = U**T*D*U.

 2. If the leading i-by-i principal minor is not positive definite,
    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.

 3. The system of equations is solved for X using the factored form
    of A.

 4. Iterative refinement is applied to improve the computed solution
    matrix and calculate error bounds and backward error estimates
    for it.
Parameters
[in]FACT
          FACT is CHARACTER*1
          Specifies whether or not the factored form of A has been
          supplied on entry.
          = 'F':  On entry, DF and EF contain the factored form of A.
                  D, E, DF, and EF will not be modified.
          = 'N':  The matrix A will be copied to DF and EF and
                  factored.
[in]N
          N is INTEGER
          The order of the matrix A.  N >= 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]D
          D is DOUBLE PRECISION array, dimension (N)
          The n diagonal elements of the tridiagonal matrix A.
[in]E
          E is DOUBLE PRECISION array, dimension (N-1)
          The (n-1) subdiagonal elements of the tridiagonal matrix A.
[in,out]DF
          DF is DOUBLE PRECISION array, dimension (N)
          If FACT = 'F', then DF is an input argument and on entry
          contains the n diagonal elements of the diagonal matrix D
          from the L*D*L**T factorization of A.
          If FACT = 'N', then DF is an output argument and on exit
          contains the n diagonal elements of the diagonal matrix D
          from the L*D*L**T factorization of A.
[in,out]EF
          EF is DOUBLE PRECISION array, dimension (N-1)
          If FACT = 'F', then EF is an input argument and on entry
          contains the (n-1) subdiagonal elements of the unit
          bidiagonal factor L from the L*D*L**T factorization of A.
          If FACT = 'N', then EF is an output argument and on exit
          contains the (n-1) subdiagonal elements of the unit
          bidiagonal factor L from the L*D*L**T factorization of A.
[in]B
          B is DOUBLE PRECISION array, dimension (LDB,NRHS)
          The N-by-NRHS right hand side matrix 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 of INFO = N+1, the N-by-NRHS solution matrix X.
[in]LDX
          LDX is INTEGER
          The leading dimension of the array X.  LDX >= max(1,N).
[out]RCOND
          RCOND is DOUBLE PRECISION
          The reciprocal condition number of the matrix A.  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 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).
[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 (2*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 minor of order i of A is
                       not positive definite, 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.
Date
September 2012

Definition at line 230 of file dptsvx.f.

230 *
231 * -- LAPACK driver routine (version 3.4.2) --
232 * -- LAPACK is a software package provided by Univ. of Tennessee, --
233 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
234 * September 2012
235 *
236 * .. Scalar Arguments ..
237  CHARACTER fact
238  INTEGER info, ldb, ldx, n, nrhs
239  DOUBLE PRECISION rcond
240 * ..
241 * .. Array Arguments ..
242  DOUBLE PRECISION b( ldb, * ), berr( * ), d( * ), df( * ),
243  $ e( * ), ef( * ), ferr( * ), work( * ),
244  $ x( ldx, * )
245 * ..
246 *
247 * =====================================================================
248 *
249 * .. Parameters ..
250  DOUBLE PRECISION zero
251  parameter ( zero = 0.0d+0 )
252 * ..
253 * .. Local Scalars ..
254  LOGICAL nofact
255  DOUBLE PRECISION anorm
256 * ..
257 * .. External Functions ..
258  LOGICAL lsame
259  DOUBLE PRECISION dlamch, dlanst
260  EXTERNAL lsame, dlamch, dlanst
261 * ..
262 * .. External Subroutines ..
263  EXTERNAL dcopy, dlacpy, dptcon, dptrfs, dpttrf, dpttrs,
264  $ xerbla
265 * ..
266 * .. Intrinsic Functions ..
267  INTRINSIC max
268 * ..
269 * .. Executable Statements ..
270 *
271 * Test the input parameters.
272 *
273  info = 0
274  nofact = lsame( fact, 'N' )
275  IF( .NOT.nofact .AND. .NOT.lsame( fact, 'F' ) ) THEN
276  info = -1
277  ELSE IF( n.LT.0 ) THEN
278  info = -2
279  ELSE IF( nrhs.LT.0 ) THEN
280  info = -3
281  ELSE IF( ldb.LT.max( 1, n ) ) THEN
282  info = -9
283  ELSE IF( ldx.LT.max( 1, n ) ) THEN
284  info = -11
285  END IF
286  IF( info.NE.0 ) THEN
287  CALL xerbla( 'DPTSVX', -info )
288  RETURN
289  END IF
290 *
291  IF( nofact ) THEN
292 *
293 * Compute the L*D*L**T (or U**T*D*U) factorization of A.
294 *
295  CALL dcopy( n, d, 1, df, 1 )
296  IF( n.GT.1 )
297  $ CALL dcopy( n-1, e, 1, ef, 1 )
298  CALL dpttrf( n, df, ef, info )
299 *
300 * Return if INFO is non-zero.
301 *
302  IF( info.GT.0 )THEN
303  rcond = zero
304  RETURN
305  END IF
306  END IF
307 *
308 * Compute the norm of the matrix A.
309 *
310  anorm = dlanst( '1', n, d, e )
311 *
312 * Compute the reciprocal of the condition number of A.
313 *
314  CALL dptcon( n, df, ef, anorm, rcond, work, info )
315 *
316 * Compute the solution vectors X.
317 *
318  CALL dlacpy( 'Full', n, nrhs, b, ldb, x, ldx )
319  CALL dpttrs( n, nrhs, df, ef, x, ldx, info )
320 *
321 * Use iterative refinement to improve the computed solutions and
322 * compute error bounds and backward error estimates for them.
323 *
324  CALL dptrfs( n, nrhs, d, e, df, ef, b, ldb, x, ldx, ferr, berr,
325  $ work, info )
326 *
327 * Set INFO = N+1 if the matrix is singular to working precision.
328 *
329  IF( rcond.LT.dlamch( 'Epsilon' ) )
330  $ info = n + 1
331 *
332  RETURN
333 *
334 * End of DPTSVX
335 *
subroutine dcopy(N, DX, INCX, DY, INCY)
DCOPY
Definition: dcopy.f:53
double precision function dlamch(CMACH)
DLAMCH
Definition: dlamch.f:65
subroutine dpttrf(N, D, E, INFO)
DPTTRF
Definition: dpttrf.f:93
subroutine dlacpy(UPLO, M, N, A, LDA, B, LDB)
DLACPY copies all or part of one two-dimensional array to another.
Definition: dlacpy.f:105
subroutine dptcon(N, D, E, ANORM, RCOND, WORK, INFO)
DPTCON
Definition: dptcon.f:120
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine dptrfs(N, NRHS, D, E, DF, EF, B, LDB, X, LDX, FERR, BERR, WORK, INFO)
DPTRFS
Definition: dptrfs.f:165
double precision function dlanst(NORM, N, D, E)
DLANST returns the value of the 1-norm, or the Frobenius norm, or the infinity norm, or the element of largest absolute value of a real symmetric tridiagonal matrix.
Definition: dlanst.f:102
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55
subroutine dpttrs(N, NRHS, D, E, B, LDB, INFO)
DPTTRS
Definition: dpttrs.f:111

Here is the call graph for this function:

Here is the caller graph for this function: