LAPACK  3.6.0
LAPACK: Linear Algebra PACKage
Collaboration diagram for double:

Functions

subroutine dposv (UPLO, N, NRHS, A, LDA, B, LDB, INFO)
  DPOSV computes the solution to system of linear equations A * X = B for PO matrices More...
 
subroutine dposvx (FACT, UPLO, N, NRHS, A, LDA, AF, LDAF, EQUED, S, B, LDB, X, LDX, RCOND, FERR, BERR, WORK, IWORK, INFO)
  DPOSVX computes the solution to system of linear equations A * X = B for PO matrices More...
 
subroutine dposvxx (FACT, UPLO, N, NRHS, A, LDA, AF, LDAF, EQUED, S, B, LDB, X, LDX, RCOND, RPVGRW, BERR, N_ERR_BNDS, ERR_BNDS_NORM, ERR_BNDS_COMP, NPARAMS, PARAMS, WORK, IWORK, INFO)
  DPOSVXX computes the solution to system of linear equations A * X = B for PO matrices More...
 
subroutine dsposv (UPLO, N, NRHS, A, LDA, B, LDB, X, LDX, WORK, SWORK, ITER, INFO)
  DSPOSV computes the solution to system of linear equations A * X = B for PO matrices More...
 

Detailed Description

This is the group of double solve driver functions for PO matrices

Function Documentation

subroutine dposv ( character  UPLO,
integer  N,
integer  NRHS,
double precision, dimension( lda, * )  A,
integer  LDA,
double precision, dimension( ldb, * )  B,
integer  LDB,
integer  INFO 
)

DPOSV computes the solution to system of linear equations A * X = B for PO matrices

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

Purpose:
 DPOSV computes the solution to a real system of linear equations
    A * X = B,
 where A is an N-by-N symmetric positive definite matrix and X and B
 are N-by-NRHS matrices.

 The Cholesky decomposition is used to factor A 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.  The factored form of A is then used to solve the system of
 equations A * X = B.
Parameters
[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 matrix B.  NRHS >= 0.
[in,out]A
          A is DOUBLE PRECISION array, dimension (LDA,N)
          On entry, the symmetric matrix A.  If UPLO = 'U', the leading
          N-by-N upper triangular part of A contains the upper
          triangular part of the matrix A, and the strictly lower
          triangular part of A is not referenced.  If UPLO = 'L', the
          leading N-by-N lower triangular part of A contains the lower
          triangular part of the matrix A, and the strictly upper
          triangular part of A is not referenced.

          On exit, if INFO = 0, the factor U or L from the Cholesky
          factorization A = U**T*U or A = L*L**T.
[in]LDA
          LDA is INTEGER
          The leading dimension of the array A.  LDA >= max(1,N).
[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 INFO = 0, the N-by-NRHS solution matrix X.
[in]LDB
          LDB is INTEGER
          The leading dimension of the array B.  LDB >= max(1,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, 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.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
November 2011

Definition at line 132 of file dposv.f.

132 *
133 * -- LAPACK driver routine (version 3.4.0) --
134 * -- LAPACK is a software package provided by Univ. of Tennessee, --
135 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
136 * November 2011
137 *
138 * .. Scalar Arguments ..
139  CHARACTER uplo
140  INTEGER info, lda, ldb, n, nrhs
141 * ..
142 * .. Array Arguments ..
143  DOUBLE PRECISION a( lda, * ), b( ldb, * )
144 * ..
145 *
146 * =====================================================================
147 *
148 * .. External Functions ..
149  LOGICAL lsame
150  EXTERNAL lsame
151 * ..
152 * .. External Subroutines ..
153  EXTERNAL dpotrf, dpotrs, xerbla
154 * ..
155 * .. Intrinsic Functions ..
156  INTRINSIC max
157 * ..
158 * .. Executable Statements ..
159 *
160 * Test the input parameters.
161 *
162  info = 0
163  IF( .NOT.lsame( uplo, 'U' ) .AND. .NOT.lsame( uplo, 'L' ) ) THEN
164  info = -1
165  ELSE IF( n.LT.0 ) THEN
166  info = -2
167  ELSE IF( nrhs.LT.0 ) THEN
168  info = -3
169  ELSE IF( lda.LT.max( 1, n ) ) THEN
170  info = -5
171  ELSE IF( ldb.LT.max( 1, n ) ) THEN
172  info = -7
173  END IF
174  IF( info.NE.0 ) THEN
175  CALL xerbla( 'DPOSV ', -info )
176  RETURN
177  END IF
178 *
179 * Compute the Cholesky factorization A = U**T*U or A = L*L**T.
180 *
181  CALL dpotrf( uplo, n, a, lda, info )
182  IF( info.EQ.0 ) THEN
183 *
184 * Solve the system A*X = B, overwriting B with X.
185 *
186  CALL dpotrs( uplo, n, nrhs, a, lda, b, ldb, info )
187 *
188  END IF
189  RETURN
190 *
191 * End of DPOSV
192 *
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55
subroutine dpotrf(UPLO, N, A, LDA, INFO)
DPOTRF
Definition: dpotrf.f:109
subroutine dpotrs(UPLO, N, NRHS, A, LDA, B, LDB, INFO)
DPOTRS
Definition: dpotrs.f:112

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine dposvx ( character  FACT,
character  UPLO,
integer  N,
integer  NRHS,
double precision, dimension( lda, * )  A,
integer  LDA,
double precision, dimension( ldaf, * )  AF,
integer  LDAF,
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 
)

DPOSVX computes the solution to system of linear equations A * X = B for PO matrices

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

Purpose:
 DPOSVX 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 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, AF contains the factored form of A.
                  If EQUED = 'Y', the matrix A has been equilibrated
                  with scaling factors given by S.  A and AF will not
                  be modified.
          = 'N':  The matrix A will be copied to AF and factored.
          = 'E':  The matrix A will be equilibrated if necessary, then
                  copied to AF 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]A
          A is DOUBLE PRECISION array, dimension (LDA,N)
          On entry, the symmetric matrix A, except if FACT = 'F' and
          EQUED = 'Y', then A must contain the equilibrated matrix
          diag(S)*A*diag(S).  If UPLO = 'U', the leading
          N-by-N upper triangular part of A contains the upper
          triangular part of the matrix A, and the strictly lower
          triangular part of A is not referenced.  If UPLO = 'L', the
          leading N-by-N lower triangular part of A contains the lower
          triangular part of the matrix A, and the strictly upper
          triangular part of A is not referenced.  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]LDA
          LDA is INTEGER
          The leading dimension of the array A.  LDA >= max(1,N).
[in,out]AF
          AF is DOUBLE PRECISION array, dimension (LDAF,N)
          If FACT = 'F', then AF 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 AF is the factored form
          of the equilibrated matrix diag(S)*A*diag(S).

          If FACT = 'N', then AF 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 AF 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 A for the form of the
          equilibrated matrix).
[in]LDAF
          LDAF is INTEGER
          The leading dimension of the array AF.  LDAF >= max(1,N).
[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

Definition at line 309 of file dposvx.f.

309 *
310 * -- LAPACK driver routine (version 3.4.1) --
311 * -- LAPACK is a software package provided by Univ. of Tennessee, --
312 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
313 * April 2012
314 *
315 * .. Scalar Arguments ..
316  CHARACTER equed, fact, uplo
317  INTEGER info, lda, ldaf, ldb, ldx, n, nrhs
318  DOUBLE PRECISION rcond
319 * ..
320 * .. Array Arguments ..
321  INTEGER iwork( * )
322  DOUBLE PRECISION a( lda, * ), af( ldaf, * ), b( ldb, * ),
323  $ berr( * ), ferr( * ), s( * ), work( * ),
324  $ 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, dlansy
341  EXTERNAL lsame, dlamch, dlansy
342 * ..
343 * .. External Subroutines ..
344  EXTERNAL dlacpy, dlaqsy, dpocon, dpoequ, dporfs, dpotrf,
345  $ dpotrs, xerbla
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( lda.LT.max( 1, n ) ) THEN
377  info = -6
378  ELSE IF( ldaf.LT.max( 1, n ) ) THEN
379  info = -8
380  ELSE IF( lsame( fact, 'F' ) .AND. .NOT.
381  $ ( rcequ .OR. lsame( equed, 'N' ) ) ) THEN
382  info = -9
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 = -10
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 = -12
402  ELSE IF( ldx.LT.max( 1, n ) ) THEN
403  info = -14
404  END IF
405  END IF
406  END IF
407 *
408  IF( info.NE.0 ) THEN
409  CALL xerbla( 'DPOSVX', -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 dpoequ( n, a, lda, s, scond, amax, infequ )
418  IF( infequ.EQ.0 ) THEN
419 *
420 * Equilibrate the matrix.
421 *
422  CALL dlaqsy( uplo, n, a, lda, 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 dlacpy( uplo, n, n, a, lda, af, ldaf )
442  CALL dpotrf( uplo, n, af, ldaf, 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 = dlansy( '1', uplo, n, a, lda, work )
455 *
456 * Compute the reciprocal of the condition number of A.
457 *
458  CALL dpocon( uplo, n, af, ldaf, 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 dpotrs( uplo, n, nrhs, af, ldaf, 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 dporfs( uplo, n, nrhs, a, lda, af, ldaf, b, ldb, x, ldx,
469  $ ferr, berr, 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 DPOSVX
493 *
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 dlaqsy(UPLO, N, A, LDA, S, SCOND, AMAX, EQUED)
DLAQSY scales a symmetric/Hermitian matrix, using scaling factors computed by spoequ.
Definition: dlaqsy.f:135
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
double precision function dlamch(CMACH)
DLAMCH
Definition: dlamch.f:65
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55
subroutine dpoequ(N, A, LDA, S, SCOND, AMAX, INFO)
DPOEQU
Definition: dpoequ.f:114
subroutine dpotrf(UPLO, N, A, LDA, INFO)
DPOTRF
Definition: dpotrf.f:109
subroutine dpotrs(UPLO, N, NRHS, A, LDA, B, LDB, INFO)
DPOTRS
Definition: dpotrs.f:112
subroutine dpocon(UPLO, N, A, LDA, ANORM, RCOND, WORK, IWORK, INFO)
DPOCON
Definition: dpocon.f:123
subroutine dporfs(UPLO, N, NRHS, A, LDA, AF, LDAF, B, LDB, X, LDX, FERR, BERR, WORK, IWORK, INFO)
DPORFS
Definition: dporfs.f:185
double precision function dlansy(NORM, UPLO, N, A, LDA, WORK)
DLANSY 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 matrix.
Definition: dlansy.f:124

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine dposvxx ( character  FACT,
character  UPLO,
integer  N,
integer  NRHS,
double precision, dimension( lda, * )  A,
integer  LDA,
double precision, dimension( ldaf, * )  AF,
integer  LDAF,
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  RPVGRW,
double precision, dimension( * )  BERR,
integer  N_ERR_BNDS,
double precision, dimension( nrhs, * )  ERR_BNDS_NORM,
double precision, dimension( nrhs, * )  ERR_BNDS_COMP,
integer  NPARAMS,
double precision, dimension( * )  PARAMS,
double precision, dimension( * )  WORK,
integer, dimension( * )  IWORK,
integer  INFO 
)

DPOSVXX computes the solution to system of linear equations A * X = B for PO matrices

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

Purpose:
    DPOSVXX uses the Cholesky factorization A = U**T*U or A = L*L**T
    to compute the solution to a double precision system of linear equations
    A * X = B, where A is an N-by-N symmetric positive definite matrix
    and X and B are N-by-NRHS matrices.

    If requested, both normwise and maximum componentwise error bounds
    are returned. DPOSVXX will return a solution with a tiny
    guaranteed error (O(eps) where eps is the working machine
    precision) unless the matrix is very ill-conditioned, in which
    case a warning is returned. Relevant condition numbers also are
    calculated and returned.

    DPOSVXX accepts user-provided factorizations and equilibration
    factors; see the definitions of the FACT and EQUED options.
    Solving with refinement and using a factorization from a previous
    DPOSVXX call will also produce a solution with either O(eps)
    errors or warnings, but we cannot make that claim for general
    user-provided factorizations and equilibration factors if they
    differ from what DPOSVXX would itself produce.
Description:
    The following steps are performed:

    1. If FACT = 'E', double precision 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 (see argument RCOND).  If the reciprocal of the condition number
    is less than machine precision, 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. By default (unless PARAMS(LA_LINRX_ITREF_I) is set to zero),
    the routine will use iterative refinement to try to get a small
    error and error bounds.  Refinement calculates the residual to at
    least twice the working precision.

    6. If equilibration was used, the matrix X is premultiplied by
    diag(S) so that it solves the original system before
    equilibration.
     Some optional parameters are bundled in the PARAMS array.  These
     settings determine how refinement is performed, but often the
     defaults are acceptable.  If the defaults are acceptable, users
     can pass NPARAMS = 0 which prevents the source code from accessing
     the PARAMS argument.
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, AF contains the factored form of A.
               If EQUED is not 'N', the matrix A has been
               equilibrated with scaling factors given by S.
               A and AF are not modified.
       = 'N':  The matrix A will be copied to AF and factored.
       = 'E':  The matrix A will be equilibrated if necessary, then
               copied to AF 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]A
          A is DOUBLE PRECISION array, dimension (LDA,N)
     On entry, the symmetric matrix A, except if FACT = 'F' and EQUED =
     'Y', then A must contain the equilibrated matrix
     diag(S)*A*diag(S).  If UPLO = 'U', the leading N-by-N upper
     triangular part of A contains the upper triangular part of the
     matrix A, and the strictly lower triangular part of A is not
     referenced.  If UPLO = 'L', the leading N-by-N lower triangular
     part of A contains the lower triangular part of the matrix A, and
     the strictly upper triangular part of A is not referenced.  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]LDA
          LDA is INTEGER
     The leading dimension of the array A.  LDA >= max(1,N).
[in,out]AF
          AF is DOUBLE PRECISION array, dimension (LDAF,N)
     If FACT = 'F', then AF 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 AF is the factored
     form of the equilibrated matrix diag(S)*A*diag(S).

     If FACT = 'N', then AF 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 AF 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 A for the form of the
     equilibrated matrix).
[in]LDAF
          LDAF is INTEGER
     The leading dimension of the array AF.  LDAF >= max(1,N).
[in,out]EQUED
          EQUED is CHARACTER*1
     Specifies the form of equilibration that was done.
       = 'N':  No equilibration (always true if FACT = 'N').
       = 'Y':  Both row and column equilibration, 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 row scale factors for A.  If EQUED = 'Y', A is multiplied on
     the left and right by diag(S).  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.  If S is output, each
     element of S is a power of the radix. If S is input, each element
     of S should be a power of the radix to ensure a reliable solution
     and error estimates. Scaling by powers of the radix does not cause
     rounding errors unless the result underflows or overflows.
     Rounding errors during scaling lead to refining with a matrix that
     is not equivalent to the input matrix, producing error estimates
     that may not be reliable.
[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, the N-by-NRHS solution matrix X to the original
     system of equations.  Note that A and B are modified on exit if
     EQUED .ne. 'N', 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
     Reciprocal scaled condition number.  This is an estimate of the
     reciprocal Skeel condition number of the matrix A after
     equilibration (if done).  If this is less than the machine
     precision (in particular, if it is zero), the matrix is singular
     to working precision.  Note that the error may still be small even
     if this number is very small and the matrix appears ill-
     conditioned.
[out]RPVGRW
          RPVGRW is DOUBLE PRECISION
     Reciprocal pivot growth.  On exit, this contains the reciprocal
     pivot growth factor norm(A)/norm(U). The "max absolute element"
     norm is used.  If this is much less than 1, then the stability of
     the LU factorization of the (equilibrated) matrix A could be poor.
     This also means that the solution X, estimated condition numbers,
     and error bounds could be unreliable. If factorization fails with
     0<INFO<=N, then this contains the reciprocal pivot growth factor
     for the leading INFO columns of A.
[out]BERR
          BERR is DOUBLE PRECISION array, dimension (NRHS)
     Componentwise relative backward error.  This is 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).
[in]N_ERR_BNDS
          N_ERR_BNDS is INTEGER
     Number of error bounds to return for each right hand side
     and each type (normwise or componentwise).  See ERR_BNDS_NORM and
     ERR_BNDS_COMP below.
[out]ERR_BNDS_NORM
          ERR_BNDS_NORM is DOUBLE PRECISION array, dimension (NRHS, N_ERR_BNDS)
     For each right-hand side, this array contains information about
     various error bounds and condition numbers corresponding to the
     normwise relative error, which is defined as follows:

     Normwise relative error in the ith solution vector:
             max_j (abs(XTRUE(j,i) - X(j,i)))
            ------------------------------
                  max_j abs(X(j,i))

     The array is indexed by the type of error information as described
     below. There currently are up to three pieces of information
     returned.

     The first index in ERR_BNDS_NORM(i,:) corresponds to the ith
     right-hand side.

     The second index in ERR_BNDS_NORM(:,err) contains the following
     three fields:
     err = 1 "Trust/don't trust" boolean. Trust the answer if the
              reciprocal condition number is less than the threshold
              sqrt(n) * dlamch('Epsilon').

     err = 2 "Guaranteed" error bound: The estimated forward error,
              almost certainly within a factor of 10 of the true error
              so long as the next entry is greater than the threshold
              sqrt(n) * dlamch('Epsilon'). This error bound should only
              be trusted if the previous boolean is true.

     err = 3  Reciprocal condition number: Estimated normwise
              reciprocal condition number.  Compared with the threshold
              sqrt(n) * dlamch('Epsilon') to determine if the error
              estimate is "guaranteed". These reciprocal condition
              numbers are 1 / (norm(Z^{-1},inf) * norm(Z,inf)) for some
              appropriately scaled matrix Z.
              Let Z = S*A, where S scales each row by a power of the
              radix so all absolute row sums of Z are approximately 1.

     See Lapack Working Note 165 for further details and extra
     cautions.
[out]ERR_BNDS_COMP
          ERR_BNDS_COMP is DOUBLE PRECISION array, dimension (NRHS, N_ERR_BNDS)
     For each right-hand side, this array contains information about
     various error bounds and condition numbers corresponding to the
     componentwise relative error, which is defined as follows:

     Componentwise relative error in the ith solution vector:
                    abs(XTRUE(j,i) - X(j,i))
             max_j ----------------------
                         abs(X(j,i))

     The array is indexed by the right-hand side i (on which the
     componentwise relative error depends), and the type of error
     information as described below. There currently are up to three
     pieces of information returned for each right-hand side. If
     componentwise accuracy is not requested (PARAMS(3) = 0.0), then
     ERR_BNDS_COMP is not accessed.  If N_ERR_BNDS .LT. 3, then at most
     the first (:,N_ERR_BNDS) entries are returned.

     The first index in ERR_BNDS_COMP(i,:) corresponds to the ith
     right-hand side.

     The second index in ERR_BNDS_COMP(:,err) contains the following
     three fields:
     err = 1 "Trust/don't trust" boolean. Trust the answer if the
              reciprocal condition number is less than the threshold
              sqrt(n) * dlamch('Epsilon').

     err = 2 "Guaranteed" error bound: The estimated forward error,
              almost certainly within a factor of 10 of the true error
              so long as the next entry is greater than the threshold
              sqrt(n) * dlamch('Epsilon'). This error bound should only
              be trusted if the previous boolean is true.

     err = 3  Reciprocal condition number: Estimated componentwise
              reciprocal condition number.  Compared with the threshold
              sqrt(n) * dlamch('Epsilon') to determine if the error
              estimate is "guaranteed". These reciprocal condition
              numbers are 1 / (norm(Z^{-1},inf) * norm(Z,inf)) for some
              appropriately scaled matrix Z.
              Let Z = S*(A*diag(x)), where x is the solution for the
              current right-hand side and S scales each row of
              A*diag(x) by a power of the radix so all absolute row
              sums of Z are approximately 1.

     See Lapack Working Note 165 for further details and extra
     cautions.
[in]NPARAMS
          NPARAMS is INTEGER
     Specifies the number of parameters set in PARAMS.  If .LE. 0, the
     PARAMS array is never referenced and default values are used.
[in,out]PARAMS
          PARAMS is DOUBLE PRECISION array, dimension NPARAMS
     Specifies algorithm parameters.  If an entry is .LT. 0.0, then
     that entry will be filled with default value used for that
     parameter.  Only positions up to NPARAMS are accessed; defaults
     are used for higher-numbered parameters.

       PARAMS(LA_LINRX_ITREF_I = 1) : Whether to perform iterative
            refinement or not.
         Default: 1.0D+0
            = 0.0 : No refinement is performed, and no error bounds are
                    computed.
            = 1.0 : Use the extra-precise refinement algorithm.
              (other values are reserved for future use)

       PARAMS(LA_LINRX_ITHRESH_I = 2) : Maximum number of residual
            computations allowed for refinement.
         Default: 10
         Aggressive: Set to 100 to permit convergence using approximate
                     factorizations or factorizations other than LU. If
                     the factorization uses a technique other than
                     Gaussian elimination, the guarantees in
                     err_bnds_norm and err_bnds_comp may no longer be
                     trustworthy.

       PARAMS(LA_LINRX_CWISE_I = 3) : Flag determining if the code
            will attempt to find a solution with small componentwise
            relative error in the double-precision algorithm.  Positive
            is true, 0.0 is false.
         Default: 1.0 (attempt componentwise convergence)
[out]WORK
          WORK is DOUBLE PRECISION array, dimension (4*N)
[out]IWORK
          IWORK is INTEGER array, dimension (N)
[out]INFO
          INFO is INTEGER
       = 0:  Successful exit. The solution to every right-hand side is
         guaranteed.
       < 0:  If INFO = -i, the i-th argument had an illegal value
       > 0 and <= N:  U(INFO,INFO) is exactly zero.  The factorization
         has been completed, but the factor U is exactly singular, so
         the solution and error bounds could not be computed. RCOND = 0
         is returned.
       = N+J: The solution corresponding to the Jth right-hand side is
         not guaranteed. The solutions corresponding to other right-
         hand sides K with K > J may not be guaranteed as well, but
         only the first such right-hand side is reported. If a small
         componentwise error is not requested (PARAMS(3) = 0.0) then
         the Jth right-hand side is the first with a normwise error
         bound that is not guaranteed (the smallest J such
         that ERR_BNDS_NORM(J,1) = 0.0). By default (PARAMS(3) = 1.0)
         the Jth right-hand side is the first with either a normwise or
         componentwise error bound that is not guaranteed (the smallest
         J such that either ERR_BNDS_NORM(J,1) = 0.0 or
         ERR_BNDS_COMP(J,1) = 0.0). See the definition of
         ERR_BNDS_NORM(:,1) and ERR_BNDS_COMP(:,1). To get information
         about all of the right-hand sides check ERR_BNDS_NORM or
         ERR_BNDS_COMP.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
April 2012

Definition at line 496 of file dposvxx.f.

496 *
497 * -- LAPACK driver routine (version 3.4.1) --
498 * -- LAPACK is a software package provided by Univ. of Tennessee, --
499 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
500 * April 2012
501 *
502 * .. Scalar Arguments ..
503  CHARACTER equed, fact, uplo
504  INTEGER info, lda, ldaf, ldb, ldx, n, nrhs, nparams,
505  $ n_err_bnds
506  DOUBLE PRECISION rcond, rpvgrw
507 * ..
508 * .. Array Arguments ..
509  INTEGER iwork( * )
510  DOUBLE PRECISION a( lda, * ), af( ldaf, * ), b( ldb, * ),
511  $ x( ldx, * ), work( * )
512  DOUBLE PRECISION s( * ), params( * ), berr( * ),
513  $ err_bnds_norm( nrhs, * ),
514  $ err_bnds_comp( nrhs, * )
515 * ..
516 *
517 * ==================================================================
518 *
519 * .. Parameters ..
520  DOUBLE PRECISION zero, one
521  parameter( zero = 0.0d+0, one = 1.0d+0 )
522  INTEGER final_nrm_err_i, final_cmp_err_i, berr_i
523  INTEGER rcond_i, nrm_rcond_i, nrm_err_i, cmp_rcond_i
524  INTEGER cmp_err_i, piv_growth_i
525  parameter( final_nrm_err_i = 1, final_cmp_err_i = 2,
526  $ berr_i = 3 )
527  parameter( rcond_i = 4, nrm_rcond_i = 5, nrm_err_i = 6 )
528  parameter( cmp_rcond_i = 7, cmp_err_i = 8,
529  $ piv_growth_i = 9 )
530 * ..
531 * .. Local Scalars ..
532  LOGICAL equil, nofact, rcequ
533  INTEGER infequ, j
534  DOUBLE PRECISION amax, bignum, smin, smax,
535  $ scond, smlnum
536 * ..
537 * .. External Functions ..
538  EXTERNAL lsame, dlamch, dla_porpvgrw
539  LOGICAL lsame
540  DOUBLE PRECISION dlamch, dla_porpvgrw
541 * ..
542 * .. External Subroutines ..
543  EXTERNAL dpoequb, dpotrf, dpotrs, dlacpy, dlaqsy,
545 * ..
546 * .. Intrinsic Functions ..
547  INTRINSIC max, min
548 * ..
549 * .. Executable Statements ..
550 *
551  info = 0
552  nofact = lsame( fact, 'N' )
553  equil = lsame( fact, 'E' )
554  smlnum = dlamch( 'Safe minimum' )
555  bignum = one / smlnum
556  IF( nofact .OR. equil ) THEN
557  equed = 'N'
558  rcequ = .false.
559  ELSE
560  rcequ = lsame( equed, 'Y' )
561  ENDIF
562 *
563 * Default is failure. If an input parameter is wrong or
564 * factorization fails, make everything look horrible. Only the
565 * pivot growth is set here, the rest is initialized in DPORFSX.
566 *
567  rpvgrw = zero
568 *
569 * Test the input parameters. PARAMS is not tested until DPORFSX.
570 *
571  IF( .NOT.nofact .AND. .NOT.equil .AND. .NOT.
572  $ lsame( fact, 'F' ) ) THEN
573  info = -1
574  ELSE IF( .NOT.lsame( uplo, 'U' ) .AND.
575  $ .NOT.lsame( uplo, 'L' ) ) THEN
576  info = -2
577  ELSE IF( n.LT.0 ) THEN
578  info = -3
579  ELSE IF( nrhs.LT.0 ) THEN
580  info = -4
581  ELSE IF( lda.LT.max( 1, n ) ) THEN
582  info = -6
583  ELSE IF( ldaf.LT.max( 1, n ) ) THEN
584  info = -8
585  ELSE IF( lsame( fact, 'F' ) .AND. .NOT.
586  $ ( rcequ .OR. lsame( equed, 'N' ) ) ) THEN
587  info = -9
588  ELSE
589  IF ( rcequ ) THEN
590  smin = bignum
591  smax = zero
592  DO 10 j = 1, n
593  smin = min( smin, s( j ) )
594  smax = max( smax, s( j ) )
595  10 CONTINUE
596  IF( smin.LE.zero ) THEN
597  info = -10
598  ELSE IF( n.GT.0 ) THEN
599  scond = max( smin, smlnum ) / min( smax, bignum )
600  ELSE
601  scond = one
602  END IF
603  END IF
604  IF( info.EQ.0 ) THEN
605  IF( ldb.LT.max( 1, n ) ) THEN
606  info = -12
607  ELSE IF( ldx.LT.max( 1, n ) ) THEN
608  info = -14
609  END IF
610  END IF
611  END IF
612 *
613  IF( info.NE.0 ) THEN
614  CALL xerbla( 'DPOSVXX', -info )
615  RETURN
616  END IF
617 *
618  IF( equil ) THEN
619 *
620 * Compute row and column scalings to equilibrate the matrix A.
621 *
622  CALL dpoequb( n, a, lda, s, scond, amax, infequ )
623  IF( infequ.EQ.0 ) THEN
624 *
625 * Equilibrate the matrix.
626 *
627  CALL dlaqsy( uplo, n, a, lda, s, scond, amax, equed )
628  rcequ = lsame( equed, 'Y' )
629  END IF
630  END IF
631 *
632 * Scale the right-hand side.
633 *
634  IF( rcequ ) CALL dlascl2( n, nrhs, s, b, ldb )
635 *
636  IF( nofact .OR. equil ) THEN
637 *
638 * Compute the Cholesky factorization of A.
639 *
640  CALL dlacpy( uplo, n, n, a, lda, af, ldaf )
641  CALL dpotrf( uplo, n, af, ldaf, info )
642 *
643 * Return if INFO is non-zero.
644 *
645  IF( info.NE.0 ) THEN
646 *
647 * Pivot in column INFO is exactly 0
648 * Compute the reciprocal pivot growth factor of the
649 * leading rank-deficient INFO columns of A.
650 *
651  rpvgrw = dla_porpvgrw( uplo, info, a, lda, af, ldaf, work )
652  RETURN
653  ENDIF
654  END IF
655 *
656 * Compute the reciprocal growth factor RPVGRW.
657 *
658  rpvgrw = dla_porpvgrw( uplo, n, a, lda, af, ldaf, work )
659 *
660 * Compute the solution matrix X.
661 *
662  CALL dlacpy( 'Full', n, nrhs, b, ldb, x, ldx )
663  CALL dpotrs( uplo, n, nrhs, af, ldaf, x, ldx, info )
664 *
665 * Use iterative refinement to improve the computed solution and
666 * compute error bounds and backward error estimates for it.
667 *
668  CALL dporfsx( uplo, equed, n, nrhs, a, lda, af, ldaf,
669  $ s, b, ldb, x, ldx, rcond, berr, n_err_bnds, err_bnds_norm,
670  $ err_bnds_comp, nparams, params, work, iwork, info )
671 
672 *
673 * Scale solutions.
674 *
675  IF ( rcequ ) THEN
676  CALL dlascl2 ( n, nrhs, s, x, ldx )
677  END IF
678 *
679  RETURN
680 *
681 * End of DPOSVXX
682 *
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 dlaqsy(UPLO, N, A, LDA, S, SCOND, AMAX, EQUED)
DLAQSY scales a symmetric/Hermitian matrix, using scaling factors computed by spoequ.
Definition: dlaqsy.f:135
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine dlascl2(M, N, D, X, LDX)
DLASCL2 performs diagonal scaling on a vector.
Definition: dlascl2.f:92
subroutine dporfsx(UPLO, EQUED, N, NRHS, A, LDA, AF, LDAF, S, B, LDB, X, LDX, RCOND, BERR, N_ERR_BNDS, ERR_BNDS_NORM, ERR_BNDS_COMP, NPARAMS, PARAMS, WORK, IWORK, INFO)
DPORFSX
Definition: dporfsx.f:396
double precision function dlamch(CMACH)
DLAMCH
Definition: dlamch.f:65
subroutine dpoequb(N, A, LDA, S, SCOND, AMAX, INFO)
DPOEQUB
Definition: dpoequb.f:114
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55
subroutine dpotrf(UPLO, N, A, LDA, INFO)
DPOTRF
Definition: dpotrf.f:109
subroutine dpotrs(UPLO, N, NRHS, A, LDA, B, LDB, INFO)
DPOTRS
Definition: dpotrs.f:112
double precision function dla_porpvgrw(UPLO, NCOLS, A, LDA, AF, LDAF, WORK)
DLA_PORPVGRW computes the reciprocal pivot growth factor norm(A)/norm(U) for a symmetric or Hermitian...
Definition: dla_porpvgrw.f:108

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine dsposv ( character  UPLO,
integer  N,
integer  NRHS,
double precision, dimension( lda, * )  A,
integer  LDA,
double precision, dimension( ldb, * )  B,
integer  LDB,
double precision, dimension( ldx, * )  X,
integer  LDX,
double precision, dimension( n, * )  WORK,
real, dimension( * )  SWORK,
integer  ITER,
integer  INFO 
)

DSPOSV computes the solution to system of linear equations A * X = B for PO matrices

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

Purpose:
 DSPOSV computes the solution to a real system of linear equations
    A * X = B,
 where A is an N-by-N symmetric positive definite matrix and X and B
 are N-by-NRHS matrices.

 DSPOSV first attempts to factorize the matrix in SINGLE PRECISION
 and use this factorization within an iterative refinement procedure
 to produce a solution with DOUBLE PRECISION normwise backward error
 quality (see below). If the approach fails the method switches to a
 DOUBLE PRECISION factorization and solve.

 The iterative refinement is not going to be a winning strategy if
 the ratio SINGLE PRECISION performance over DOUBLE PRECISION
 performance is too small. A reasonable strategy should take the
 number of right-hand sides and the size of the matrix into account.
 This might be done with a call to ILAENV in the future. Up to now, we
 always try iterative refinement.

 The iterative refinement process is stopped if
     ITER > ITERMAX
 or for all the RHS we have:
     RNRM < SQRT(N)*XNRM*ANRM*EPS*BWDMAX
 where
     o ITER is the number of the current iteration in the iterative
       refinement process
     o RNRM is the infinity-norm of the residual
     o XNRM is the infinity-norm of the solution
     o ANRM is the infinity-operator-norm of the matrix A
     o EPS is the machine epsilon returned by DLAMCH('Epsilon')
 The value ITERMAX and BWDMAX are fixed to 30 and 1.0D+00
 respectively.
Parameters
[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 matrix B.  NRHS >= 0.
[in,out]A
          A is DOUBLE PRECISION array,
          dimension (LDA,N)
          On entry, the symmetric matrix A.  If UPLO = 'U', the leading
          N-by-N upper triangular part of A contains the upper
          triangular part of the matrix A, and the strictly lower
          triangular part of A is not referenced.  If UPLO = 'L', the
          leading N-by-N lower triangular part of A contains the lower
          triangular part of the matrix A, and the strictly upper
          triangular part of A is not referenced.
          On exit, if iterative refinement has been successfully used
          (INFO.EQ.0 and ITER.GE.0, see description below), then A is
          unchanged, if double precision factorization has been used
          (INFO.EQ.0 and ITER.LT.0, see description below), then the
          array A contains the factor U or L from the Cholesky
          factorization A = U**T*U or A = L*L**T.
[in]LDA
          LDA is INTEGER
          The leading dimension of the array A.  LDA >= max(1,N).
[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, the N-by-NRHS solution matrix X.
[in]LDX
          LDX is INTEGER
          The leading dimension of the array X.  LDX >= max(1,N).
[out]WORK
          WORK is DOUBLE PRECISION array, dimension (N,NRHS)
          This array is used to hold the residual vectors.
[out]SWORK
          SWORK is REAL array, dimension (N*(N+NRHS))
          This array is used to use the single precision matrix and the
          right-hand sides or solutions in single precision.
[out]ITER
          ITER is INTEGER
          < 0: iterative refinement has failed, double precision
               factorization has been performed
               -1 : the routine fell back to full precision for
                    implementation- or machine-specific reasons
               -2 : narrowing the precision induced an overflow,
                    the routine fell back to full precision
               -3 : failure of SPOTRF
               -31: stop the iterative refinement after the 30th
                    iterations
          > 0: iterative refinement has been sucessfully used.
               Returns the number of iterations
[out]INFO
          INFO is INTEGER
          = 0:  successful exit
          < 0:  if INFO = -i, the i-th argument had an illegal value
          > 0:  if INFO = i, the leading minor of order i of (DOUBLE
                PRECISION) A is not positive definite, so the
                factorization could not be completed, and the solution
                has not been computed.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
November 2011

Definition at line 201 of file dsposv.f.

201 *
202 * -- LAPACK driver routine (version 3.4.0) --
203 * -- LAPACK is a software package provided by Univ. of Tennessee, --
204 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
205 * November 2011
206 *
207 * .. Scalar Arguments ..
208  CHARACTER uplo
209  INTEGER info, iter, lda, ldb, ldx, n, nrhs
210 * ..
211 * .. Array Arguments ..
212  REAL swork( * )
213  DOUBLE PRECISION a( lda, * ), b( ldb, * ), work( n, * ),
214  $ x( ldx, * )
215 * ..
216 *
217 * =====================================================================
218 *
219 * .. Parameters ..
220  LOGICAL doitref
221  parameter( doitref = .true. )
222 *
223  INTEGER itermax
224  parameter( itermax = 30 )
225 *
226  DOUBLE PRECISION bwdmax
227  parameter( bwdmax = 1.0e+00 )
228 *
229  DOUBLE PRECISION negone, one
230  parameter( negone = -1.0d+0, one = 1.0d+0 )
231 *
232 * .. Local Scalars ..
233  INTEGER i, iiter, ptsa, ptsx
234  DOUBLE PRECISION anrm, cte, eps, rnrm, xnrm
235 *
236 * .. External Subroutines ..
237  EXTERNAL daxpy, dsymm, dlacpy, dlat2s, dlag2s, slag2d,
238  $ spotrf, spotrs, xerbla
239 * ..
240 * .. External Functions ..
241  INTEGER idamax
242  DOUBLE PRECISION dlamch, dlansy
243  LOGICAL lsame
244  EXTERNAL idamax, dlamch, dlansy, lsame
245 * ..
246 * .. Intrinsic Functions ..
247  INTRINSIC abs, dble, max, sqrt
248 * ..
249 * .. Executable Statements ..
250 *
251  info = 0
252  iter = 0
253 *
254 * Test the input parameters.
255 *
256  IF( .NOT.lsame( uplo, 'U' ) .AND. .NOT.lsame( uplo, 'L' ) ) THEN
257  info = -1
258  ELSE IF( n.LT.0 ) THEN
259  info = -2
260  ELSE IF( nrhs.LT.0 ) THEN
261  info = -3
262  ELSE IF( lda.LT.max( 1, n ) ) THEN
263  info = -5
264  ELSE IF( ldb.LT.max( 1, n ) ) THEN
265  info = -7
266  ELSE IF( ldx.LT.max( 1, n ) ) THEN
267  info = -9
268  END IF
269  IF( info.NE.0 ) THEN
270  CALL xerbla( 'DSPOSV', -info )
271  RETURN
272  END IF
273 *
274 * Quick return if (N.EQ.0).
275 *
276  IF( n.EQ.0 )
277  $ RETURN
278 *
279 * Skip single precision iterative refinement if a priori slower
280 * than double precision factorization.
281 *
282  IF( .NOT.doitref ) THEN
283  iter = -1
284  GO TO 40
285  END IF
286 *
287 * Compute some constants.
288 *
289  anrm = dlansy( 'I', uplo, n, a, lda, work )
290  eps = dlamch( 'Epsilon' )
291  cte = anrm*eps*sqrt( dble( n ) )*bwdmax
292 *
293 * Set the indices PTSA, PTSX for referencing SA and SX in SWORK.
294 *
295  ptsa = 1
296  ptsx = ptsa + n*n
297 *
298 * Convert B from double precision to single precision and store the
299 * result in SX.
300 *
301  CALL dlag2s( n, nrhs, b, ldb, swork( ptsx ), n, info )
302 *
303  IF( info.NE.0 ) THEN
304  iter = -2
305  GO TO 40
306  END IF
307 *
308 * Convert A from double precision to single precision and store the
309 * result in SA.
310 *
311  CALL dlat2s( uplo, n, a, lda, swork( ptsa ), n, info )
312 *
313  IF( info.NE.0 ) THEN
314  iter = -2
315  GO TO 40
316  END IF
317 *
318 * Compute the Cholesky factorization of SA.
319 *
320  CALL spotrf( uplo, n, swork( ptsa ), n, info )
321 *
322  IF( info.NE.0 ) THEN
323  iter = -3
324  GO TO 40
325  END IF
326 *
327 * Solve the system SA*SX = SB.
328 *
329  CALL spotrs( uplo, n, nrhs, swork( ptsa ), n, swork( ptsx ), n,
330  $ info )
331 *
332 * Convert SX back to double precision
333 *
334  CALL slag2d( n, nrhs, swork( ptsx ), n, x, ldx, info )
335 *
336 * Compute R = B - AX (R is WORK).
337 *
338  CALL dlacpy( 'All', n, nrhs, b, ldb, work, n )
339 *
340  CALL dsymm( 'Left', uplo, n, nrhs, negone, a, lda, x, ldx, one,
341  $ work, n )
342 *
343 * Check whether the NRHS normwise backward errors satisfy the
344 * stopping criterion. If yes, set ITER=0 and return.
345 *
346  DO i = 1, nrhs
347  xnrm = abs( x( idamax( n, x( 1, i ), 1 ), i ) )
348  rnrm = abs( work( idamax( n, work( 1, i ), 1 ), i ) )
349  IF( rnrm.GT.xnrm*cte )
350  $ GO TO 10
351  END DO
352 *
353 * If we are here, the NRHS normwise backward errors satisfy the
354 * stopping criterion. We are good to exit.
355 *
356  iter = 0
357  RETURN
358 *
359  10 CONTINUE
360 *
361  DO 30 iiter = 1, itermax
362 *
363 * Convert R (in WORK) from double precision to single precision
364 * and store the result in SX.
365 *
366  CALL dlag2s( n, nrhs, work, n, swork( ptsx ), n, info )
367 *
368  IF( info.NE.0 ) THEN
369  iter = -2
370  GO TO 40
371  END IF
372 *
373 * Solve the system SA*SX = SR.
374 *
375  CALL spotrs( uplo, n, nrhs, swork( ptsa ), n, swork( ptsx ), n,
376  $ info )
377 *
378 * Convert SX back to double precision and update the current
379 * iterate.
380 *
381  CALL slag2d( n, nrhs, swork( ptsx ), n, work, n, info )
382 *
383  DO i = 1, nrhs
384  CALL daxpy( n, one, work( 1, i ), 1, x( 1, i ), 1 )
385  END DO
386 *
387 * Compute R = B - AX (R is WORK).
388 *
389  CALL dlacpy( 'All', n, nrhs, b, ldb, work, n )
390 *
391  CALL dsymm( 'L', uplo, n, nrhs, negone, a, lda, x, ldx, one,
392  $ work, n )
393 *
394 * Check whether the NRHS normwise backward errors satisfy the
395 * stopping criterion. If yes, set ITER=IITER>0 and return.
396 *
397  DO i = 1, nrhs
398  xnrm = abs( x( idamax( n, x( 1, i ), 1 ), i ) )
399  rnrm = abs( work( idamax( n, work( 1, i ), 1 ), i ) )
400  IF( rnrm.GT.xnrm*cte )
401  $ GO TO 20
402  END DO
403 *
404 * If we are here, the NRHS normwise backward errors satisfy the
405 * stopping criterion, we are good to exit.
406 *
407  iter = iiter
408 *
409  RETURN
410 *
411  20 CONTINUE
412 *
413  30 CONTINUE
414 *
415 * If we are at this place of the code, this is because we have
416 * performed ITER=ITERMAX iterations and never satisified the
417 * stopping criterion, set up the ITER flag accordingly and follow
418 * up on double precision routine.
419 *
420  iter = -itermax - 1
421 *
422  40 CONTINUE
423 *
424 * Single-precision iterative refinement failed to converge to a
425 * satisfactory solution, so we resort to double precision.
426 *
427  CALL dpotrf( uplo, n, a, lda, info )
428 *
429  IF( info.NE.0 )
430  $ RETURN
431 *
432  CALL dlacpy( 'All', n, nrhs, b, ldb, x, ldx )
433  CALL dpotrs( uplo, n, nrhs, a, lda, x, ldx, info )
434 *
435  RETURN
436 *
437 * End of DSPOSV.
438 *
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 spotrs(UPLO, N, NRHS, A, LDA, B, LDB, INFO)
SPOTRS
Definition: spotrs.f:112
subroutine dsymm(SIDE, UPLO, M, N, ALPHA, A, LDA, B, LDB, BETA, C, LDC)
DSYMM
Definition: dsymm.f:191
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine dlag2s(M, N, A, LDA, SA, LDSA, INFO)
DLAG2S converts a double precision matrix to a single precision matrix.
Definition: dlag2s.f:110
subroutine spotrf(UPLO, N, A, LDA, INFO)
SPOTRF
Definition: spotrf.f:109
double precision function dlamch(CMACH)
DLAMCH
Definition: dlamch.f:65
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55
subroutine dpotrf(UPLO, N, A, LDA, INFO)
DPOTRF
Definition: dpotrf.f:109
integer function idamax(N, DX, INCX)
IDAMAX
Definition: idamax.f:53
subroutine dpotrs(UPLO, N, NRHS, A, LDA, B, LDB, INFO)
DPOTRS
Definition: dpotrs.f:112
subroutine daxpy(N, DA, DX, INCX, DY, INCY)
DAXPY
Definition: daxpy.f:54
subroutine slag2d(M, N, SA, LDSA, A, LDA, INFO)
SLAG2D converts a single precision matrix to a double precision matrix.
Definition: slag2d.f:106
double precision function dlansy(NORM, UPLO, N, A, LDA, WORK)
DLANSY 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 matrix.
Definition: dlansy.f:124
subroutine dlat2s(UPLO, N, A, LDA, SA, LDSA, INFO)
DLAT2S converts a double-precision triangular matrix to a single-precision triangular matrix...
Definition: dlat2s.f:113

Here is the call graph for this function:

Here is the caller graph for this function: