LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine dppsvx ( character  FACT,
character  UPLO,
integer  N,
integer  NRHS,
double precision, dimension( * )  AP,
double precision, dimension( * )  AFP,
character  EQUED,
double precision, dimension( * )  S,
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, dimension( * )  IWORK,
integer  INFO 
)

DPPSVX computes the solution to system of linear equations A * X = B for OTHER matrices

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

Purpose:
 DPPSVX uses the Cholesky factorization A = U**T*U or A = L*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 matrix stored in
 packed format 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 = 'E', real scaling factors are computed to equilibrate
    the system:
       diag(S) * A * diag(S) * inv(diag(S)) * X = diag(S) * B
    Whether or not the system will be equilibrated depends on the
    scaling of the matrix A, but if equilibration is used, A is
    overwritten by diag(S)*A*diag(S) and B by diag(S)*B.

 2. If FACT = 'N' or 'E', the Cholesky decomposition is used to
    factor the matrix A (after equilibration if FACT = 'E') as
       A = U**T* U,  if UPLO = 'U', or
       A = L * L**T,  if UPLO = 'L',
    where U is an upper triangular matrix and L is a lower triangular
    matrix.

 3. 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.

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

 5. Iterative refinement is applied to improve the computed solution
    matrix and calculate error bounds and backward error estimates
    for it.

 6. If equilibration was used, the matrix X is premultiplied by
    diag(S) so that it solves the original system before
    equilibration.
Parameters
[in]FACT
          FACT is CHARACTER*1
          Specifies whether or not the factored form of the matrix A is
          supplied on entry, and if not, whether the matrix A should be
          equilibrated before it is factored.
          = 'F':  On entry, AFP contains the factored form of A.
                  If EQUED = 'Y', the matrix A has been equilibrated
                  with scaling factors given by S.  AP and AFP will not
                  be modified.
          = 'N':  The matrix A will be copied to AFP and factored.
          = 'E':  The matrix A will be equilibrated if necessary, then
                  copied to AFP and factored.
[in]UPLO
          UPLO is CHARACTER*1
          = 'U':  Upper triangle of A is stored;
          = 'L':  Lower triangle of A is stored.
[in]N
          N is INTEGER
          The number of linear equations, i.e., 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,out]AP
          AP is DOUBLE PRECISION array, dimension (N*(N+1)/2)
          On entry, the upper or lower triangle of the symmetric matrix
          A, packed columnwise in a linear array, except if FACT = 'F'
          and EQUED = 'Y', then A must contain the equilibrated matrix
          diag(S)*A*diag(S).  The j-th column of A is stored in the
          array AP as follows:
          if UPLO = 'U', AP(i + (j-1)*j/2) = A(i,j) for 1<=i<=j;
          if UPLO = 'L', AP(i + (j-1)*(2n-j)/2) = A(i,j) for j<=i<=n.
          See below for further details.  A is not modified if
          FACT = 'F' or 'N', or if FACT = 'E' and EQUED = 'N' on exit.

          On exit, if FACT = 'E' and EQUED = 'Y', A is overwritten by
          diag(S)*A*diag(S).
[in,out]AFP
          AFP is DOUBLE PRECISION array, dimension
                            (N*(N+1)/2)
          If FACT = 'F', then AFP is an input argument and on entry
          contains the triangular factor U or L from the Cholesky
          factorization A = U**T*U or A = L*L**T, in the same storage
          format as A.  If EQUED .ne. 'N', then AFP is the factored
          form of the equilibrated matrix A.

          If FACT = 'N', then AFP is an output argument and on exit
          returns the triangular factor U or L from the Cholesky
          factorization A = U**T * U or A = L * L**T of the original
          matrix A.

          If FACT = 'E', then AFP is an output argument and on exit
          returns the triangular factor U or L from the Cholesky
          factorization A = U**T * U or A = L * L**T of the equilibrated
          matrix A (see the description of AP for the form of the
          equilibrated matrix).
[in,out]EQUED
          EQUED is CHARACTER*1
          Specifies the form of equilibration that was done.
          = 'N':  No equilibration (always true if FACT = 'N').
          = 'Y':  Equilibration was done, i.e., A has been replaced by
                  diag(S) * A * diag(S).
          EQUED is an input argument if FACT = 'F'; otherwise, it is an
          output argument.
[in,out]S
          S is DOUBLE PRECISION array, dimension (N)
          The scale factors for A; not accessed if EQUED = 'N'.  S is
          an input argument if FACT = 'F'; otherwise, S is an output
          argument.  If FACT = 'F' and EQUED = 'Y', each element of S
          must be positive.
[in,out]B
          B is DOUBLE PRECISION array, dimension (LDB,NRHS)
          On entry, the N-by-NRHS right hand side matrix B.
          On exit, if EQUED = 'N', B is not modified; if EQUED = 'Y',
          B is overwritten by diag(S) * 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 or INFO = N+1, the N-by-NRHS solution matrix X to
          the original system of equations.  Note that if EQUED = 'Y',
          A and B are modified on exit, and the solution to the
          equilibrated system is inv(diag(S))*X.
[in]LDX
          LDX is INTEGER
          The leading dimension of the array X.  LDX >= max(1,N).
[out]RCOND
          RCOND is DOUBLE PRECISION
          The estimate of the reciprocal condition number of the matrix
          A after equilibration (if done).  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 estimated 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).  The estimate is as reliable as
          the estimate for RCOND, and is almost always a slight
          overestimate of the true error.
[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 (3*N)
[out]IWORK
          IWORK is INTEGER array, dimension (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
April 2012
Further Details:
  The packed storage scheme is illustrated by the following example
  when N = 4, UPLO = 'U':

  Two-dimensional storage of the symmetric matrix A:

     a11 a12 a13 a14
         a22 a23 a24
             a33 a34     (aij = conjg(aji))
                 a44

  Packed storage of the upper triangle of A:

  AP = [ a11, a12, a22, a13, a23, a33, a14, a24, a34, a44 ]

Definition at line 314 of file dppsvx.f.

314 *
315 * -- LAPACK driver routine (version 3.4.1) --
316 * -- LAPACK is a software package provided by Univ. of Tennessee, --
317 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
318 * April 2012
319 *
320 * .. Scalar Arguments ..
321  CHARACTER equed, fact, uplo
322  INTEGER info, ldb, ldx, n, nrhs
323  DOUBLE PRECISION rcond
324 * ..
325 * .. Array Arguments ..
326  INTEGER iwork( * )
327  DOUBLE PRECISION afp( * ), ap( * ), b( ldb, * ), berr( * ),
328  $ ferr( * ), s( * ), work( * ), x( ldx, * )
329 * ..
330 *
331 * =====================================================================
332 *
333 * .. Parameters ..
334  DOUBLE PRECISION zero, one
335  parameter ( zero = 0.0d+0, one = 1.0d+0 )
336 * ..
337 * .. Local Scalars ..
338  LOGICAL equil, nofact, rcequ
339  INTEGER i, infequ, j
340  DOUBLE PRECISION amax, anorm, bignum, scond, smax, smin, smlnum
341 * ..
342 * .. External Functions ..
343  LOGICAL lsame
344  DOUBLE PRECISION dlamch, dlansp
345  EXTERNAL lsame, dlamch, dlansp
346 * ..
347 * .. External Subroutines ..
348  EXTERNAL dcopy, dlacpy, dlaqsp, dppcon, dppequ, dpprfs,
349  $ dpptrf, dpptrs, xerbla
350 * ..
351 * .. Intrinsic Functions ..
352  INTRINSIC max, min
353 * ..
354 * .. Executable Statements ..
355 *
356  info = 0
357  nofact = lsame( fact, 'N' )
358  equil = lsame( fact, 'E' )
359  IF( nofact .OR. equil ) THEN
360  equed = 'N'
361  rcequ = .false.
362  ELSE
363  rcequ = lsame( equed, 'Y' )
364  smlnum = dlamch( 'Safe minimum' )
365  bignum = one / smlnum
366  END IF
367 *
368 * Test the input parameters.
369 *
370  IF( .NOT.nofact .AND. .NOT.equil .AND. .NOT.lsame( fact, 'F' ) )
371  $ THEN
372  info = -1
373  ELSE IF( .NOT.lsame( uplo, 'U' ) .AND. .NOT.lsame( uplo, 'L' ) )
374  $ THEN
375  info = -2
376  ELSE IF( n.LT.0 ) THEN
377  info = -3
378  ELSE IF( nrhs.LT.0 ) THEN
379  info = -4
380  ELSE IF( lsame( fact, 'F' ) .AND. .NOT.
381  $ ( rcequ .OR. lsame( equed, 'N' ) ) ) THEN
382  info = -7
383  ELSE
384  IF( rcequ ) THEN
385  smin = bignum
386  smax = zero
387  DO 10 j = 1, n
388  smin = min( smin, s( j ) )
389  smax = max( smax, s( j ) )
390  10 CONTINUE
391  IF( smin.LE.zero ) THEN
392  info = -8
393  ELSE IF( n.GT.0 ) THEN
394  scond = max( smin, smlnum ) / min( smax, bignum )
395  ELSE
396  scond = one
397  END IF
398  END IF
399  IF( info.EQ.0 ) THEN
400  IF( ldb.LT.max( 1, n ) ) THEN
401  info = -10
402  ELSE IF( ldx.LT.max( 1, n ) ) THEN
403  info = -12
404  END IF
405  END IF
406  END IF
407 *
408  IF( info.NE.0 ) THEN
409  CALL xerbla( 'DPPSVX', -info )
410  RETURN
411  END IF
412 *
413  IF( equil ) THEN
414 *
415 * Compute row and column scalings to equilibrate the matrix A.
416 *
417  CALL dppequ( uplo, n, ap, s, scond, amax, infequ )
418  IF( infequ.EQ.0 ) THEN
419 *
420 * Equilibrate the matrix.
421 *
422  CALL dlaqsp( uplo, n, ap, s, scond, amax, equed )
423  rcequ = lsame( equed, 'Y' )
424  END IF
425  END IF
426 *
427 * Scale the right-hand side.
428 *
429  IF( rcequ ) THEN
430  DO 30 j = 1, nrhs
431  DO 20 i = 1, n
432  b( i, j ) = s( i )*b( i, j )
433  20 CONTINUE
434  30 CONTINUE
435  END IF
436 *
437  IF( nofact .OR. equil ) THEN
438 *
439 * Compute the Cholesky factorization A = U**T * U or A = L * L**T.
440 *
441  CALL dcopy( n*( n+1 ) / 2, ap, 1, afp, 1 )
442  CALL dpptrf( uplo, n, afp, info )
443 *
444 * Return if INFO is non-zero.
445 *
446  IF( info.GT.0 )THEN
447  rcond = zero
448  RETURN
449  END IF
450  END IF
451 *
452 * Compute the norm of the matrix A.
453 *
454  anorm = dlansp( 'I', uplo, n, ap, work )
455 *
456 * Compute the reciprocal of the condition number of A.
457 *
458  CALL dppcon( uplo, n, afp, anorm, rcond, work, iwork, info )
459 *
460 * Compute the solution matrix X.
461 *
462  CALL dlacpy( 'Full', n, nrhs, b, ldb, x, ldx )
463  CALL dpptrs( uplo, n, nrhs, afp, x, ldx, info )
464 *
465 * Use iterative refinement to improve the computed solution and
466 * compute error bounds and backward error estimates for it.
467 *
468  CALL dpprfs( uplo, n, nrhs, ap, afp, b, ldb, x, ldx, ferr, berr,
469  $ work, iwork, info )
470 *
471 * Transform the solution matrix X to a solution of the original
472 * system.
473 *
474  IF( rcequ ) THEN
475  DO 50 j = 1, nrhs
476  DO 40 i = 1, n
477  x( i, j ) = s( i )*x( i, j )
478  40 CONTINUE
479  50 CONTINUE
480  DO 60 j = 1, nrhs
481  ferr( j ) = ferr( j ) / scond
482  60 CONTINUE
483  END IF
484 *
485 * Set INFO = N+1 if the matrix is singular to working precision.
486 *
487  IF( rcond.LT.dlamch( 'Epsilon' ) )
488  $ info = n + 1
489 *
490  RETURN
491 *
492 * End of DPPSVX
493 *
subroutine dpptrs(UPLO, N, NRHS, AP, B, LDB, INFO)
DPPTRS
Definition: dpptrs.f:110
subroutine dcopy(N, DX, INCX, DY, INCY)
DCOPY
Definition: dcopy.f:53
double precision function dlamch(CMACH)
DLAMCH
Definition: dlamch.f:65
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 dpprfs(UPLO, N, NRHS, AP, AFP, B, LDB, X, LDX, FERR, BERR, WORK, IWORK, INFO)
DPPRFS
Definition: dpprfs.f:173
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine dpptrf(UPLO, N, AP, INFO)
DPPTRF
Definition: dpptrf.f:121
subroutine dppequ(UPLO, N, AP, S, SCOND, AMAX, INFO)
DPPEQU
Definition: dppequ.f:118
subroutine dppcon(UPLO, N, AP, ANORM, RCOND, WORK, IWORK, INFO)
DPPCON
Definition: dppcon.f:120
subroutine dlaqsp(UPLO, N, AP, S, SCOND, AMAX, EQUED)
DLAQSP scales a symmetric/Hermitian matrix in packed storage, using scaling factors computed by sppeq...
Definition: dlaqsp.f:127
double precision function dlansp(NORM, UPLO, N, AP, WORK)
DLANSP returns the value of the 1-norm, or the Frobenius norm, or the infinity norm, or the element of largest absolute value of a symmetric matrix supplied in packed form.
Definition: dlansp.f:116
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55

Here is the call graph for this function:

Here is the caller graph for this function: