.TH "doublePOcomputational" 3 "Sun Nov 27 2022" "Version 3.11.0" "LAPACK" \" -*- nroff -*-
.ad l
.nh
.SH NAME
doublePOcomputational \- double
.SH SYNOPSIS
.br
.PP
.SS "Functions"

.in +1c
.ti -1c
.RI "double precision function \fBdla_porcond\fP (UPLO, N, A, LDA, AF, LDAF, CMODE, C, INFO, WORK, IWORK)"
.br
.RI "\fBDLA_PORCOND\fP estimates the Skeel condition number for a symmetric positive-definite matrix\&. "
.ti -1c
.RI "subroutine \fBdla_porfsx_extended\fP (PREC_TYPE, UPLO, N, NRHS, A, LDA, AF, LDAF, COLEQU, C, B, LDB, Y, LDY, BERR_OUT, N_NORMS, ERR_BNDS_NORM, ERR_BNDS_COMP, RES, AYB, DY, Y_TAIL, RCOND, ITHRESH, RTHRESH, DZ_UB, IGNORE_CWISE, INFO)"
.br
.RI "\fBDLA_PORFSX_EXTENDED\fP improves the computed solution to a system of linear equations for symmetric or Hermitian positive-definite matrices by performing extra-precise iterative refinement and provides error bounds and backward error estimates for the solution\&. "
.ti -1c
.RI "double precision function \fBdla_porpvgrw\fP (UPLO, NCOLS, A, LDA, AF, LDAF, WORK)"
.br
.RI "\fBDLA_PORPVGRW\fP computes the reciprocal pivot growth factor norm(A)/norm(U) for a symmetric or Hermitian positive-definite matrix\&. "
.ti -1c
.RI "subroutine \fBdpocon\fP (UPLO, N, A, LDA, ANORM, RCOND, WORK, IWORK, INFO)"
.br
.RI "\fBDPOCON\fP "
.ti -1c
.RI "subroutine \fBdpoequ\fP (N, A, LDA, S, SCOND, AMAX, INFO)"
.br
.RI "\fBDPOEQU\fP "
.ti -1c
.RI "subroutine \fBdpoequb\fP (N, A, LDA, S, SCOND, AMAX, INFO)"
.br
.RI "\fBDPOEQUB\fP "
.ti -1c
.RI "subroutine \fBdporfs\fP (UPLO, N, NRHS, A, LDA, AF, LDAF, B, LDB, X, LDX, FERR, BERR, WORK, IWORK, INFO)"
.br
.RI "\fBDPORFS\fP "
.ti -1c
.RI "subroutine \fBdporfsx\fP (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)"
.br
.RI "\fBDPORFSX\fP "
.ti -1c
.RI "subroutine \fBdpotf2\fP (UPLO, N, A, LDA, INFO)"
.br
.RI "\fBDPOTF2\fP computes the Cholesky factorization of a symmetric/Hermitian positive definite matrix (unblocked algorithm)\&. "
.ti -1c
.RI "subroutine \fBdpotrf\fP (UPLO, N, A, LDA, INFO)"
.br
.RI "\fBDPOTRF\fP "
.ti -1c
.RI "recursive subroutine \fBdpotrf2\fP (UPLO, N, A, LDA, INFO)"
.br
.RI "\fBDPOTRF2\fP "
.ti -1c
.RI "subroutine \fBdpotri\fP (UPLO, N, A, LDA, INFO)"
.br
.RI "\fBDPOTRI\fP "
.ti -1c
.RI "subroutine \fBdpotrs\fP (UPLO, N, NRHS, A, LDA, B, LDB, INFO)"
.br
.RI "\fBDPOTRS\fP "
.in -1c
.SH "Detailed Description"
.PP 
This is the group of double computational functions for PO matrices 
.SH "Function Documentation"
.PP 
.SS "double precision function dla_porcond (character UPLO, integer N, double precision, dimension( lda, * ) A, integer LDA, double precision, dimension( ldaf, * ) AF, integer LDAF, integer CMODE, double precision, dimension( * ) C, integer INFO, double precision, dimension( * ) WORK, integer, dimension( * ) IWORK)"

.PP
\fBDLA_PORCOND\fP estimates the Skeel condition number for a symmetric positive-definite matrix\&.  
.PP
\fBPurpose:\fP
.RS 4

.PP
.nf
    DLA_PORCOND Estimates the Skeel condition number of  op(A) * op2(C)
    where op2 is determined by CMODE as follows
    CMODE =  1    op2(C) = C
    CMODE =  0    op2(C) = I
    CMODE = -1    op2(C) = inv(C)
    The Skeel condition number  cond(A) = norminf( |inv(A)||A| )
    is computed by computing scaling factors R such that
    diag(R)*A*op2(C) is row equilibrated and computing the standard
    infinity-norm condition number\&.
.fi
.PP
 
.RE
.PP
\fBParameters\fP
.RS 4
\fIUPLO\fP 
.PP
.nf
          UPLO is CHARACTER*1
       = 'U':  Upper triangle of A is stored;
       = 'L':  Lower triangle of A is stored\&.
.fi
.PP
.br
\fIN\fP 
.PP
.nf
          N is INTEGER
     The number of linear equations, i\&.e\&., the order of the
     matrix A\&.  N >= 0\&.
.fi
.PP
.br
\fIA\fP 
.PP
.nf
          A is DOUBLE PRECISION array, dimension (LDA,N)
     On entry, the N-by-N matrix A\&.
.fi
.PP
.br
\fILDA\fP 
.PP
.nf
          LDA is INTEGER
     The leading dimension of the array A\&.  LDA >= max(1,N)\&.
.fi
.PP
.br
\fIAF\fP 
.PP
.nf
          AF is DOUBLE PRECISION array, dimension (LDAF,N)
     The triangular factor U or L from the Cholesky factorization
     A = U**T*U or A = L*L**T, as computed by DPOTRF\&.
.fi
.PP
.br
\fILDAF\fP 
.PP
.nf
          LDAF is INTEGER
     The leading dimension of the array AF\&.  LDAF >= max(1,N)\&.
.fi
.PP
.br
\fICMODE\fP 
.PP
.nf
          CMODE is INTEGER
     Determines op2(C) in the formula op(A) * op2(C) as follows:
     CMODE =  1    op2(C) = C
     CMODE =  0    op2(C) = I
     CMODE = -1    op2(C) = inv(C)
.fi
.PP
.br
\fIC\fP 
.PP
.nf
          C is DOUBLE PRECISION array, dimension (N)
     The vector C in the formula op(A) * op2(C)\&.
.fi
.PP
.br
\fIINFO\fP 
.PP
.nf
          INFO is INTEGER
       = 0:  Successful exit\&.
     i > 0:  The ith argument is invalid\&.
.fi
.PP
.br
\fIWORK\fP 
.PP
.nf
          WORK is DOUBLE PRECISION array, dimension (3*N)\&.
     Workspace\&.
.fi
.PP
.br
\fIIWORK\fP 
.PP
.nf
          IWORK is INTEGER array, dimension (N)\&.
     Workspace\&.
.fi
.PP
 
.RE
.PP
\fBAuthor\fP
.RS 4
Univ\&. of Tennessee 
.PP
Univ\&. of California Berkeley 
.PP
Univ\&. of Colorado Denver 
.PP
NAG Ltd\&. 
.RE
.PP

.SS "subroutine dla_porfsx_extended (integer PREC_TYPE, character UPLO, integer N, integer NRHS, double precision, dimension( lda, * ) A, integer LDA, double precision, dimension( ldaf, * ) AF, integer LDAF, logical COLEQU, double precision, dimension( * ) C, double precision, dimension( ldb, * ) B, integer LDB, double precision, dimension( ldy, * ) Y, integer LDY, double precision, dimension( * ) BERR_OUT, integer N_NORMS, double precision, dimension( nrhs, * ) ERR_BNDS_NORM, double precision, dimension( nrhs, * ) ERR_BNDS_COMP, double precision, dimension( * ) RES, double precision, dimension(*) AYB, double precision, dimension( * ) DY, double precision, dimension( * ) Y_TAIL, double precision RCOND, integer ITHRESH, double precision RTHRESH, double precision DZ_UB, logical IGNORE_CWISE, integer INFO)"

.PP
\fBDLA_PORFSX_EXTENDED\fP improves the computed solution to a system of linear equations for symmetric or Hermitian positive-definite matrices by performing extra-precise iterative refinement and provides error bounds and backward error estimates for the solution\&.  
.PP
\fBPurpose:\fP
.RS 4

.PP
.nf
 DLA_PORFSX_EXTENDED improves the computed solution to a system of
 linear equations by performing extra-precise iterative refinement
 and provides error bounds and backward error estimates for the solution\&.
 This subroutine is called by DPORFSX to perform iterative refinement\&.
 In addition to normwise error bound, the code provides maximum
 componentwise error bound if possible\&. See comments for ERR_BNDS_NORM
 and ERR_BNDS_COMP for details of the error bounds\&. Note that this
 subroutine is only responsible for setting the second fields of
 ERR_BNDS_NORM and ERR_BNDS_COMP\&.
.fi
.PP
 
.RE
.PP
\fBParameters\fP
.RS 4
\fIPREC_TYPE\fP 
.PP
.nf
          PREC_TYPE is INTEGER
     Specifies the intermediate precision to be used in refinement\&.
     The value is defined by ILAPREC(P) where P is a CHARACTER and P
          = 'S':  Single
          = 'D':  Double
          = 'I':  Indigenous
          = 'X' or 'E':  Extra
.fi
.PP
.br
\fIUPLO\fP 
.PP
.nf
          UPLO is CHARACTER*1
       = 'U':  Upper triangle of A is stored;
       = 'L':  Lower triangle of A is stored\&.
.fi
.PP
.br
\fIN\fP 
.PP
.nf
          N is INTEGER
     The number of linear equations, i\&.e\&., the order of the
     matrix A\&.  N >= 0\&.
.fi
.PP
.br
\fINRHS\fP 
.PP
.nf
          NRHS is INTEGER
     The number of right-hand-sides, i\&.e\&., the number of columns of the
     matrix B\&.
.fi
.PP
.br
\fIA\fP 
.PP
.nf
          A is DOUBLE PRECISION array, dimension (LDA,N)
     On entry, the N-by-N matrix A\&.
.fi
.PP
.br
\fILDA\fP 
.PP
.nf
          LDA is INTEGER
     The leading dimension of the array A\&.  LDA >= max(1,N)\&.
.fi
.PP
.br
\fIAF\fP 
.PP
.nf
          AF is DOUBLE PRECISION array, dimension (LDAF,N)
     The triangular factor U or L from the Cholesky factorization
     A = U**T*U or A = L*L**T, as computed by DPOTRF\&.
.fi
.PP
.br
\fILDAF\fP 
.PP
.nf
          LDAF is INTEGER
     The leading dimension of the array AF\&.  LDAF >= max(1,N)\&.
.fi
.PP
.br
\fICOLEQU\fP 
.PP
.nf
          COLEQU is LOGICAL
     If \&.TRUE\&. then column equilibration was done to A before calling
     this routine\&. This is needed to compute the solution and error
     bounds correctly\&.
.fi
.PP
.br
\fIC\fP 
.PP
.nf
          C is DOUBLE PRECISION array, dimension (N)
     The column scale factors for A\&. If COLEQU = \&.FALSE\&., C
     is not accessed\&. If C is input, each element of C 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\&.
.fi
.PP
.br
\fIB\fP 
.PP
.nf
          B is DOUBLE PRECISION array, dimension (LDB,NRHS)
     The right-hand-side matrix B\&.
.fi
.PP
.br
\fILDB\fP 
.PP
.nf
          LDB is INTEGER
     The leading dimension of the array B\&.  LDB >= max(1,N)\&.
.fi
.PP
.br
\fIY\fP 
.PP
.nf
          Y is DOUBLE PRECISION array, dimension (LDY,NRHS)
     On entry, the solution matrix X, as computed by DPOTRS\&.
     On exit, the improved solution matrix Y\&.
.fi
.PP
.br
\fILDY\fP 
.PP
.nf
          LDY is INTEGER
     The leading dimension of the array Y\&.  LDY >= max(1,N)\&.
.fi
.PP
.br
\fIBERR_OUT\fP 
.PP
.nf
          BERR_OUT is DOUBLE PRECISION array, dimension (NRHS)
     On exit, BERR_OUT(j) contains the componentwise relative backward
     error for right-hand-side j from the formula
         max(i) ( abs(RES(i)) / ( abs(op(A_s))*abs(Y) + abs(B_s) )(i) )
     where abs(Z) is the componentwise absolute value of the matrix
     or vector Z\&. This is computed by DLA_LIN_BERR\&.
.fi
.PP
.br
\fIN_NORMS\fP 
.PP
.nf
          N_NORMS is INTEGER
     Determines which error bounds to return (see ERR_BNDS_NORM
     and ERR_BNDS_COMP)\&.
     If N_NORMS >= 1 return normwise error bounds\&.
     If N_NORMS >= 2 return componentwise error bounds\&.
.fi
.PP
.br
\fIERR_BNDS_NORM\fP 
.PP
.nf
          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) * slamch('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) * slamch('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) * slamch('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\&.

     This subroutine is only responsible for setting the second field
     above\&.
     See Lapack Working Note 165 for further details and extra
     cautions\&.
.fi
.PP
.br
\fIERR_BNDS_COMP\fP 
.PP
.nf
          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 < 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) * slamch('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) * slamch('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) * slamch('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\&.

     This subroutine is only responsible for setting the second field
     above\&.
     See Lapack Working Note 165 for further details and extra
     cautions\&.
.fi
.PP
.br
\fIRES\fP 
.PP
.nf
          RES is DOUBLE PRECISION array, dimension (N)
     Workspace to hold the intermediate residual\&.
.fi
.PP
.br
\fIAYB\fP 
.PP
.nf
          AYB is DOUBLE PRECISION array, dimension (N)
     Workspace\&. This can be the same workspace passed for Y_TAIL\&.
.fi
.PP
.br
\fIDY\fP 
.PP
.nf
          DY is DOUBLE PRECISION array, dimension (N)
     Workspace to hold the intermediate solution\&.
.fi
.PP
.br
\fIY_TAIL\fP 
.PP
.nf
          Y_TAIL is DOUBLE PRECISION array, dimension (N)
     Workspace to hold the trailing bits of the intermediate solution\&.
.fi
.PP
.br
\fIRCOND\fP 
.PP
.nf
          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\&.
.fi
.PP
.br
\fIITHRESH\fP 
.PP
.nf
          ITHRESH is INTEGER
     The maximum number of residual computations allowed for
     refinement\&. The default is 10\&. For '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\&.
.fi
.PP
.br
\fIRTHRESH\fP 
.PP
.nf
          RTHRESH is DOUBLE PRECISION
     Determines when to stop refinement if the error estimate stops
     decreasing\&. Refinement will stop when the next solution no longer
     satisfies norm(dx_{i+1}) < RTHRESH * norm(dx_i) where norm(Z) is
     the infinity norm of Z\&. RTHRESH satisfies 0 < RTHRESH <= 1\&. The
     default value is 0\&.5\&. For 'aggressive' set to 0\&.9 to permit
     convergence on extremely ill-conditioned matrices\&. See LAWN 165
     for more details\&.
.fi
.PP
.br
\fIDZ_UB\fP 
.PP
.nf
          DZ_UB is DOUBLE PRECISION
     Determines when to start considering componentwise convergence\&.
     Componentwise convergence is only considered after each component
     of the solution Y is stable, which we define as the relative
     change in each component being less than DZ_UB\&. The default value
     is 0\&.25, requiring the first bit to be stable\&. See LAWN 165 for
     more details\&.
.fi
.PP
.br
\fIIGNORE_CWISE\fP 
.PP
.nf
          IGNORE_CWISE is LOGICAL
     If \&.TRUE\&. then ignore componentwise convergence\&. Default value
     is \&.FALSE\&.\&.
.fi
.PP
.br
\fIINFO\fP 
.PP
.nf
          INFO is INTEGER
       = 0:  Successful exit\&.
       < 0:  if INFO = -i, the ith argument to DPOTRS had an illegal
             value
.fi
.PP
 
.RE
.PP
\fBAuthor\fP
.RS 4
Univ\&. of Tennessee 
.PP
Univ\&. of California Berkeley 
.PP
Univ\&. of Colorado Denver 
.PP
NAG Ltd\&. 
.RE
.PP

.SS "double precision function dla_porpvgrw (character*1 UPLO, integer NCOLS, double precision, dimension( lda, * ) A, integer LDA, double precision, dimension( ldaf, * ) AF, integer LDAF, double precision, dimension( * ) WORK)"

.PP
\fBDLA_PORPVGRW\fP computes the reciprocal pivot growth factor norm(A)/norm(U) for a symmetric or Hermitian positive-definite matrix\&.  
.PP
\fBPurpose:\fP
.RS 4

.PP
.nf
 DLA_PORPVGRW computes the reciprocal pivot growth factor
 norm(A)/norm(U)\&. The 'max absolute element' norm is used\&. If this is
 much less than 1, 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\&.
.fi
.PP
 
.RE
.PP
\fBParameters\fP
.RS 4
\fIUPLO\fP 
.PP
.nf
          UPLO is CHARACTER*1
       = 'U':  Upper triangle of A is stored;
       = 'L':  Lower triangle of A is stored\&.
.fi
.PP
.br
\fINCOLS\fP 
.PP
.nf
          NCOLS is INTEGER
     The number of columns of the matrix A\&. NCOLS >= 0\&.
.fi
.PP
.br
\fIA\fP 
.PP
.nf
          A is DOUBLE PRECISION array, dimension (LDA,N)
     On entry, the N-by-N matrix A\&.
.fi
.PP
.br
\fILDA\fP 
.PP
.nf
          LDA is INTEGER
     The leading dimension of the array A\&.  LDA >= max(1,N)\&.
.fi
.PP
.br
\fIAF\fP 
.PP
.nf
          AF is DOUBLE PRECISION array, dimension (LDAF,N)
     The triangular factor U or L from the Cholesky factorization
     A = U**T*U or A = L*L**T, as computed by DPOTRF\&.
.fi
.PP
.br
\fILDAF\fP 
.PP
.nf
          LDAF is INTEGER
     The leading dimension of the array AF\&.  LDAF >= max(1,N)\&.
.fi
.PP
.br
\fIWORK\fP 
.PP
.nf
          WORK is DOUBLE PRECISION array, dimension (2*N)
.fi
.PP
 
.RE
.PP
\fBAuthor\fP
.RS 4
Univ\&. of Tennessee 
.PP
Univ\&. of California Berkeley 
.PP
Univ\&. of Colorado Denver 
.PP
NAG Ltd\&. 
.RE
.PP

.SS "subroutine dpocon (character UPLO, integer N, double precision, dimension( lda, * ) A, integer LDA, double precision ANORM, double precision RCOND, double precision, dimension( * ) WORK, integer, dimension( * ) IWORK, integer INFO)"

.PP
\fBDPOCON\fP  
.PP
\fBPurpose:\fP
.RS 4

.PP
.nf
 DPOCON estimates the reciprocal of the condition number (in the
 1-norm) of a real symmetric positive definite matrix using the
 Cholesky factorization A = U**T*U or A = L*L**T computed by DPOTRF\&.

 An estimate is obtained for norm(inv(A)), and the reciprocal of the
 condition number is computed as RCOND = 1 / (ANORM * norm(inv(A)))\&.
.fi
.PP
 
.RE
.PP
\fBParameters\fP
.RS 4
\fIUPLO\fP 
.PP
.nf
          UPLO is CHARACTER*1
          = 'U':  Upper triangle of A is stored;
          = 'L':  Lower triangle of A is stored\&.
.fi
.PP
.br
\fIN\fP 
.PP
.nf
          N is INTEGER
          The order of the matrix A\&.  N >= 0\&.
.fi
.PP
.br
\fIA\fP 
.PP
.nf
          A is DOUBLE PRECISION array, dimension (LDA,N)
          The triangular factor U or L from the Cholesky factorization
          A = U**T*U or A = L*L**T, as computed by DPOTRF\&.
.fi
.PP
.br
\fILDA\fP 
.PP
.nf
          LDA is INTEGER
          The leading dimension of the array A\&.  LDA >= max(1,N)\&.
.fi
.PP
.br
\fIANORM\fP 
.PP
.nf
          ANORM is DOUBLE PRECISION
          The 1-norm (or infinity-norm) of the symmetric matrix A\&.
.fi
.PP
.br
\fIRCOND\fP 
.PP
.nf
          RCOND is DOUBLE PRECISION
          The reciprocal of the condition number of the matrix A,
          computed as RCOND = 1/(ANORM * AINVNM), where AINVNM is an
          estimate of the 1-norm of inv(A) computed in this routine\&.
.fi
.PP
.br
\fIWORK\fP 
.PP
.nf
          WORK is DOUBLE PRECISION array, dimension (3*N)
.fi
.PP
.br
\fIIWORK\fP 
.PP
.nf
          IWORK is INTEGER array, dimension (N)
.fi
.PP
.br
\fIINFO\fP 
.PP
.nf
          INFO is INTEGER
          = 0:  successful exit
          < 0:  if INFO = -i, the i-th argument had an illegal value
.fi
.PP
 
.RE
.PP
\fBAuthor\fP
.RS 4
Univ\&. of Tennessee 
.PP
Univ\&. of California Berkeley 
.PP
Univ\&. of Colorado Denver 
.PP
NAG Ltd\&. 
.RE
.PP

.SS "subroutine dpoequ (integer N, double precision, dimension( lda, * ) A, integer LDA, double precision, dimension( * ) S, double precision SCOND, double precision AMAX, integer INFO)"

.PP
\fBDPOEQU\fP  
.PP
\fBPurpose:\fP
.RS 4

.PP
.nf
 DPOEQU computes row and column scalings intended to equilibrate a
 symmetric positive definite matrix A and reduce its condition number
 (with respect to the two-norm)\&.  S contains the scale factors,
 S(i) = 1/sqrt(A(i,i)), chosen so that the scaled matrix B with
 elements B(i,j) = S(i)*A(i,j)*S(j) has ones on the diagonal\&.  This
 choice of S puts the condition number of B within a factor N of the
 smallest possible condition number over all possible diagonal
 scalings\&.
.fi
.PP
 
.RE
.PP
\fBParameters\fP
.RS 4
\fIN\fP 
.PP
.nf
          N is INTEGER
          The order of the matrix A\&.  N >= 0\&.
.fi
.PP
.br
\fIA\fP 
.PP
.nf
          A is DOUBLE PRECISION array, dimension (LDA,N)
          The N-by-N symmetric positive definite matrix whose scaling
          factors are to be computed\&.  Only the diagonal elements of A
          are referenced\&.
.fi
.PP
.br
\fILDA\fP 
.PP
.nf
          LDA is INTEGER
          The leading dimension of the array A\&.  LDA >= max(1,N)\&.
.fi
.PP
.br
\fIS\fP 
.PP
.nf
          S is DOUBLE PRECISION array, dimension (N)
          If INFO = 0, S contains the scale factors for A\&.
.fi
.PP
.br
\fISCOND\fP 
.PP
.nf
          SCOND is DOUBLE PRECISION
          If INFO = 0, S contains the ratio of the smallest S(i) to
          the largest S(i)\&.  If SCOND >= 0\&.1 and AMAX is neither too
          large nor too small, it is not worth scaling by S\&.
.fi
.PP
.br
\fIAMAX\fP 
.PP
.nf
          AMAX is DOUBLE PRECISION
          Absolute value of largest matrix element\&.  If AMAX is very
          close to overflow or very close to underflow, the matrix
          should be scaled\&.
.fi
.PP
.br
\fIINFO\fP 
.PP
.nf
          INFO is INTEGER
          = 0:  successful exit
          < 0:  if INFO = -i, the i-th argument had an illegal value
          > 0:  if INFO = i, the i-th diagonal element is nonpositive\&.
.fi
.PP
 
.RE
.PP
\fBAuthor\fP
.RS 4
Univ\&. of Tennessee 
.PP
Univ\&. of California Berkeley 
.PP
Univ\&. of Colorado Denver 
.PP
NAG Ltd\&. 
.RE
.PP

.SS "subroutine dpoequb (integer N, double precision, dimension( lda, * ) A, integer LDA, double precision, dimension( * ) S, double precision SCOND, double precision AMAX, integer INFO)"

.PP
\fBDPOEQUB\fP  
.PP
\fBPurpose:\fP
.RS 4

.PP
.nf
 DPOEQUB computes row and column scalings intended to equilibrate a
 symmetric positive definite matrix A and reduce its condition number
 (with respect to the two-norm)\&.  S contains the scale factors,
 S(i) = 1/sqrt(A(i,i)), chosen so that the scaled matrix B with
 elements B(i,j) = S(i)*A(i,j)*S(j) has ones on the diagonal\&.  This
 choice of S puts the condition number of B within a factor N of the
 smallest possible condition number over all possible diagonal
 scalings\&.

 This routine differs from DPOEQU by restricting the scaling factors
 to a power of the radix\&.  Barring over- and underflow, scaling by
 these factors introduces no additional rounding errors\&.  However, the
 scaled diagonal entries are no longer approximately 1 but lie
 between sqrt(radix) and 1/sqrt(radix)\&.
.fi
.PP
 
.RE
.PP
\fBParameters\fP
.RS 4
\fIN\fP 
.PP
.nf
          N is INTEGER
          The order of the matrix A\&.  N >= 0\&.
.fi
.PP
.br
\fIA\fP 
.PP
.nf
          A is DOUBLE PRECISION array, dimension (LDA,N)
          The N-by-N symmetric positive definite matrix whose scaling
          factors are to be computed\&.  Only the diagonal elements of A
          are referenced\&.
.fi
.PP
.br
\fILDA\fP 
.PP
.nf
          LDA is INTEGER
          The leading dimension of the array A\&.  LDA >= max(1,N)\&.
.fi
.PP
.br
\fIS\fP 
.PP
.nf
          S is DOUBLE PRECISION array, dimension (N)
          If INFO = 0, S contains the scale factors for A\&.
.fi
.PP
.br
\fISCOND\fP 
.PP
.nf
          SCOND is DOUBLE PRECISION
          If INFO = 0, S contains the ratio of the smallest S(i) to
          the largest S(i)\&.  If SCOND >= 0\&.1 and AMAX is neither too
          large nor too small, it is not worth scaling by S\&.
.fi
.PP
.br
\fIAMAX\fP 
.PP
.nf
          AMAX is DOUBLE PRECISION
          Absolute value of largest matrix element\&.  If AMAX is very
          close to overflow or very close to underflow, the matrix
          should be scaled\&.
.fi
.PP
.br
\fIINFO\fP 
.PP
.nf
          INFO is INTEGER
          = 0:  successful exit
          < 0:  if INFO = -i, the i-th argument had an illegal value
          > 0:  if INFO = i, the i-th diagonal element is nonpositive\&.
.fi
.PP
 
.RE
.PP
\fBAuthor\fP
.RS 4
Univ\&. of Tennessee 
.PP
Univ\&. of California Berkeley 
.PP
Univ\&. of Colorado Denver 
.PP
NAG Ltd\&. 
.RE
.PP

.SS "subroutine dporfs (character UPLO, integer N, integer NRHS, double precision, dimension( lda, * ) A, integer LDA, double precision, dimension( ldaf, * ) AF, integer LDAF, double precision, dimension( ldb, * ) B, integer LDB, double precision, dimension( ldx, * ) X, integer LDX, double precision, dimension( * ) FERR, double precision, dimension( * ) BERR, double precision, dimension( * ) WORK, integer, dimension( * ) IWORK, integer INFO)"

.PP
\fBDPORFS\fP  
.PP
\fBPurpose:\fP
.RS 4

.PP
.nf
 DPORFS improves the computed solution to a system of linear
 equations when the coefficient matrix is symmetric positive definite,
 and provides error bounds and backward error estimates for the
 solution\&.
.fi
.PP
 
.RE
.PP
\fBParameters\fP
.RS 4
\fIUPLO\fP 
.PP
.nf
          UPLO is CHARACTER*1
          = 'U':  Upper triangle of A is stored;
          = 'L':  Lower triangle of A is stored\&.
.fi
.PP
.br
\fIN\fP 
.PP
.nf
          N is INTEGER
          The order of the matrix A\&.  N >= 0\&.
.fi
.PP
.br
\fINRHS\fP 
.PP
.nf
          NRHS is INTEGER
          The number of right hand sides, i\&.e\&., the number of columns
          of the matrices B and X\&.  NRHS >= 0\&.
.fi
.PP
.br
\fIA\fP 
.PP
.nf
          A is DOUBLE PRECISION array, dimension (LDA,N)
          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\&.
.fi
.PP
.br
\fILDA\fP 
.PP
.nf
          LDA is INTEGER
          The leading dimension of the array A\&.  LDA >= max(1,N)\&.
.fi
.PP
.br
\fIAF\fP 
.PP
.nf
          AF is DOUBLE PRECISION array, dimension (LDAF,N)
          The triangular factor U or L from the Cholesky factorization
          A = U**T*U or A = L*L**T, as computed by DPOTRF\&.
.fi
.PP
.br
\fILDAF\fP 
.PP
.nf
          LDAF is INTEGER
          The leading dimension of the array AF\&.  LDAF >= max(1,N)\&.
.fi
.PP
.br
\fIB\fP 
.PP
.nf
          B is DOUBLE PRECISION array, dimension (LDB,NRHS)
          The right hand side matrix B\&.
.fi
.PP
.br
\fILDB\fP 
.PP
.nf
          LDB is INTEGER
          The leading dimension of the array B\&.  LDB >= max(1,N)\&.
.fi
.PP
.br
\fIX\fP 
.PP
.nf
          X is DOUBLE PRECISION array, dimension (LDX,NRHS)
          On entry, the solution matrix X, as computed by DPOTRS\&.
          On exit, the improved solution matrix X\&.
.fi
.PP
.br
\fILDX\fP 
.PP
.nf
          LDX is INTEGER
          The leading dimension of the array X\&.  LDX >= max(1,N)\&.
.fi
.PP
.br
\fIFERR\fP 
.PP
.nf
          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\&.
.fi
.PP
.br
\fIBERR\fP 
.PP
.nf
          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)\&.
.fi
.PP
.br
\fIWORK\fP 
.PP
.nf
          WORK is DOUBLE PRECISION array, dimension (3*N)
.fi
.PP
.br
\fIIWORK\fP 
.PP
.nf
          IWORK is INTEGER array, dimension (N)
.fi
.PP
.br
\fIINFO\fP 
.PP
.nf
          INFO is INTEGER
          = 0:  successful exit
          < 0:  if INFO = -i, the i-th argument had an illegal value
.fi
.PP
 
.RE
.PP
\fBInternal Parameters:\fP
.RS 4

.PP
.nf
  ITMAX is the maximum number of steps of iterative refinement\&.
.fi
.PP
 
.RE
.PP
\fBAuthor\fP
.RS 4
Univ\&. of Tennessee 
.PP
Univ\&. of California Berkeley 
.PP
Univ\&. of Colorado Denver 
.PP
NAG Ltd\&. 
.RE
.PP

.SS "subroutine dporfsx (character UPLO, character EQUED, integer N, integer NRHS, double precision, dimension( lda, * ) A, integer LDA, double precision, dimension( ldaf, * ) AF, integer LDAF, double precision, dimension( * ) S, double precision, dimension( ldb, * ) B, integer LDB, double precision, dimension( ldx, * ) X, integer LDX, double precision RCOND, 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)"

.PP
\fBDPORFSX\fP  
.PP
\fBPurpose:\fP
.RS 4

.PP
.nf
    DPORFSX improves the computed solution to a system of linear
    equations when the coefficient matrix is symmetric positive
    definite, and provides error bounds and backward error estimates
    for the solution\&.  In addition to normwise error bound, the code
    provides maximum componentwise error bound if possible\&.  See
    comments for ERR_BNDS_NORM and ERR_BNDS_COMP for details of the
    error bounds\&.

    The original system of linear equations may have been equilibrated
    before calling this routine, as described by arguments EQUED and S
    below\&. In this case, the solution and error bounds returned are
    for the original unequilibrated system\&.
.fi
.PP
 
.PP
.nf
     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\&.
.fi
.PP
.RE
.PP
\fBParameters\fP
.RS 4
\fIUPLO\fP 
.PP
.nf
          UPLO is CHARACTER*1
       = 'U':  Upper triangle of A is stored;
       = 'L':  Lower triangle of A is stored\&.
.fi
.PP
.br
\fIEQUED\fP 
.PP
.nf
          EQUED is CHARACTER*1
     Specifies the form of equilibration that was done to A
     before calling this routine\&. This is needed to compute
     the solution and error bounds correctly\&.
       = 'N':  No equilibration
       = 'Y':  Both row and column equilibration, i\&.e\&., A has been
               replaced by diag(S) * A * diag(S)\&.
               The right hand side B has been changed accordingly\&.
.fi
.PP
.br
\fIN\fP 
.PP
.nf
          N is INTEGER
     The order of the matrix A\&.  N >= 0\&.
.fi
.PP
.br
\fINRHS\fP 
.PP
.nf
          NRHS is INTEGER
     The number of right hand sides, i\&.e\&., the number of columns
     of the matrices B and X\&.  NRHS >= 0\&.
.fi
.PP
.br
\fIA\fP 
.PP
.nf
          A is DOUBLE PRECISION array, dimension (LDA,N)
     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\&.
.fi
.PP
.br
\fILDA\fP 
.PP
.nf
          LDA is INTEGER
     The leading dimension of the array A\&.  LDA >= max(1,N)\&.
.fi
.PP
.br
\fIAF\fP 
.PP
.nf
          AF is DOUBLE PRECISION array, dimension (LDAF,N)
     The triangular factor U or L from the Cholesky factorization
     A = U**T*U or A = L*L**T, as computed by DPOTRF\&.
.fi
.PP
.br
\fILDAF\fP 
.PP
.nf
          LDAF is INTEGER
     The leading dimension of the array AF\&.  LDAF >= max(1,N)\&.
.fi
.PP
.br
\fIS\fP 
.PP
.nf
          S is DOUBLE PRECISION array, dimension (N)
     The 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\&.
.fi
.PP
.br
\fIB\fP 
.PP
.nf
          B is DOUBLE PRECISION array, dimension (LDB,NRHS)
     The right hand side matrix B\&.
.fi
.PP
.br
\fILDB\fP 
.PP
.nf
          LDB is INTEGER
     The leading dimension of the array B\&.  LDB >= max(1,N)\&.
.fi
.PP
.br
\fIX\fP 
.PP
.nf
          X is DOUBLE PRECISION array, dimension (LDX,NRHS)
     On entry, the solution matrix X, as computed by DGETRS\&.
     On exit, the improved solution matrix X\&.
.fi
.PP
.br
\fILDX\fP 
.PP
.nf
          LDX is INTEGER
     The leading dimension of the array X\&.  LDX >= max(1,N)\&.
.fi
.PP
.br
\fIRCOND\fP 
.PP
.nf
          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\&.
.fi
.PP
.br
\fIBERR\fP 
.PP
.nf
          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)\&.
.fi
.PP
.br
\fIN_ERR_BNDS\fP 
.PP
.nf
          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\&.
.fi
.PP
.br
\fIERR_BNDS_NORM\fP 
.PP
.nf
          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\&.
.fi
.PP
.br
\fIERR_BNDS_COMP\fP 
.PP
.nf
          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 < 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\&.
.fi
.PP
.br
\fINPARAMS\fP 
.PP
.nf
          NPARAMS is INTEGER
     Specifies the number of parameters set in PARAMS\&.  If <= 0, the
     PARAMS array is never referenced and default values are used\&.
.fi
.PP
.br
\fIPARAMS\fP 
.PP
.nf
          PARAMS is DOUBLE PRECISION array, dimension (NPARAMS)
     Specifies algorithm parameters\&.  If an entry is < 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 double-precision refinement algorithm,
                    possibly with doubled-single computations if the
                    compilation environment does not support DOUBLE
                    PRECISION\&.
              (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)
.fi
.PP
.br
\fIWORK\fP 
.PP
.nf
          WORK is DOUBLE PRECISION array, dimension (4*N)
.fi
.PP
.br
\fIIWORK\fP 
.PP
.nf
          IWORK is INTEGER array, dimension (N)
.fi
.PP
.br
\fIINFO\fP 
.PP
.nf
          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\&.
.fi
.PP
 
.RE
.PP
\fBAuthor\fP
.RS 4
Univ\&. of Tennessee 
.PP
Univ\&. of California Berkeley 
.PP
Univ\&. of Colorado Denver 
.PP
NAG Ltd\&. 
.RE
.PP

.SS "subroutine dpotf2 (character UPLO, integer N, double precision, dimension( lda, * ) A, integer LDA, integer INFO)"

.PP
\fBDPOTF2\fP computes the Cholesky factorization of a symmetric/Hermitian positive definite matrix (unblocked algorithm)\&.  
.PP
\fBPurpose:\fP
.RS 4

.PP
.nf
 DPOTF2 computes the Cholesky factorization of a real symmetric
 positive definite matrix A\&.

 The factorization has the form
    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 lower triangular\&.

 This is the unblocked version of the algorithm, calling Level 2 BLAS\&.
.fi
.PP
 
.RE
.PP
\fBParameters\fP
.RS 4
\fIUPLO\fP 
.PP
.nf
          UPLO is CHARACTER*1
          Specifies whether the upper or lower triangular part of the
          symmetric matrix A is stored\&.
          = 'U':  Upper triangular
          = 'L':  Lower triangular
.fi
.PP
.br
\fIN\fP 
.PP
.nf
          N is INTEGER
          The order of the matrix A\&.  N >= 0\&.
.fi
.PP
.br
\fIA\fP 
.PP
.nf
          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\&.
.fi
.PP
.br
\fILDA\fP 
.PP
.nf
          LDA is INTEGER
          The leading dimension of the array A\&.  LDA >= max(1,N)\&.
.fi
.PP
.br
\fIINFO\fP 
.PP
.nf
          INFO is INTEGER
          = 0: successful exit
          < 0: if INFO = -k, the k-th argument had an illegal value
          > 0: if INFO = k, the leading minor of order k is not
               positive definite, and the factorization could not be
               completed\&.
.fi
.PP
 
.RE
.PP
\fBAuthor\fP
.RS 4
Univ\&. of Tennessee 
.PP
Univ\&. of California Berkeley 
.PP
Univ\&. of Colorado Denver 
.PP
NAG Ltd\&. 
.RE
.PP

.SS "subroutine dpotrf (character UPLO, integer N, double precision, dimension( lda, * ) A, integer LDA, integer INFO)"

.PP
\fBDPOTRF\fP \fBDPOTRF\fP VARIANT: top-looking block version of the algorithm, calling Level 3 BLAS\&.
.PP
 
.PP
\fBPurpose:\fP
.RS 4

.PP
.nf
 DPOTRF computes the Cholesky factorization of a real symmetric
 positive definite matrix A\&.

 The factorization has the form
    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 lower triangular\&.

 This is the block version of the algorithm, calling Level 3 BLAS\&.
.fi
.PP
 
.RE
.PP
\fBParameters\fP
.RS 4
\fIUPLO\fP 
.PP
.nf
          UPLO is CHARACTER*1
          = 'U':  Upper triangle of A is stored;
          = 'L':  Lower triangle of A is stored\&.
.fi
.PP
.br
\fIN\fP 
.PP
.nf
          N is INTEGER
          The order of the matrix A\&.  N >= 0\&.
.fi
.PP
.br
\fIA\fP 
.PP
.nf
          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\&.
.fi
.PP
.br
\fILDA\fP 
.PP
.nf
          LDA is INTEGER
          The leading dimension of the array A\&.  LDA >= max(1,N)\&.
.fi
.PP
.br
\fIINFO\fP 
.PP
.nf
          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 is not
                positive definite, and the factorization could not be
                completed\&.
.fi
.PP
 
.RE
.PP
\fBAuthor\fP
.RS 4
Univ\&. of Tennessee 
.PP
Univ\&. of California Berkeley 
.PP
Univ\&. of Colorado Denver 
.PP
NAG Ltd\&.
.RE
.PP
\fBPurpose:\fP 
.PP
.nf
 DPOTRF computes the Cholesky factorization of a real symmetric
 positive definite matrix A\&.

 The factorization has the form
    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 lower triangular\&.

 This is the top-looking block version of the algorithm, calling Level 3 BLAS\&.
.fi
.PP
 
.PP
\fBParameters\fP
.RS 4
\fIUPLO\fP 
.PP
.nf
          UPLO is CHARACTER*1
          = 'U':  Upper triangle of A is stored;
          = 'L':  Lower triangle of A is stored\&.
.fi
.PP
.br
\fIN\fP 
.PP
.nf
          N is INTEGER
          The order of the matrix A\&.  N >= 0\&.
.fi
.PP
.br
\fIA\fP 
.PP
.nf
          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\&.
.fi
.PP
 
.PP
.nf
          On exit, if INFO = 0, the factor U or L from the Cholesky
          factorization A = U**T*U or A = L*L**T\&.
.fi
.PP
.br
\fILDA\fP 
.PP
.nf
          LDA is INTEGER
          The leading dimension of the array A\&.  LDA >= max(1,N)\&.
.fi
.PP
.br
\fIINFO\fP 
.PP
.nf
          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 is not
                positive definite, and the factorization could not be
                completed\&.
.fi
.PP
 
.RE
.PP
\fBAuthor\fP
.RS 4
Univ\&. of Tennessee 
.PP
Univ\&. of California Berkeley 
.PP
Univ\&. of Colorado Denver 
.PP
NAG Ltd\&. 
.RE
.PP
\fBDate\fP
.RS 4
December 2016 
.RE
.PP

.SS "recursive subroutine dpotrf2 (character UPLO, integer N, double precision, dimension( lda, * ) A, integer LDA, integer INFO)"

.PP
\fBDPOTRF2\fP 
.PP
\fBPurpose:\fP
.RS 4

.PP
.nf
 DPOTRF2 computes the Cholesky factorization of a real symmetric
 positive definite matrix A using the recursive algorithm\&.

 The factorization has the form
    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 lower triangular\&.

 This is the recursive version of the algorithm\&. It divides
 the matrix into four submatrices:

        [  A11 | A12  ]  where A11 is n1 by n1 and A22 is n2 by n2
    A = [ -----|----- ]  with n1 = n/2
        [  A21 | A22  ]       n2 = n-n1

 The subroutine calls itself to factor A11\&. Update and scale A21
 or A12, update A22 then calls itself to factor A22\&.
.fi
.PP
 
.RE
.PP
\fBParameters\fP
.RS 4
\fIUPLO\fP 
.PP
.nf
          UPLO is CHARACTER*1
          = 'U':  Upper triangle of A is stored;
          = 'L':  Lower triangle of A is stored\&.
.fi
.PP
.br
\fIN\fP 
.PP
.nf
          N is INTEGER
          The order of the matrix A\&.  N >= 0\&.
.fi
.PP
.br
\fIA\fP 
.PP
.nf
          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\&.
.fi
.PP
.br
\fILDA\fP 
.PP
.nf
          LDA is INTEGER
          The leading dimension of the array A\&.  LDA >= max(1,N)\&.
.fi
.PP
.br
\fIINFO\fP 
.PP
.nf
          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 is not
                positive definite, and the factorization could not be
                completed\&.
.fi
.PP
 
.RE
.PP
\fBAuthor\fP
.RS 4
Univ\&. of Tennessee 
.PP
Univ\&. of California Berkeley 
.PP
Univ\&. of Colorado Denver 
.PP
NAG Ltd\&. 
.RE
.PP

.SS "subroutine dpotri (character UPLO, integer N, double precision, dimension( lda, * ) A, integer LDA, integer INFO)"

.PP
\fBDPOTRI\fP  
.PP
\fBPurpose:\fP
.RS 4

.PP
.nf
 DPOTRI computes the inverse of a real symmetric positive definite
 matrix A using the Cholesky factorization A = U**T*U or A = L*L**T
 computed by DPOTRF\&.
.fi
.PP
 
.RE
.PP
\fBParameters\fP
.RS 4
\fIUPLO\fP 
.PP
.nf
          UPLO is CHARACTER*1
          = 'U':  Upper triangle of A is stored;
          = 'L':  Lower triangle of A is stored\&.
.fi
.PP
.br
\fIN\fP 
.PP
.nf
          N is INTEGER
          The order of the matrix A\&.  N >= 0\&.
.fi
.PP
.br
\fIA\fP 
.PP
.nf
          A is DOUBLE PRECISION array, dimension (LDA,N)
          On entry, the triangular factor U or L from the Cholesky
          factorization A = U**T*U or A = L*L**T, as computed by
          DPOTRF\&.
          On exit, the upper or lower triangle of the (symmetric)
          inverse of A, overwriting the input factor U or L\&.
.fi
.PP
.br
\fILDA\fP 
.PP
.nf
          LDA is INTEGER
          The leading dimension of the array A\&.  LDA >= max(1,N)\&.
.fi
.PP
.br
\fIINFO\fP 
.PP
.nf
          INFO is INTEGER
          = 0:  successful exit
          < 0:  if INFO = -i, the i-th argument had an illegal value
          > 0:  if INFO = i, the (i,i) element of the factor U or L is
                zero, and the inverse could not be computed\&.
.fi
.PP
 
.RE
.PP
\fBAuthor\fP
.RS 4
Univ\&. of Tennessee 
.PP
Univ\&. of California Berkeley 
.PP
Univ\&. of Colorado Denver 
.PP
NAG Ltd\&. 
.RE
.PP

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

.PP
\fBDPOTRS\fP  
.PP
\fBPurpose:\fP
.RS 4

.PP
.nf
 DPOTRS solves a system of linear equations A*X = B with a symmetric
 positive definite matrix A using the Cholesky factorization
 A = U**T*U or A = L*L**T computed by DPOTRF\&.
.fi
.PP
 
.RE
.PP
\fBParameters\fP
.RS 4
\fIUPLO\fP 
.PP
.nf
          UPLO is CHARACTER*1
          = 'U':  Upper triangle of A is stored;
          = 'L':  Lower triangle of A is stored\&.
.fi
.PP
.br
\fIN\fP 
.PP
.nf
          N is INTEGER
          The order of the matrix A\&.  N >= 0\&.
.fi
.PP
.br
\fINRHS\fP 
.PP
.nf
          NRHS is INTEGER
          The number of right hand sides, i\&.e\&., the number of columns
          of the matrix B\&.  NRHS >= 0\&.
.fi
.PP
.br
\fIA\fP 
.PP
.nf
          A is DOUBLE PRECISION array, dimension (LDA,N)
          The triangular factor U or L from the Cholesky factorization
          A = U**T*U or A = L*L**T, as computed by DPOTRF\&.
.fi
.PP
.br
\fILDA\fP 
.PP
.nf
          LDA is INTEGER
          The leading dimension of the array A\&.  LDA >= max(1,N)\&.
.fi
.PP
.br
\fIB\fP 
.PP
.nf
          B is DOUBLE PRECISION array, dimension (LDB,NRHS)
          On entry, the right hand side matrix B\&.
          On exit, the solution matrix X\&.
.fi
.PP
.br
\fILDB\fP 
.PP
.nf
          LDB is INTEGER
          The leading dimension of the array B\&.  LDB >= max(1,N)\&.
.fi
.PP
.br
\fIINFO\fP 
.PP
.nf
          INFO is INTEGER
          = 0:  successful exit
          < 0:  if INFO = -i, the i-th argument had an illegal value
.fi
.PP
 
.RE
.PP
\fBAuthor\fP
.RS 4
Univ\&. of Tennessee 
.PP
Univ\&. of California Berkeley 
.PP
Univ\&. of Colorado Denver 
.PP
NAG Ltd\&. 
.RE
.PP

.SH "Author"
.PP 
Generated automatically by Doxygen for LAPACK from the source code\&.