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

◆ dppsvx()

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

 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 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.
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 309 of file dppsvx.f.

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