COMPLEX routines for symmetric or Hermitian positive definite matrix
- cpocon : Estimates the reciprocal of the condition number of aHermitian positive definite matrix, using theCholesky factorization computed by CPOTRF.
- cpoequ : Computes row and column scalings to equilibrate a Hermitianpositive definite matrix and reduce its condition number.
- cporfs : Improves the computed solution to a Hermitian positivedefinite system of linear equations AX=B, and provides forwardand backward error bounds for the solution.
- cposv : Solves a Hermitian positive definite system of linearequations AX=B.
- cposvx : Solves a Hermitian positive definite system of linearequations AX=B, and provides an estimate of the condition numberand error bounds on the solution.
- cpotf2 :
- cpotrf : Computes the Cholesky factorization of a Hermitianpositive definite matrix.
- cpotri : Computes the inverse of a Hermitian positive definitematrix, using the Cholesky factorization computed by CPOTRF.
- cpotrs : Solves a Hermitian positive definite system of linearequations AX=B, using the Cholesky factorization computed byCPOTRF.
cpocon
Estimates the reciprocal of the condition number of aHermitian positive definite matrix, using theCholesky factorization computed by CPOTRF.
USAGE:
rcond, info = NumRu::Lapack.cpocon( uplo, a, anorm)
or
NumRu::Lapack.cpocon # print help
FORTRAN MANUAL
SUBROUTINE CPOCON( UPLO, N, A, LDA, ANORM, RCOND, WORK, RWORK, INFO )
* Purpose
* =======
*
* CPOCON estimates the reciprocal of the condition number (in the
* 1-norm) of a complex Hermitian positive definite matrix using the
* Cholesky factorization A = U**H*U or A = L*L**H computed by CPOTRF.
*
* An estimate is obtained for norm(inv(A)), and the reciprocal of the
* condition number is computed as RCOND = 1 / (ANORM * norm(inv(A))).
*
* Arguments
* =========
*
* UPLO (input) CHARACTER*1
* = 'U': Upper triangle of A is stored;
* = 'L': Lower triangle of A is stored.
*
* N (input) INTEGER
* The order of the matrix A. N >= 0.
*
* A (input) COMPLEX array, dimension (LDA,N)
* The triangular factor U or L from the Cholesky factorization
* A = U**H*U or A = L*L**H, as computed by CPOTRF.
*
* LDA (input) INTEGER
* The leading dimension of the array A. LDA >= max(1,N).
*
* ANORM (input) REAL
* The 1-norm (or infinity-norm) of the Hermitian matrix A.
*
* RCOND (output) REAL
* 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.
*
* WORK (workspace) COMPLEX array, dimension (2*N)
*
* RWORK (workspace) REAL array, dimension (N)
*
* INFO (output) INTEGER
* = 0: successful exit
* < 0: if INFO = -i, the i-th argument had an illegal value
*
* =====================================================================
*
go to the page top
cpoequ
Computes row and column scalings to equilibrate a Hermitianpositive definite matrix and reduce its condition number.
USAGE:
s, scond, amax, info = NumRu::Lapack.cpoequ( a)
or
NumRu::Lapack.cpoequ # print help
FORTRAN MANUAL
SUBROUTINE CPOEQU( N, A, LDA, S, SCOND, AMAX, INFO )
* Purpose
* =======
*
* CPOEQU computes row and column scalings intended to equilibrate a
* Hermitian 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.
*
* Arguments
* =========
*
* N (input) INTEGER
* The order of the matrix A. N >= 0.
*
* A (input) COMPLEX array, dimension (LDA,N)
* The N-by-N Hermitian positive definite matrix whose scaling
* factors are to be computed. Only the diagonal elements of A
* are referenced.
*
* LDA (input) INTEGER
* The leading dimension of the array A. LDA >= max(1,N).
*
* S (output) REAL array, dimension (N)
* If INFO = 0, S contains the scale factors for A.
*
* SCOND (output) REAL
* 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.
*
* AMAX (output) REAL
* Absolute value of largest matrix element. If AMAX is very
* close to overflow or very close to underflow, the matrix
* should be scaled.
*
* INFO (output) 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.
*
* =====================================================================
*
go to the page top
cporfs
Improves the computed solution to a Hermitian positivedefinite system of linear equations AX=B, and provides forwardand backward error bounds for the solution.
USAGE:
ferr, berr, info, x = NumRu::Lapack.cporfs( uplo, a, af, b, x)
or
NumRu::Lapack.cporfs # print help
FORTRAN MANUAL
SUBROUTINE CPORFS( UPLO, N, NRHS, A, LDA, AF, LDAF, B, LDB, X, LDX, FERR, BERR, WORK, RWORK, INFO )
* Purpose
* =======
*
* CPORFS improves the computed solution to a system of linear
* equations when the coefficient matrix is Hermitian positive definite,
* and provides error bounds and backward error estimates for the
* solution.
*
* Arguments
* =========
*
* UPLO (input) CHARACTER*1
* = 'U': Upper triangle of A is stored;
* = 'L': Lower triangle of A is stored.
*
* N (input) INTEGER
* The order of the matrix A. N >= 0.
*
* NRHS (input) INTEGER
* The number of right hand sides, i.e., the number of columns
* of the matrices B and X. NRHS >= 0.
*
* A (input) COMPLEX array, dimension (LDA,N)
* The Hermitian 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.
*
* LDA (input) INTEGER
* The leading dimension of the array A. LDA >= max(1,N).
*
* AF (input) COMPLEX array, dimension (LDAF,N)
* The triangular factor U or L from the Cholesky factorization
* A = U**H*U or A = L*L**H, as computed by CPOTRF.
*
* LDAF (input) INTEGER
* The leading dimension of the array AF. LDAF >= max(1,N).
*
* B (input) COMPLEX array, dimension (LDB,NRHS)
* The right hand side matrix B.
*
* LDB (input) INTEGER
* The leading dimension of the array B. LDB >= max(1,N).
*
* X (input/output) COMPLEX array, dimension (LDX,NRHS)
* On entry, the solution matrix X, as computed by CPOTRS.
* On exit, the improved solution matrix X.
*
* LDX (input) INTEGER
* The leading dimension of the array X. LDX >= max(1,N).
*
* FERR (output) REAL 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.
*
* BERR (output) REAL 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).
*
* WORK (workspace) COMPLEX array, dimension (2*N)
*
* RWORK (workspace) REAL array, dimension (N)
*
* INFO (output) INTEGER
* = 0: successful exit
* < 0: if INFO = -i, the i-th argument had an illegal value
*
* Internal Parameters
* ===================
*
* ITMAX is the maximum number of steps of iterative refinement.
*
* ====================================================================
*
go to the page top
cposv
Solves a Hermitian positive definite system of linearequations AX=B.
USAGE:
info, a, b = NumRu::Lapack.cposv( uplo, a, b)
or
NumRu::Lapack.cposv # print help
FORTRAN MANUAL
SUBROUTINE CPOSV( UPLO, N, NRHS, A, LDA, B, LDB, INFO )
* Purpose
* =======
*
* CPOSV computes the solution to a complex system of linear equations
* A * X = B,
* where A is an N-by-N Hermitian positive definite matrix and X and B
* are N-by-NRHS matrices.
*
* The Cholesky decomposition is used to factor A as
* A = U**H* U, if UPLO = 'U', or
* A = L * L**H, 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.
*
* Arguments
* =========
*
* UPLO (input) CHARACTER*1
* = 'U': Upper triangle of A is stored;
* = 'L': Lower triangle of A is stored.
*
* N (input) INTEGER
* The number of linear equations, i.e., the order of the
* matrix A. N >= 0.
*
* NRHS (input) INTEGER
* The number of right hand sides, i.e., the number of columns
* of the matrix B. NRHS >= 0.
*
* A (input/output) COMPLEX array, dimension (LDA,N)
* On entry, the Hermitian 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**H*U or A = L*L**H.
*
* LDA (input) INTEGER
* The leading dimension of the array A. LDA >= max(1,N).
*
* B (input/output) COMPLEX 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.
*
* LDB (input) INTEGER
* The leading dimension of the array B. LDB >= max(1,N).
*
* INFO (output) 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.
*
* =====================================================================
*
* .. External Functions ..
LOGICAL LSAME
EXTERNAL LSAME
* ..
* .. External Subroutines ..
EXTERNAL CPOTRF, CPOTRS, XERBLA
* ..
* .. Intrinsic Functions ..
INTRINSIC MAX
* ..
go to the page top
cposvx
Solves a Hermitian positive definite system of linearequations AX=B, and provides an estimate of the condition numberand error bounds on the solution.
USAGE:
x, rcond, ferr, berr, info, a, af, equed, s, b = NumRu::Lapack.cposvx( fact, uplo, a, af, equed, s, b)
or
NumRu::Lapack.cposvx # print help
FORTRAN MANUAL
SUBROUTINE CPOSVX( FACT, UPLO, N, NRHS, A, LDA, AF, LDAF, EQUED, S, B, LDB, X, LDX, RCOND, FERR, BERR, WORK, RWORK, INFO )
* Purpose
* =======
*
* CPOSVX uses the Cholesky factorization A = U**H*U or A = L*L**H to
* compute the solution to a complex system of linear equations
* A * X = B,
* where A is an N-by-N Hermitian 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**H* U, if UPLO = 'U', or
* A = L * L**H, 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.
*
* Arguments
* =========
*
* FACT (input) 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.
*
* UPLO (input) CHARACTER*1
* = 'U': Upper triangle of A is stored;
* = 'L': Lower triangle of A is stored.
*
* N (input) INTEGER
* The number of linear equations, i.e., the order of the
* matrix A. N >= 0.
*
* NRHS (input) INTEGER
* The number of right hand sides, i.e., the number of columns
* of the matrices B and X. NRHS >= 0.
*
* A (input/output) COMPLEX array, dimension (LDA,N)
* On entry, the Hermitian 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).
*
* LDA (input) INTEGER
* The leading dimension of the array A. LDA >= max(1,N).
*
* AF (input or output) COMPLEX 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**H*U or A = L*L**H, 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**H*U or A = L*L**H 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**H*U or A = L*L**H of the equilibrated
* matrix A (see the description of A for the form of the
* equilibrated matrix).
*
* LDAF (input) INTEGER
* The leading dimension of the array AF. LDAF >= max(1,N).
*
* EQUED (input or output) 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.
*
* S (input or output) REAL 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.
*
* B (input/output) COMPLEX array, dimension (LDB,NRHS)
* On entry, the N-by-NRHS righthand side matrix B.
* On exit, if EQUED = 'N', B is not modified; if EQUED = 'Y',
* B is overwritten by diag(S) * B.
*
* LDB (input) INTEGER
* The leading dimension of the array B. LDB >= max(1,N).
*
* X (output) COMPLEX 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.
*
* LDX (input) INTEGER
* The leading dimension of the array X. LDX >= max(1,N).
*
* RCOND (output) REAL
* 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.
*
* FERR (output) REAL 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.
*
* BERR (output) REAL 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).
*
* WORK (workspace) COMPLEX array, dimension (2*N)
*
* RWORK (workspace) REAL array, dimension (N)
*
* INFO (output) 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.
*
* =====================================================================
*
go to the page top
cpotf2
USAGE:
info, a = NumRu::Lapack.cpotf2( uplo, a)
or
NumRu::Lapack.cpotf2 # print help
FORTRAN MANUAL
SUBROUTINE CPOTF2( UPLO, N, A, LDA, INFO )
* Purpose
* =======
*
* CPOTF2 computes the Cholesky factorization of a complex Hermitian
* positive definite matrix A.
*
* The factorization has the form
* A = U' * U , if UPLO = 'U', or
* A = L * L', 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.
*
* Arguments
* =========
*
* UPLO (input) CHARACTER*1
* Specifies whether the upper or lower triangular part of the
* Hermitian matrix A is stored.
* = 'U': Upper triangular
* = 'L': Lower triangular
*
* N (input) INTEGER
* The order of the matrix A. N >= 0.
*
* A (input/output) COMPLEX array, dimension (LDA,N)
* On entry, the Hermitian 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'*U or A = L*L'.
*
* LDA (input) INTEGER
* The leading dimension of the array A. LDA >= max(1,N).
*
* INFO (output) 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.
*
* =====================================================================
*
go to the page top
cpotrf
Computes the Cholesky factorization of a Hermitianpositive definite matrix.
USAGE:
info, a = NumRu::Lapack.cpotrf( uplo, a)
or
NumRu::Lapack.cpotrf # print help
FORTRAN MANUAL
SUBROUTINE CPOTRF( UPLO, N, A, LDA, INFO )
* Purpose
* =======
*
* CPOTRF computes the Cholesky factorization of a complex Hermitian
* positive definite matrix A.
*
* The factorization has the form
* A = U**H * U, if UPLO = 'U', or
* A = L * L**H, 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.
*
* Arguments
* =========
*
* UPLO (input) CHARACTER*1
* = 'U': Upper triangle of A is stored;
* = 'L': Lower triangle of A is stored.
*
* N (input) INTEGER
* The order of the matrix A. N >= 0.
*
* A (input/output) COMPLEX array, dimension (LDA,N)
* On entry, the Hermitian 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**H*U or A = L*L**H.
*
* LDA (input) INTEGER
* The leading dimension of the array A. LDA >= max(1,N).
*
* INFO (output) 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.
*
* =====================================================================
*
go to the page top
cpotri
Computes the inverse of a Hermitian positive definitematrix, using the Cholesky factorization computed by CPOTRF.
USAGE:
info, a = NumRu::Lapack.cpotri( uplo, a)
or
NumRu::Lapack.cpotri # print help
FORTRAN MANUAL
SUBROUTINE CPOTRI( UPLO, N, A, LDA, INFO )
* Purpose
* =======
*
* CPOTRI computes the inverse of a complex Hermitian positive definite
* matrix A using the Cholesky factorization A = U**H*U or A = L*L**H
* computed by CPOTRF.
*
* Arguments
* =========
*
* UPLO (input) CHARACTER*1
* = 'U': Upper triangle of A is stored;
* = 'L': Lower triangle of A is stored.
*
* N (input) INTEGER
* The order of the matrix A. N >= 0.
*
* A (input/output) COMPLEX array, dimension (LDA,N)
* On entry, the triangular factor U or L from the Cholesky
* factorization A = U**H*U or A = L*L**H, as computed by
* CPOTRF.
* On exit, the upper or lower triangle of the (Hermitian)
* inverse of A, overwriting the input factor U or L.
*
* LDA (input) INTEGER
* The leading dimension of the array A. LDA >= max(1,N).
*
* INFO (output) 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.
*
* =====================================================================
*
* .. External Functions ..
LOGICAL LSAME
EXTERNAL LSAME
* ..
* .. External Subroutines ..
EXTERNAL CLAUUM, CTRTRI, XERBLA
* ..
* .. Intrinsic Functions ..
INTRINSIC MAX
* ..
go to the page top
cpotrs
Solves a Hermitian positive definite system of linearequations AX=B, using the Cholesky factorization computed byCPOTRF.
USAGE:
info, b = NumRu::Lapack.cpotrs( uplo, a, b)
or
NumRu::Lapack.cpotrs # print help
FORTRAN MANUAL
SUBROUTINE CPOTRS( UPLO, N, NRHS, A, LDA, B, LDB, INFO )
* Purpose
* =======
*
* CPOTRS solves a system of linear equations A*X = B with a Hermitian
* positive definite matrix A using the Cholesky factorization
* A = U**H*U or A = L*L**H computed by CPOTRF.
*
* Arguments
* =========
*
* UPLO (input) CHARACTER*1
* = 'U': Upper triangle of A is stored;
* = 'L': Lower triangle of A is stored.
*
* N (input) INTEGER
* The order of the matrix A. N >= 0.
*
* NRHS (input) INTEGER
* The number of right hand sides, i.e., the number of columns
* of the matrix B. NRHS >= 0.
*
* A (input) COMPLEX array, dimension (LDA,N)
* The triangular factor U or L from the Cholesky factorization
* A = U**H*U or A = L*L**H, as computed by CPOTRF.
*
* LDA (input) INTEGER
* The leading dimension of the array A. LDA >= max(1,N).
*
* B (input/output) COMPLEX array, dimension (LDB,NRHS)
* On entry, the right hand side matrix B.
* On exit, the solution matrix X.
*
* LDB (input) INTEGER
* The leading dimension of the array B. LDB >= max(1,N).
*
* INFO (output) INTEGER
* = 0: successful exit
* < 0: if INFO = -i, the i-th argument had an illegal value
*
* =====================================================================
*
go to the page top
back to matrix types