PLASMA
Parallel Linear Algebra Software for Multicore Architectures
posv: Solves Ax = b using Cholesky factorization (driver)

Functions

int plasma_cposv (plasma_enum_t uplo, int n, int nrhs, plasma_complex32_t *pA, int lda, plasma_complex32_t *pB, int ldb)
 
void plasma_omp_cposv (plasma_enum_t uplo, plasma_desc_t A, plasma_desc_t B, plasma_sequence_t *sequence, plasma_request_t *request)
 
int plasma_dposv (plasma_enum_t uplo, int n, int nrhs, double *pA, int lda, double *pB, int ldb)
 
void plasma_omp_dposv (plasma_enum_t uplo, plasma_desc_t A, plasma_desc_t B, plasma_sequence_t *sequence, plasma_request_t *request)
 
int plasma_dsposv (plasma_enum_t uplo, int n, int nrhs, double *pA, int lda, double *pB, int ldb, double *pX, int ldx, int *iter)
 
void plasma_omp_dsposv (plasma_enum_t uplo, plasma_desc_t A, plasma_desc_t B, plasma_desc_t X, plasma_desc_t As, plasma_desc_t Xs, plasma_desc_t R, double *work, double *Rnorm, double *Xnorm, int *iter, plasma_sequence_t *sequence, plasma_request_t *request)
 
int plasma_sposv (plasma_enum_t uplo, int n, int nrhs, float *pA, int lda, float *pB, int ldb)
 
void plasma_omp_sposv (plasma_enum_t uplo, plasma_desc_t A, plasma_desc_t B, plasma_sequence_t *sequence, plasma_request_t *request)
 
int plasma_zcposv (plasma_enum_t uplo, int n, int nrhs, plasma_complex64_t *pA, int lda, plasma_complex64_t *pB, int ldb, plasma_complex64_t *pX, int ldx, int *iter)
 
void plasma_omp_zcposv (plasma_enum_t uplo, plasma_desc_t A, plasma_desc_t B, plasma_desc_t X, plasma_desc_t As, plasma_desc_t Xs, plasma_desc_t R, double *work, double *Rnorm, double *Xnorm, int *iter, plasma_sequence_t *sequence, plasma_request_t *request)
 
int plasma_zposv (plasma_enum_t uplo, int n, int nrhs, plasma_complex64_t *pA, int lda, plasma_complex64_t *pB, int ldb)
 
void plasma_omp_zposv (plasma_enum_t uplo, plasma_desc_t A, plasma_desc_t B, plasma_sequence_t *sequence, plasma_request_t *request)
 

Detailed Description

Function Documentation

int plasma_cposv ( plasma_enum_t  uplo,
int  n,
int  nrhs,
plasma_complex32_t *  pA,
int  lda,
plasma_complex32_t *  pB,
int  ldb 
)

Computes the solution to a 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 = L\times L^H, \]

if uplo = PlasmaLower, or

\[ A = U^H\times U, \]

if uplo = PlasmaUpper,

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
  • PlasmaUpper: Upper triangle of A is stored;
  • PlasmaLower: Lower triangle of A is stored.
[in]nThe number of linear equations, i.e., the order of the matrix A. n >= 0.
[in]nrhsThe number of right hand sides, i.e., the number of columns of the matrix B. nrhs >= 0.
[in,out]pAOn entry, the Hermitian positive definite matrix A. If uplo = PlasmaUpper, 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 return value = 0, the factor U or L from the Cholesky factorization A = U^H*U or A = L*L^H.
[in]ldaThe leading dimension of the array A. lda >= max(1,n).
[in,out]pBOn entry, the n-by-nrhs right hand side matrix B. On exit, if return value = 0, the n-by-nrhs solution matrix X.
[in]ldbThe leading dimension of the array B. ldb >= max(1,n).
Return values
PlasmaSuccesssuccessful exit
<0 if -i, the i-th argument had an illegal value
>0 if 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.
See also
plasma_omp_cposv
plasma_cposv
plasma_dposv
plasma_sposv
void plasma_omp_cposv ( plasma_enum_t  uplo,
plasma_desc_t  A,
plasma_desc_t  B,
plasma_sequence_t *  sequence,
plasma_request_t *  request 
)

Solves a Hermitian positive definite system of linear equations using Cholesky factorization. Non-blocking tile version of plasma_cposv(). Operates on matrices stored by tiles. All matrices are passed through descriptors. All dimensions are taken from the descriptors. Allows for pipelining of operations at runtime.

Parameters
[in]uplo
  • PlasmaUpper: Upper triangle of A is stored;
  • PlasmaLower: Lower triangle of A is stored.
[in,out]AOn entry, the Hermitian positive definite matrix A. If uplo = PlasmaUpper, 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 return value = 0, the factor U or L from the Cholesky factorization A = U^H*U or A = L*L^H.
[in,out]BOn entry, the n-by-nrhs right hand side matrix B. On exit, if return value = 0, the n-by-nrhs solution matrix X.
[in]sequenceIdentifies the sequence of function calls that this call belongs to (for completion checks and exception handling purposes). Check the sequence->status for errors.
[out]requestIdentifies this function call (for exception handling purposes).
Return values
voidErrors are returned by setting sequence->status and request->status to error values. The sequence->status and request->status should never be set to PlasmaSuccess (the initial values) since another async call may be setting a failure value at the same time.
See also
plasma_cposv
plasma_omp_cposv
plasma_omp_dposv
plasma_omp_sposv
int plasma_dposv ( plasma_enum_t  uplo,
int  n,
int  nrhs,
double *  pA,
int  lda,
double *  pB,
int  ldb 
)

Computes the solution to a 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 = L\times L^T, \]

if uplo = PlasmaLower, or

\[ A = U^T\times U, \]

if uplo = PlasmaUpper,

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
  • PlasmaUpper: Upper triangle of A is stored;
  • PlasmaLower: Lower triangle of A is stored.
[in]nThe number of linear equations, i.e., the order of the matrix A. n >= 0.
[in]nrhsThe number of right hand sides, i.e., the number of columns of the matrix B. nrhs >= 0.
[in,out]pAOn entry, the symmetric positive definite matrix A. If uplo = PlasmaUpper, 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 return value = 0, the factor U or L from the Cholesky factorization A = U^T*U or A = L*L^T.
[in]ldaThe leading dimension of the array A. lda >= max(1,n).
[in,out]pBOn entry, the n-by-nrhs right hand side matrix B. On exit, if return value = 0, the n-by-nrhs solution matrix X.
[in]ldbThe leading dimension of the array B. ldb >= max(1,n).
Return values
PlasmaSuccesssuccessful exit
<0 if -i, the i-th argument had an illegal value
>0 if 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.
See also
plasma_omp_dposv
plasma_cposv
plasma_dposv
plasma_sposv
void plasma_omp_dposv ( plasma_enum_t  uplo,
plasma_desc_t  A,
plasma_desc_t  B,
plasma_sequence_t *  sequence,
plasma_request_t *  request 
)

Solves a symmetric positive definite system of linear equations using Cholesky factorization. Non-blocking tile version of plasma_dposv(). Operates on matrices stored by tiles. All matrices are passed through descriptors. All dimensions are taken from the descriptors. Allows for pipelining of operations at runtime.

Parameters
[in]uplo
  • PlasmaUpper: Upper triangle of A is stored;
  • PlasmaLower: Lower triangle of A is stored.
[in,out]AOn entry, the symmetric positive definite matrix A. If uplo = PlasmaUpper, 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 return value = 0, the factor U or L from the Cholesky factorization A = U^T*U or A = L*L^T.
[in,out]BOn entry, the n-by-nrhs right hand side matrix B. On exit, if return value = 0, the n-by-nrhs solution matrix X.
[in]sequenceIdentifies the sequence of function calls that this call belongs to (for completion checks and exception handling purposes). Check the sequence->status for errors.
[out]requestIdentifies this function call (for exception handling purposes).
Return values
voidErrors are returned by setting sequence->status and request->status to error values. The sequence->status and request->status should never be set to PlasmaSuccess (the initial values) since another async call may be setting a failure value at the same time.
See also
plasma_dposv
plasma_omp_cposv
plasma_omp_dposv
plasma_omp_sposv
int plasma_dsposv ( plasma_enum_t  uplo,
int  n,
int  nrhs,
double *  pA,
int  lda,
double *  pB,
int  ldb,
double *  pX,
int  ldx,
int *  iter 
)

Computes the solution to a 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.

plasma_dsposv first factorizes the matrix using plasma_spotrf and uses this factorization within an iterative refinement procedure to produce a solution with COMPLEX*16 normwise backward error quality (see below). If the approach fails the method falls back to a COMPLEX*16 factorization and solve.

The iterative refinement is not going to be a winning strategy if the ratio COMPLEX performance over COMPLEX*16 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: Rnorm < sqrt(n)*Xnorm*Anorm*eps*BWDmax where:

  • iter is the number of the current iteration in the iterative refinement process
  • Rnorm is the Infinity-norm of the residual
  • Xnorm is the Infinity-norm of the solution
  • Anorm is the Infinity-operator-norm of the matrix A
  • eps is the machine epsilon returned by DLAMCH('Epsilon'). The values itermax and BWDmax are fixed to 30 and 1.0D+00 respectively.
Parameters
[in]uploSpecifies whether the matrix A is upper or lower triangular:
  • PlasmaUpper: Upper triangle of A is stored;
  • PlasmaLower: Lower triangle of A is stored.
[in]nThe number of linear equations, i.e., the order of the matrix A. n >= 0.
[in]nrhsThe number of right hand sides, i.e., the number of columns of the matrix B. nrhs >= 0.
[in,out]pAThe n-by-n symmetric positive definite coefficient matrix A. If uplo = PlasmaUpper, 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 = PlasmaLower, 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, contains the lower Cholesky factor matrix L, if uplo == PlasmaLower and upper Cholesky factor (L^T), otherwise.
[in]ldaThe leading dimension of the array A. lda >= max(1,n).
[in]pBThe n-by-nrhs matrix of right hand side matrix B. This matrix remains unchanged.
[in]ldbThe leading dimension of the array B. ldb >= max(1,n).
[out]pXIf return value = 0, the n-by-nrhs solution matrix X.
[in]ldxThe leading dimension of the array X. ldx >= max(1,n).
[out]iterThe number of the iterations in the iterative refinement process, needed for the convergence. If failed, it is set to be -(1+itermax), where itermax = 30.
Return values
PlasmaSuccesssuccessful exit
See also
plasma_omp_dsposv
plasma_dsposv
plasma_zposv
void plasma_omp_dsposv ( plasma_enum_t  uplo,
plasma_desc_t  A,
plasma_desc_t  B,
plasma_desc_t  X,
plasma_desc_t  As,
plasma_desc_t  Xs,
plasma_desc_t  R,
double *  work,
double *  Rnorm,
double *  Xnorm,
int *  iter,
plasma_sequence_t *  sequence,
plasma_request_t *  request 
)

Solves a symmetric positive definite system using iterative refinement with the Cholesky factor computed using plasma_spotrf. Non-blocking tile version of plasma_dsposv(). Operates on matrices stored by tiles. All matrices are passed through descriptors. All dimensions are taken from the descriptors. Allows for pipelining of operations at runtime.

Parameters
[in]uploSpecifies whether the matrix A is upper or lower triangular:
  • PlasmaUpper: Upper triangle of A is stored;
  • PlasmaLower: Lower triangle of A is stored.
[in]ADescriptor of matrix A.
[in]BDescriptor of matrix B.
[in,out]XDescriptor of matrix X.
[out]AsDescriptor of auxiliary matrix A in single complex precision.
[out]XsDescriptor of auxiliary matrix X in single complex precision.
[out]RDescriptor of auxiliary remainder matrix R.
[out]workWorkspace needed to compute infinity norm of the matrix A.
[out]RnormWorkspace needed to store the max value in each of resudual vectors.
[out]XnormWorkspace needed to store the max value in each of currenct solution vectors.
[out]iterThe number of the iterations in the iterative refinement process, needed for the convergence. If failed, it is set to be -(1+itermax), where itermax = 30.
[in]sequenceIdentifies the sequence of function calls that this call belongs to (for completion checks and exception handling purposes).
[out]requestIdentifies this function call (for exception handling purposes).
Return values
voidErrors are returned by setting sequence->status and request->status to error values. The sequence->status and request->status should never be set to PLASMA_SUCCESS (the initial values) since another async call may be setting a failure value at the same time.
See also
plasma_dsposv
plasma_omp_dsposv
plasma_omp_zposv
int plasma_sposv ( plasma_enum_t  uplo,
int  n,
int  nrhs,
float *  pA,
int  lda,
float *  pB,
int  ldb 
)

Computes the solution to a 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 = L\times L^T, \]

if uplo = PlasmaLower, or

\[ A = U^T\times U, \]

if uplo = PlasmaUpper,

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
  • PlasmaUpper: Upper triangle of A is stored;
  • PlasmaLower: Lower triangle of A is stored.
[in]nThe number of linear equations, i.e., the order of the matrix A. n >= 0.
[in]nrhsThe number of right hand sides, i.e., the number of columns of the matrix B. nrhs >= 0.
[in,out]pAOn entry, the symmetric positive definite matrix A. If uplo = PlasmaUpper, 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 return value = 0, the factor U or L from the Cholesky factorization A = U^T*U or A = L*L^T.
[in]ldaThe leading dimension of the array A. lda >= max(1,n).
[in,out]pBOn entry, the n-by-nrhs right hand side matrix B. On exit, if return value = 0, the n-by-nrhs solution matrix X.
[in]ldbThe leading dimension of the array B. ldb >= max(1,n).
Return values
PlasmaSuccesssuccessful exit
<0 if -i, the i-th argument had an illegal value
>0 if 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.
See also
plasma_omp_sposv
plasma_cposv
plasma_dposv
plasma_sposv
void plasma_omp_sposv ( plasma_enum_t  uplo,
plasma_desc_t  A,
plasma_desc_t  B,
plasma_sequence_t *  sequence,
plasma_request_t *  request 
)

Solves a symmetric positive definite system of linear equations using Cholesky factorization. Non-blocking tile version of plasma_sposv(). Operates on matrices stored by tiles. All matrices are passed through descriptors. All dimensions are taken from the descriptors. Allows for pipelining of operations at runtime.

Parameters
[in]uplo
  • PlasmaUpper: Upper triangle of A is stored;
  • PlasmaLower: Lower triangle of A is stored.
[in,out]AOn entry, the symmetric positive definite matrix A. If uplo = PlasmaUpper, 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 return value = 0, the factor U or L from the Cholesky factorization A = U^T*U or A = L*L^T.
[in,out]BOn entry, the n-by-nrhs right hand side matrix B. On exit, if return value = 0, the n-by-nrhs solution matrix X.
[in]sequenceIdentifies the sequence of function calls that this call belongs to (for completion checks and exception handling purposes). Check the sequence->status for errors.
[out]requestIdentifies this function call (for exception handling purposes).
Return values
voidErrors are returned by setting sequence->status and request->status to error values. The sequence->status and request->status should never be set to PlasmaSuccess (the initial values) since another async call may be setting a failure value at the same time.
See also
plasma_sposv
plasma_omp_cposv
plasma_omp_dposv
plasma_omp_sposv
int plasma_zcposv ( plasma_enum_t  uplo,
int  n,
int  nrhs,
plasma_complex64_t *  pA,
int  lda,
plasma_complex64_t *  pB,
int  ldb,
plasma_complex64_t *  pX,
int  ldx,
int *  iter 
)

Computes the solution to a 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.

plasma_zcposv first factorizes the matrix using plasma_cpotrf and uses this factorization within an iterative refinement procedure to produce a solution with COMPLEX*16 normwise backward error quality (see below). If the approach fails the method falls back to a COMPLEX*16 factorization and solve.

The iterative refinement is not going to be a winning strategy if the ratio COMPLEX performance over COMPLEX*16 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: Rnorm < sqrt(n)*Xnorm*Anorm*eps*BWDmax where:

  • iter is the number of the current iteration in the iterative refinement process
  • Rnorm is the Infinity-norm of the residual
  • Xnorm is the Infinity-norm of the solution
  • Anorm is the Infinity-operator-norm of the matrix A
  • eps is the machine epsilon returned by DLAMCH('Epsilon'). The values itermax and BWDmax are fixed to 30 and 1.0D+00 respectively.
Parameters
[in]uploSpecifies whether the matrix A is upper or lower triangular:
  • PlasmaUpper: Upper triangle of A is stored;
  • PlasmaLower: Lower triangle of A is stored.
[in]nThe number of linear equations, i.e., the order of the matrix A. n >= 0.
[in]nrhsThe number of right hand sides, i.e., the number of columns of the matrix B. nrhs >= 0.
[in,out]pAThe n-by-n Hermitian positive definite coefficient matrix A. If uplo = PlasmaUpper, 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 = PlasmaLower, 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, contains the lower Cholesky factor matrix L, if uplo == PlasmaLower and upper Cholesky factor conj(L^T), otherwise.
[in]ldaThe leading dimension of the array A. lda >= max(1,n).
[in]pBThe n-by-nrhs matrix of right hand side matrix B. This matrix remains unchanged.
[in]ldbThe leading dimension of the array B. ldb >= max(1,n).
[out]pXIf return value = 0, the n-by-nrhs solution matrix X.
[in]ldxThe leading dimension of the array X. ldx >= max(1,n).
[out]iterThe number of the iterations in the iterative refinement process, needed for the convergence. If failed, it is set to be -(1+itermax), where itermax = 30.
Return values
PlasmaSuccesssuccessful exit
See also
plasma_omp_zcposv
plasma_dsposv
plasma_zposv
void plasma_omp_zcposv ( plasma_enum_t  uplo,
plasma_desc_t  A,
plasma_desc_t  B,
plasma_desc_t  X,
plasma_desc_t  As,
plasma_desc_t  Xs,
plasma_desc_t  R,
double *  work,
double *  Rnorm,
double *  Xnorm,
int *  iter,
plasma_sequence_t *  sequence,
plasma_request_t *  request 
)

Solves a Hermitian positive definite system using iterative refinement with the Cholesky factor computed using plasma_cpotrf. Non-blocking tile version of plasma_zcposv(). Operates on matrices stored by tiles. All matrices are passed through descriptors. All dimensions are taken from the descriptors. Allows for pipelining of operations at runtime.

Parameters
[in]uploSpecifies whether the matrix A is upper or lower triangular:
  • PlasmaUpper: Upper triangle of A is stored;
  • PlasmaLower: Lower triangle of A is stored.
[in]ADescriptor of matrix A.
[in]BDescriptor of matrix B.
[in,out]XDescriptor of matrix X.
[out]AsDescriptor of auxiliary matrix A in single complex precision.
[out]XsDescriptor of auxiliary matrix X in single complex precision.
[out]RDescriptor of auxiliary remainder matrix R.
[out]workWorkspace needed to compute infinity norm of the matrix A.
[out]RnormWorkspace needed to store the max value in each of resudual vectors.
[out]XnormWorkspace needed to store the max value in each of currenct solution vectors.
[out]iterThe number of the iterations in the iterative refinement process, needed for the convergence. If failed, it is set to be -(1+itermax), where itermax = 30.
[in]sequenceIdentifies the sequence of function calls that this call belongs to (for completion checks and exception handling purposes).
[out]requestIdentifies this function call (for exception handling purposes).
Return values
voidErrors are returned by setting sequence->status and request->status to error values. The sequence->status and request->status should never be set to PLASMA_SUCCESS (the initial values) since another async call may be setting a failure value at the same time.
See also
plasma_zcposv
plasma_omp_dsposv
plasma_omp_zposv
int plasma_zposv ( plasma_enum_t  uplo,
int  n,
int  nrhs,
plasma_complex64_t *  pA,
int  lda,
plasma_complex64_t *  pB,
int  ldb 
)

Computes the solution to a 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 = L\times L^H, \]

if uplo = PlasmaLower, or

\[ A = U^H\times U, \]

if uplo = PlasmaUpper,

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
  • PlasmaUpper: Upper triangle of A is stored;
  • PlasmaLower: Lower triangle of A is stored.
[in]nThe number of linear equations, i.e., the order of the matrix A. n >= 0.
[in]nrhsThe number of right hand sides, i.e., the number of columns of the matrix B. nrhs >= 0.
[in,out]pAOn entry, the Hermitian positive definite matrix A. If uplo = PlasmaUpper, 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 return value = 0, the factor U or L from the Cholesky factorization A = U^H*U or A = L*L^H.
[in]ldaThe leading dimension of the array A. lda >= max(1,n).
[in,out]pBOn entry, the n-by-nrhs right hand side matrix B. On exit, if return value = 0, the n-by-nrhs solution matrix X.
[in]ldbThe leading dimension of the array B. ldb >= max(1,n).
Return values
PlasmaSuccesssuccessful exit
<0 if -i, the i-th argument had an illegal value
>0 if 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.
See also
plasma_omp_zposv
plasma_cposv
plasma_dposv
plasma_sposv
void plasma_omp_zposv ( plasma_enum_t  uplo,
plasma_desc_t  A,
plasma_desc_t  B,
plasma_sequence_t *  sequence,
plasma_request_t *  request 
)

Solves a Hermitian positive definite system of linear equations using Cholesky factorization. Non-blocking tile version of plasma_zposv(). Operates on matrices stored by tiles. All matrices are passed through descriptors. All dimensions are taken from the descriptors. Allows for pipelining of operations at runtime.

Parameters
[in]uplo
  • PlasmaUpper: Upper triangle of A is stored;
  • PlasmaLower: Lower triangle of A is stored.
[in,out]AOn entry, the Hermitian positive definite matrix A. If uplo = PlasmaUpper, 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 return value = 0, the factor U or L from the Cholesky factorization A = U^H*U or A = L*L^H.
[in,out]BOn entry, the n-by-nrhs right hand side matrix B. On exit, if return value = 0, the n-by-nrhs solution matrix X.
[in]sequenceIdentifies the sequence of function calls that this call belongs to (for completion checks and exception handling purposes). Check the sequence->status for errors.
[out]requestIdentifies this function call (for exception handling purposes).
Return values
voidErrors are returned by setting sequence->status and request->status to error values. The sequence->status and request->status should never be set to PlasmaSuccess (the initial values) since another async call may be setting a failure value at the same time.
See also
plasma_zposv
plasma_omp_cposv
plasma_omp_dposv
plasma_omp_sposv