LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches

◆ dptsvx()

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 principal minor of order i is not positive,
    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 principal minor of order i of A
                       is not positive, 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.

Definition at line 226 of file dptsvx.f.

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